energymodule.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
29#ifndef EWOMS_ENERGY_MODULE_HH
30#define EWOMS_ENERGY_MODULE_HH
31
32#include <dune/common/fvector.hh>
33
34#include <opm/material/common/Valgrind.hpp>
35
38
40
41#include <string>
42
43namespace Opm {
49template <class TypeTag, bool enableEnergy>
51
55template <class TypeTag>
56class EnergyModule<TypeTag, /*enableEnergy=*/false>
57{
66
67 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
68
69 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
70
71public:
75 static void registerParameters()
76 {}
77
83 static std::string primaryVarName(unsigned)
84 { return ""; }
85
91 static std::string eqName(unsigned)
92 { return ""; }
93
98 static Scalar primaryVarWeight(const Model&,
99 unsigned,
100 unsigned)
101 { return -1; }
102
106 static Scalar eqWeight(const Model&,
107 unsigned,
108 unsigned)
109 { return -1; }
110
114 template <class FluidState>
115 static void setPriVarTemperatures(PrimaryVariables&,
116 const FluidState&)
117 {}
118
123 template <class RateVector, class FluidState>
124 static void setEnthalpyRate(RateVector&,
125 const FluidState&,
126 unsigned,
127 const Evaluation&)
128 {}
129
133 static void setEnthalpyRate(RateVector&,
134 const Evaluation&)
135 {}
136
140 static void addToEnthalpyRate(RateVector&,
141 const Evaluation&)
142 {}
143
147 static Scalar thermalConductionRate(const ExtensiveQuantities&)
148 { return 0.0; }
149
154 template <class LhsEval>
155 static void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>&,
156 const IntensiveQuantities&,
157 unsigned)
158 {}
159
164 template <class LhsEval, class Scv>
165 static void addFracturePhaseStorage(Dune::FieldVector<LhsEval, numEq>&,
166 const IntensiveQuantities&,
167 const Scv&,
168 unsigned)
169 {}
170
175 template <class LhsEval>
176 static void addSolidEnergyStorage(Dune::FieldVector<LhsEval, numEq>&,
177 const IntensiveQuantities&)
178 {}
179
186 template <class Context>
187 static void addAdvectiveFlux(RateVector&,
188 const Context&,
189 unsigned,
190 unsigned)
191 {}
192
198 template <class Context>
199 static void handleFractureFlux(RateVector&,
200 const Context&,
201 unsigned,
202 unsigned)
203 {}
204
211 template <class Context>
212 static void addDiffusiveFlux(RateVector&,
213 const Context&,
214 unsigned,
215 unsigned)
216 {}
217};
218
222template <class TypeTag>
223class EnergyModule<TypeTag, /*enableEnergy=*/true>
224{
235
236 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
237 enum { numPhases = FluidSystem::numPhases };
238 enum { energyEqIdx = Indices::energyEqIdx };
239 enum { temperatureIdx = Indices::temperatureIdx };
240
241 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
242 using Toolbox = Opm::MathToolbox<Evaluation>;
243
244public:
248 static void registerParameters()
249 {}
250
256 static std::string primaryVarName(unsigned pvIdx)
257 {
258 if (pvIdx == temperatureIdx)
259 return "temperature";
260 return "";
261 }
262
268 static std::string eqName(unsigned eqIdx)
269 {
270 if (eqIdx == energyEqIdx)
271 return "energy";
272 return "";
273 }
274
279 static Scalar primaryVarWeight(const Model& model, unsigned globalDofIdx, unsigned pvIdx)
280 {
281 if (pvIdx != temperatureIdx)
282 return -1;
283
284 // make the weight of the temperature primary variable inversly proportional to its value
285 return std::max(1.0/1000, 1.0/model.solution(/*timeIdx=*/0)[globalDofIdx][temperatureIdx]);
286 }
287
291 static Scalar eqWeight(const Model&,
292 unsigned,
293 unsigned eqIdx)
294 {
295 if (eqIdx != energyEqIdx)
296 return -1;
297
298 // approximate change of internal energy of 1kg of liquid water for a temperature
299 // change of 30K
300 return 1.0 / (4.184e3 * 30.0);
301 }
302
306 static void setEnthalpyRate(RateVector& rateVec, const Evaluation& rate)
307 { rateVec[energyEqIdx] = rate; }
308
312 static void addToEnthalpyRate(RateVector& rateVec, const Evaluation& rate)
313 { rateVec[energyEqIdx] += rate; }
314
318 static Evaluation thermalConductionRate(const ExtensiveQuantities& extQuants)
319 { return -extQuants.temperatureGradNormal() * extQuants.thermalConductivity(); }
320
325 template <class RateVector, class FluidState>
326 static void setEnthalpyRate(RateVector& rateVec,
327 const FluidState& fluidState,
328 unsigned phaseIdx,
329 const Evaluation& volume)
330 {
331 rateVec[energyEqIdx] =
332 volume
333 * fluidState.density(phaseIdx)
334 * fluidState.enthalpy(phaseIdx);
335 }
336
340 template <class FluidState>
341 static void setPriVarTemperatures(PrimaryVariables& priVars,
342 const FluidState& fs)
343 {
344 priVars[temperatureIdx] = Toolbox::value(fs.temperature(/*phaseIdx=*/0));
345#ifndef NDEBUG
346 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
347 assert(std::abs(Toolbox::value(fs.temperature(/*phaseIdx=*/0))
348 - Toolbox::value(fs.temperature(phaseIdx))) < 1e-30);
349 }
350#endif
351 }
352
357 template <class LhsEval>
358 static void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
359 const IntensiveQuantities& intQuants,
360 unsigned phaseIdx)
361 {
362 const auto& fs = intQuants.fluidState();
363 storage[energyEqIdx] +=
364 Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
365 * Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx))
366 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
367 * Toolbox::template decay<LhsEval>(intQuants.porosity());
368 }
369
374 template <class Scv, class LhsEval>
375 static void addFracturePhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
376 const IntensiveQuantities& intQuants,
377 const Scv& scv,
378 unsigned phaseIdx)
379 {
380 const auto& fs = intQuants.fractureFluidState();
381 storage[energyEqIdx] +=
382 Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
383 * Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx))
384 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
385 * Toolbox::template decay<LhsEval>(intQuants.fracturePorosity())
386 * Toolbox::template decay<LhsEval>(intQuants.fractureVolume())/scv.volume();
387 }
388
393 template <class LhsEval>
394 static void addSolidEnergyStorage(Dune::FieldVector<LhsEval, numEq>& storage,
395 const IntensiveQuantities& intQuants)
396 { storage[energyEqIdx] += Opm::decay<LhsEval>(intQuants.solidInternalEnergy()); }
397
404 template <class Context>
405 static void addAdvectiveFlux(RateVector& flux,
406 const Context& context,
407 unsigned spaceIdx,
408 unsigned timeIdx)
409 {
410 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
411
412 // advective energy flux in all phases
413 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
414 if (!context.model().phaseIsConsidered(phaseIdx))
415 continue;
416
417 // intensive quantities of the upstream and the downstream DOFs
418 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
419 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
420
421 flux[energyEqIdx] +=
422 extQuants.volumeFlux(phaseIdx)
423 * up.fluidState().enthalpy(phaseIdx)
424 * up.fluidState().density(phaseIdx);
425 }
426 }
427
433 template <class Context>
434 static void handleFractureFlux(RateVector& flux,
435 const Context& context,
436 unsigned spaceIdx,
437 unsigned timeIdx)
438 {
439 const auto& scvf = context.stencil(timeIdx).interiorFace(spaceIdx);
440 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
441
442 // reduce the energy flux in the matrix by the half the width occupied by the
443 // fracture
444 flux[energyEqIdx] *=
445 1 - extQuants.fractureWidth()/(2*scvf.area());
446
447 // advective energy flux in all phases
448 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
449 if (!context.model().phaseIsConsidered(phaseIdx))
450 continue;
451
452 // intensive quantities of the upstream and the downstream DOFs
453 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
454 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
455
456 flux[energyEqIdx] +=
457 extQuants.fractureVolumeFlux(phaseIdx)
458 * up.fluidState().enthalpy(phaseIdx)
459 * up.fluidState().density(phaseIdx);
460 }
461 }
462
469 template <class Context>
470 static void addDiffusiveFlux(RateVector& flux,
471 const Context& context,
472 unsigned spaceIdx,
473 unsigned timeIdx)
474 {
475 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
476
477 // diffusive energy flux
478 flux[energyEqIdx] +=
479 - extQuants.temperatureGradNormal()
480 * extQuants.thermalConductivity();
481 }
482};
483
490template <unsigned PVOffset, bool enableEnergy>
492
496template <unsigned PVOffset>
497struct EnergyIndices<PVOffset, /*enableEnergy=*/false>
498{
500 enum { temperatureIdx = -1000 };
501
503 enum { energyEqIdx = -1000 };
504
505protected:
506 enum { numEq_ = 0 };
507};
508
512template <unsigned PVOffset>
513struct EnergyIndices<PVOffset, /*enableEnergy=*/true>
514{
516 enum { temperatureIdx = PVOffset };
517
519 enum { energyEqIdx = PVOffset };
520
521protected:
522 enum { numEq_ = 1 };
523};
524
531template <class TypeTag, bool enableEnergy>
533
537template <class TypeTag>
538class EnergyIntensiveQuantities<TypeTag, /*enableEnergy=*/false>
539{
544
545 using Toolbox = Opm::MathToolbox<Evaluation>;
546
547public:
552 Evaluation solidInternalEnergy() const
553 {
554 throw std::logic_error("solidInternalEnergy() does not make sense for isothermal models");
555 }
556
561 Evaluation thermalConductivity() const
562 {
563 throw std::logic_error("thermalConductivity() does not make sense for isothermal models");
564 }
565
566protected:
570 template <class FluidState, class Context>
571 static void updateTemperatures_(FluidState& fluidState,
572 const Context& context,
573 unsigned spaceIdx,
574 unsigned timeIdx)
575 {
576 Scalar T = context.problem().temperature(context, spaceIdx, timeIdx);
577 fluidState.setTemperature(Toolbox::createConstant(T));
578 }
579
584 template <class FluidState>
585 void update_(FluidState&,
586 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
587 const ElementContext&,
588 unsigned,
589 unsigned)
590 { }
591};
592
596template <class TypeTag>
597class EnergyIntensiveQuantities<TypeTag, /*enableEnergy=*/true>
598{
606
607 enum { numPhases = FluidSystem::numPhases };
608 enum { energyEqIdx = Indices::energyEqIdx };
609 enum { temperatureIdx = Indices::temperatureIdx };
610
611 using Toolbox = Opm::MathToolbox<Evaluation>;
612
613protected:
617 template <class FluidState, class Context>
618 static void updateTemperatures_(FluidState& fluidState,
619 const Context& context,
620 unsigned spaceIdx,
621 unsigned timeIdx)
622 {
623 const auto& priVars = context.primaryVars(spaceIdx, timeIdx);
624 Evaluation val;
625 if (std::is_same<Evaluation, Scalar>::value) // finite differences
626 val = Toolbox::createConstant(priVars[temperatureIdx]);
627 else {
628 // automatic differentiation
629 if (timeIdx == 0)
630 val = Toolbox::createVariable(priVars[temperatureIdx], temperatureIdx);
631 else
632 val = Toolbox::createConstant(priVars[temperatureIdx]);
633 }
634 fluidState.setTemperature(val);
635 }
636
641 template <class FluidState>
642 void update_(FluidState& fs,
643 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
644 const ElementContext& elemCtx,
645 unsigned dofIdx,
646 unsigned timeIdx)
647 {
648 // set the specific enthalpies of the fluids
649 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
650 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
651 continue;
652
653 fs.setEnthalpy(phaseIdx,
654 FluidSystem::enthalpy(fs, paramCache, phaseIdx));
655 }
656
657 // compute and set the volumetric internal energy of the solid phase
658 const auto& problem = elemCtx.problem();
659 const auto& solidEnergyParams = problem.solidEnergyLawParams(elemCtx, dofIdx, timeIdx);
660 const auto& thermalCondParams = problem.thermalConductionLawParams(elemCtx, dofIdx, timeIdx);
661
662 solidInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyParams, fs);
663 thermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalCondParams, fs);
664
665 Opm::Valgrind::CheckDefined(solidInternalEnergy_);
666 Opm::Valgrind::CheckDefined(thermalConductivity_);
667 }
668
669public:
674 const Evaluation& solidInternalEnergy() const
675 { return solidInternalEnergy_; }
676
681 const Evaluation& thermalConductivity() const
682 { return thermalConductivity_; }
683
684private:
685 Evaluation solidInternalEnergy_;
686 Evaluation thermalConductivity_;
687};
688
695template <class TypeTag, bool enableEnergy>
697
701template <class TypeTag>
702class EnergyExtensiveQuantities<TypeTag, /*enableEnergy=*/false>
703{
706
707protected:
712 void update_(const ElementContext&,
713 unsigned,
714 unsigned)
715 {}
716
717 template <class Context, class FluidState>
718 void updateBoundary_(const Context&,
719 unsigned,
720 unsigned,
721 const FluidState&)
722 {}
723
724public:
729 {
730 throw std::logic_error("Calling temperatureGradNormal() does not make sense "
731 "for isothermal models");
732 }
733
737 Scalar thermalConductivity() const
738 {
739 throw std::logic_error("Calling thermalConductivity() does not make sense for "
740 "isothermal models");
741 }
742};
743
747template <class TypeTag>
748class EnergyExtensiveQuantities<TypeTag, /*enableEnergy=*/true>
749{
754
755 enum { dimWorld = GridView::dimensionworld };
756 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
757 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
758
759protected:
764 void update_(const ElementContext& elemCtx, unsigned faceIdx, unsigned timeIdx)
765 {
766 const auto& gradCalc = elemCtx.gradientCalculator();
767 Opm::TemperatureCallback<TypeTag> temperatureCallback(elemCtx);
768
769 EvalDimVector temperatureGrad;
770 gradCalc.calculateGradient(temperatureGrad,
771 elemCtx,
772 faceIdx,
773 temperatureCallback);
774
775 // scalar product of temperature gradient and scvf normal
776 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(faceIdx);
777
778 temperatureGradNormal_ = 0;
779 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
780 temperatureGradNormal_ += (face.normal()[dimIdx]*temperatureGrad[dimIdx]);
781
782 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
783 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
784 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
785
786 // arithmetic mean
787 thermalConductivity_ =
788 0.5 * (intQuantsInside.thermalConductivity() + intQuantsOutside.thermalConductivity());
789 Opm::Valgrind::CheckDefined(thermalConductivity_);
790 }
791
792 template <class Context, class FluidState>
793 void updateBoundary_(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fs)
794 {
795 const auto& stencil = context.stencil(timeIdx);
796 const auto& face = stencil.boundaryFace(bfIdx);
797
798 const auto& elemCtx = context.elementContext();
799 unsigned insideScvIdx = face.interiorIndex();
800 const auto& insideScv = elemCtx.stencil(timeIdx).subControlVolume(insideScvIdx);
801
802 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
803 const auto& fsInside = intQuantsInside.fluidState();
804
805 // distance between the center of the SCV and center of the boundary face
806 DimVector distVec = face.integrationPos();
807 distVec -= insideScv.geometry().center();
808
809 Scalar tmp = 0;
810 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
811 tmp += distVec[dimIdx] * face.normal()[dimIdx];
812 Scalar dist = tmp;
813
814 // if the following assertation triggers, the center of the
815 // center of the interior SCV was not inside the element!
816 assert(dist > 0);
817
818 // calculate the temperature gradient using two-point gradient
819 // appoximation
820 temperatureGradNormal_ =
821 (fs.temperature(/*phaseIdx=*/0) - fsInside.temperature(/*phaseIdx=*/0)) / dist;
822
823 // take the value for thermal conductivity from the interior finite volume
824 thermalConductivity_ = intQuantsInside.thermalConductivity();
825 }
826
827public:
831 const Evaluation& temperatureGradNormal() const
832 { return temperatureGradNormal_; }
833
838 const Evaluation& thermalConductivity() const
839 { return thermalConductivity_; }
840
841private:
842 Evaluation temperatureGradNormal_;
843 Evaluation thermalConductivity_;
844};
845
846} // namespace Opm
847
848#endif
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: energymodule.hh:718
Scalar temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition: energymodule.hh:728
Scalar thermalConductivity() const
The total thermal conductivity at the face .
Definition: energymodule.hh:737
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:712
const Evaluation & temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition: energymodule.hh:831
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:764
const Evaluation & thermalConductivity() const
The total thermal conductivity at the face .
Definition: energymodule.hh:838
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fs)
Definition: energymodule.hh:793
Provides the quantities required to calculate energy fluxes.
Definition: energymodule.hh:696
Evaluation thermalConductivity() const
Returns the total thermal conductivity of the solid matrix in the sub-control volume.
Definition: energymodule.hh:561
static void updateTemperatures_(FluidState &fluidState, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Update the temperatures of the fluids of a fluid state.
Definition: energymodule.hh:571
void update_(FluidState &, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &, const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:585
Evaluation solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition: energymodule.hh:552
const Evaluation & solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition: energymodule.hh:674
const Evaluation & thermalConductivity() const
Returns the total conductivity capacity of the solid matrix in the sub-control volume.
Definition: energymodule.hh:681
void update_(FluidState &fs, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > &paramCache, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:642
static void updateTemperatures_(FluidState &fluidState, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Update the temperatures of the fluids of a fluid state.
Definition: energymodule.hh:618
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:532
static Scalar eqWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a equation of the residual.
Definition: energymodule.hh:106
static void setPriVarTemperatures(PrimaryVariables &, const FluidState &)
Given a fluid state, set the temperature in the primary variables.
Definition: energymodule.hh:115
static void handleFractureFlux(RateVector &, const Context &, unsigned, unsigned)
Evaluates the advective energy fluxes over a fracture which should be attributed to a face of a subco...
Definition: energymodule.hh:199
static void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &, const IntensiveQuantities &, unsigned)
Add the energy storage term for a fluid phase to an equation vector.
Definition: energymodule.hh:155
static void addToEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:140
static void addDiffusiveFlux(RateVector &, const Context &, unsigned, unsigned)
Adds the diffusive energy flux to the flux vector over the face of a sub-control volume.
Definition: energymodule.hh:212
static void addSolidEnergyStorage(Dune::FieldVector< LhsEval, numEq > &, const IntensiveQuantities &)
Add the energy storage term for the fracture part a fluid phase to an equation vector.
Definition: energymodule.hh:176
static void setEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:133
static void addAdvectiveFlux(RateVector &, const Context &, unsigned, unsigned)
Evaluates the advective energy fluxes over a face of a subcontrol volume and adds the result in the f...
Definition: energymodule.hh:187
static void addFracturePhaseStorage(Dune::FieldVector< LhsEval, numEq > &, const IntensiveQuantities &, const Scv &, unsigned)
Add the energy storage term for a fluid phase to an equation vector.
Definition: energymodule.hh:165
static std::string primaryVarName(unsigned)
Returns the name of a primary variable or an empty string if the specified primary variable index doe...
Definition: energymodule.hh:83
static std::string eqName(unsigned)
Returns the name of an equation or an empty string if the specified equation index does not belong to...
Definition: energymodule.hh:91
static Scalar thermalConductionRate(const ExtensiveQuantities &)
Add the rate of the conductive energy flux to a rate vector.
Definition: energymodule.hh:147
static void setEnthalpyRate(RateVector &, const FluidState &, unsigned, const Evaluation &)
Given a fluid state, set the enthalpy rate which emerges from a volumetric rate.
Definition: energymodule.hh:124
static Scalar primaryVarWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a primary variable for calculating relative errors.
Definition: energymodule.hh:98
static void registerParameters()
Register all run-time parameters for the energy module.
Definition: energymodule.hh:75
static void handleFractureFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Evaluates the advective energy fluxes over a fracture which should be attributed to a face of a subco...
Definition: energymodule.hh:434
static void addToEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:312
static void addFracturePhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants, const Scv &scv, unsigned phaseIdx)
Add the energy storage term for a fluid phase to an equation vector.
Definition: energymodule.hh:375
static void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants, unsigned phaseIdx)
Add the energy storage term for a fluid phase to an equation vector.
Definition: energymodule.hh:358
static void setEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Set the rate of energy flux of a rate vector.
Definition: energymodule.hh:306
static void setPriVarTemperatures(PrimaryVariables &priVars, const FluidState &fs)
Given a fluid state, set the temperature in the primary variables.
Definition: energymodule.hh:341
static void addDiffusiveFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Adds the diffusive energy flux to the flux vector over the face of a sub-control volume.
Definition: energymodule.hh:470
static void setEnthalpyRate(RateVector &rateVec, const FluidState &fluidState, unsigned phaseIdx, const Evaluation &volume)
Given a fluid state, set the enthalpy rate which emerges from a volumetric rate.
Definition: energymodule.hh:326
static std::string eqName(unsigned eqIdx)
Returns the name of an equation or an empty string if the specified equation index does not belong to...
Definition: energymodule.hh:268
static void registerParameters()
Register all run-time parameters for the energy module.
Definition: energymodule.hh:248
static void addSolidEnergyStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Add the energy storage term for a fluid phase to an equation vector.
Definition: energymodule.hh:394
static void addAdvectiveFlux(RateVector &flux, const Context &context, unsigned spaceIdx, unsigned timeIdx)
Evaluates the advective energy fluxes for a flux integration point and adds the result in the flux ve...
Definition: energymodule.hh:405
static std::string primaryVarName(unsigned pvIdx)
Returns the name of a primary variable or an empty string if the specified primary variable index doe...
Definition: energymodule.hh:256
static Scalar primaryVarWeight(const Model &model, unsigned globalDofIdx, unsigned pvIdx)
Returns the relative weight of a primary variable for calculating relative errors.
Definition: energymodule.hh:279
static Scalar eqWeight(const Model &, unsigned, unsigned eqIdx)
Returns the relative weight of a equation.
Definition: energymodule.hh:291
static Evaluation thermalConductionRate(const ExtensiveQuantities &extQuants)
Returns the conductive energy flux for a given flux integration point.
Definition: energymodule.hh:318
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Callback class for temperature.
Definition: quantitycallbacks.hh:48
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
This method contains all callback classes for quantities that are required by some extensive quantiti...
Provides the indices required for the energy equation.
Definition: energymodule.hh:491