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
34
35#include <opm/material/common/Valgrind.hpp>
36
37#include <dune/common/fvector.hh>
38
39#include <string>
40
41namespace Opm {
47template <class TypeTag, bool enableEnergy>
49
53template <class TypeTag>
54class EnergyModule<TypeTag, /*enableEnergy=*/false>
55{
64
65 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
66
67 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
68
69public:
73 static void registerParameters()
74 {}
75
81 static std::string primaryVarName(unsigned)
82 { return ""; }
83
89 static std::string eqName(unsigned)
90 { return ""; }
91
96 static Scalar primaryVarWeight(const Model&,
97 unsigned,
98 unsigned)
99 { return -1; }
100
104 static Scalar eqWeight(const Model&,
105 unsigned,
106 unsigned)
107 { return -1; }
108
112 template <class FluidState>
113 static void setPriVarTemperatures(PrimaryVariables&,
114 const FluidState&)
115 {}
116
121 template <class RateVector, class FluidState>
122 static void setEnthalpyRate(RateVector&,
123 const FluidState&,
124 unsigned,
125 const Evaluation&)
126 {}
127
131 static void setEnthalpyRate(RateVector&,
132 const Evaluation&)
133 {}
134
138 static void addToEnthalpyRate(RateVector&,
139 const Evaluation&)
140 {}
141
145 static Scalar thermalConductionRate(const ExtensiveQuantities&)
146 { return 0.0; }
147
152 template <class LhsEval>
153 static void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>&,
154 const IntensiveQuantities&,
155 unsigned)
156 {}
157
162 template <class LhsEval, class Scv>
163 static void addFracturePhaseStorage(Dune::FieldVector<LhsEval, numEq>&,
164 const IntensiveQuantities&,
165 const Scv&,
166 unsigned)
167 {}
168
173 template <class LhsEval>
174 static void addSolidEnergyStorage(Dune::FieldVector<LhsEval, numEq>&,
175 const IntensiveQuantities&)
176 {}
177
184 template <class Context>
185 static void addAdvectiveFlux(RateVector&,
186 const Context&,
187 unsigned,
188 unsigned)
189 {}
190
196 template <class Context>
197 static void handleFractureFlux(RateVector&,
198 const Context&,
199 unsigned,
200 unsigned)
201 {}
202
209 template <class Context>
210 static void addDiffusiveFlux(RateVector&,
211 const Context&,
212 unsigned,
213 unsigned)
214 {}
215};
216
220template <class TypeTag>
221class EnergyModule<TypeTag, /*enableEnergy=*/true>
222{
233
234 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
235 enum { numPhases = FluidSystem::numPhases };
236 enum { energyEqIdx = Indices::energyEqIdx };
237 enum { temperatureIdx = Indices::temperatureIdx };
238
239 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
240 using Toolbox = Opm::MathToolbox<Evaluation>;
241
242public:
246 static void registerParameters()
247 {}
248
254 static std::string primaryVarName(unsigned pvIdx)
255 {
256 if (pvIdx == temperatureIdx)
257 return "temperature";
258 return "";
259 }
260
266 static std::string eqName(unsigned eqIdx)
267 {
268 if (eqIdx == energyEqIdx)
269 return "energy";
270 return "";
271 }
272
277 static Scalar primaryVarWeight(const Model& model, unsigned globalDofIdx, unsigned pvIdx)
278 {
279 if (pvIdx != temperatureIdx)
280 return -1;
281
282 // make the weight of the temperature primary variable inversly proportional to its value
283 return std::max(1.0/1000, 1.0/model.solution(/*timeIdx=*/0)[globalDofIdx][temperatureIdx]);
284 }
285
289 static Scalar eqWeight(const Model&,
290 unsigned,
291 unsigned eqIdx)
292 {
293 if (eqIdx != energyEqIdx)
294 return -1;
295
296 // approximate change of internal energy of 1kg of liquid water for a temperature
297 // change of 30K
298 return 1.0 / (4.184e3 * 30.0);
299 }
300
304 static void setEnthalpyRate(RateVector& rateVec, const Evaluation& rate)
305 { rateVec[energyEqIdx] = rate; }
306
310 static void addToEnthalpyRate(RateVector& rateVec, const Evaluation& rate)
311 { rateVec[energyEqIdx] += rate; }
312
316 static Evaluation thermalConductionRate(const ExtensiveQuantities& extQuants)
317 { return -extQuants.temperatureGradNormal() * extQuants.thermalConductivity(); }
318
323 template <class RateVector, class FluidState>
324 static void setEnthalpyRate(RateVector& rateVec,
325 const FluidState& fluidState,
326 unsigned phaseIdx,
327 const Evaluation& volume)
328 {
329 rateVec[energyEqIdx] =
330 volume
331 * fluidState.density(phaseIdx)
332 * fluidState.enthalpy(phaseIdx);
333 }
334
338 template <class FluidState>
339 static void setPriVarTemperatures(PrimaryVariables& priVars,
340 const FluidState& fs)
341 {
342 priVars[temperatureIdx] = Toolbox::value(fs.temperature(/*phaseIdx=*/0));
343#ifndef NDEBUG
344 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
345 assert(std::abs(Toolbox::value(fs.temperature(/*phaseIdx=*/0))
346 - Toolbox::value(fs.temperature(phaseIdx))) < 1e-30);
347 }
348#endif
349 }
350
355 template <class LhsEval>
356 static void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
357 const IntensiveQuantities& intQuants,
358 unsigned phaseIdx)
359 {
360 const auto& fs = intQuants.fluidState();
361 storage[energyEqIdx] +=
362 Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
363 * Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx))
364 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
365 * Toolbox::template decay<LhsEval>(intQuants.porosity());
366 }
367
372 template <class Scv, class LhsEval>
373 static void addFracturePhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
374 const IntensiveQuantities& intQuants,
375 const Scv& scv,
376 unsigned phaseIdx)
377 {
378 const auto& fs = intQuants.fractureFluidState();
379 storage[energyEqIdx] +=
380 Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
381 * Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx))
382 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
383 * Toolbox::template decay<LhsEval>(intQuants.fracturePorosity())
384 * Toolbox::template decay<LhsEval>(intQuants.fractureVolume())/scv.volume();
385 }
386
391 template <class LhsEval>
392 static void addSolidEnergyStorage(Dune::FieldVector<LhsEval, numEq>& storage,
393 const IntensiveQuantities& intQuants)
394 { storage[energyEqIdx] += Opm::decay<LhsEval>(intQuants.solidInternalEnergy()); }
395
402 template <class Context>
403 static void addAdvectiveFlux(RateVector& flux,
404 const Context& context,
405 unsigned spaceIdx,
406 unsigned timeIdx)
407 {
408 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
409
410 // advective energy flux in all phases
411 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
412 if (!context.model().phaseIsConsidered(phaseIdx))
413 continue;
414
415 // intensive quantities of the upstream and the downstream DOFs
416 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
417 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
418
419 flux[energyEqIdx] +=
420 extQuants.volumeFlux(phaseIdx)
421 * up.fluidState().enthalpy(phaseIdx)
422 * up.fluidState().density(phaseIdx);
423 }
424 }
425
431 template <class Context>
432 static void handleFractureFlux(RateVector& flux,
433 const Context& context,
434 unsigned spaceIdx,
435 unsigned timeIdx)
436 {
437 const auto& scvf = context.stencil(timeIdx).interiorFace(spaceIdx);
438 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
439
440 // reduce the energy flux in the matrix by the half the width occupied by the
441 // fracture
442 flux[energyEqIdx] *=
443 1 - extQuants.fractureWidth()/(2*scvf.area());
444
445 // advective energy flux in all phases
446 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
447 if (!context.model().phaseIsConsidered(phaseIdx))
448 continue;
449
450 // intensive quantities of the upstream and the downstream DOFs
451 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
452 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
453
454 flux[energyEqIdx] +=
455 extQuants.fractureVolumeFlux(phaseIdx)
456 * up.fluidState().enthalpy(phaseIdx)
457 * up.fluidState().density(phaseIdx);
458 }
459 }
460
467 template <class Context>
468 static void addDiffusiveFlux(RateVector& flux,
469 const Context& context,
470 unsigned spaceIdx,
471 unsigned timeIdx)
472 {
473 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
474
475 // diffusive energy flux
476 flux[energyEqIdx] +=
477 - extQuants.temperatureGradNormal()
478 * extQuants.thermalConductivity();
479 }
480};
481
488template <unsigned PVOffset, bool enableEnergy>
490
494template <unsigned PVOffset>
495struct EnergyIndices<PVOffset, /*enableEnergy=*/false>
496{
498 enum { temperatureIdx = -1000 };
499
501 enum { energyEqIdx = -1000 };
502
503protected:
504 enum { numEq_ = 0 };
505};
506
510template <unsigned PVOffset>
511struct EnergyIndices<PVOffset, /*enableEnergy=*/true>
512{
514 enum { temperatureIdx = PVOffset };
515
517 enum { energyEqIdx = PVOffset };
518
519protected:
520 enum { numEq_ = 1 };
521};
522
529template <class TypeTag, bool enableEnergy>
531
535template <class TypeTag>
536class EnergyIntensiveQuantities<TypeTag, /*enableEnergy=*/false>
537{
542
543 using Toolbox = Opm::MathToolbox<Evaluation>;
544
545public:
550 Evaluation solidInternalEnergy() const
551 {
552 throw std::logic_error("solidInternalEnergy() does not make sense for isothermal models");
553 }
554
559 Evaluation thermalConductivity() const
560 {
561 throw std::logic_error("thermalConductivity() does not make sense for isothermal models");
562 }
563
564protected:
568 template <class FluidState, class Context>
569 static void updateTemperatures_(FluidState& fluidState,
570 const Context& context,
571 unsigned spaceIdx,
572 unsigned timeIdx)
573 {
574 Scalar T = context.problem().temperature(context, spaceIdx, timeIdx);
575 fluidState.setTemperature(Toolbox::createConstant(T));
576 }
577
582 template <class FluidState>
583 void update_(FluidState&,
584 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
585 const ElementContext&,
586 unsigned,
587 unsigned)
588 { }
589};
590
594template <class TypeTag>
595class EnergyIntensiveQuantities<TypeTag, /*enableEnergy=*/true>
596{
604
605 enum { numPhases = FluidSystem::numPhases };
606 enum { energyEqIdx = Indices::energyEqIdx };
607 enum { temperatureIdx = Indices::temperatureIdx };
608
609 using Toolbox = Opm::MathToolbox<Evaluation>;
610
611protected:
615 template <class FluidState, class Context>
616 static void updateTemperatures_(FluidState& fluidState,
617 const Context& context,
618 unsigned spaceIdx,
619 unsigned timeIdx)
620 {
621 const auto& priVars = context.primaryVars(spaceIdx, timeIdx);
622 Evaluation val;
623 if (std::is_same<Evaluation, Scalar>::value) // finite differences
624 val = Toolbox::createConstant(priVars[temperatureIdx]);
625 else {
626 // automatic differentiation
627 if (timeIdx == 0)
628 val = Toolbox::createVariable(priVars[temperatureIdx], temperatureIdx);
629 else
630 val = Toolbox::createConstant(priVars[temperatureIdx]);
631 }
632 fluidState.setTemperature(val);
633 }
634
639 template <class FluidState>
640 void update_(FluidState& fs,
641 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
642 const ElementContext& elemCtx,
643 unsigned dofIdx,
644 unsigned timeIdx)
645 {
646 // set the specific enthalpies of the fluids
647 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
648 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
649 continue;
650
651 fs.setEnthalpy(phaseIdx,
652 FluidSystem::enthalpy(fs, paramCache, phaseIdx));
653 }
654
655 // compute and set the volumetric internal energy of the solid phase
656 const auto& problem = elemCtx.problem();
657 const auto& solidEnergyParams = problem.solidEnergyLawParams(elemCtx, dofIdx, timeIdx);
658 const auto& thermalCondParams = problem.thermalConductionLawParams(elemCtx, dofIdx, timeIdx);
659
660 solidInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyParams, fs);
661 thermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalCondParams, fs);
662
663 Opm::Valgrind::CheckDefined(solidInternalEnergy_);
664 Opm::Valgrind::CheckDefined(thermalConductivity_);
665 }
666
667public:
672 const Evaluation& solidInternalEnergy() const
673 { return solidInternalEnergy_; }
674
679 const Evaluation& thermalConductivity() const
680 { return thermalConductivity_; }
681
682private:
683 Evaluation solidInternalEnergy_;
684 Evaluation thermalConductivity_;
685};
686
693template <class TypeTag, bool enableEnergy>
695
699template <class TypeTag>
700class EnergyExtensiveQuantities<TypeTag, /*enableEnergy=*/false>
701{
704
705protected:
710 void update_(const ElementContext&,
711 unsigned,
712 unsigned)
713 {}
714
715 template <class Context, class FluidState>
716 void updateBoundary_(const Context&,
717 unsigned,
718 unsigned,
719 const FluidState&)
720 {}
721
722public:
727 {
728 throw std::logic_error("Calling temperatureGradNormal() does not make sense "
729 "for isothermal models");
730 }
731
735 Scalar thermalConductivity() const
736 {
737 throw std::logic_error("Calling thermalConductivity() does not make sense for "
738 "isothermal models");
739 }
740};
741
745template <class TypeTag>
746class EnergyExtensiveQuantities<TypeTag, /*enableEnergy=*/true>
747{
752
753 enum { dimWorld = GridView::dimensionworld };
754 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
755 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
756
757protected:
762 void update_(const ElementContext& elemCtx, unsigned faceIdx, unsigned timeIdx)
763 {
764 const auto& gradCalc = elemCtx.gradientCalculator();
765 Opm::TemperatureCallback<TypeTag> temperatureCallback(elemCtx);
766
767 EvalDimVector temperatureGrad;
768 gradCalc.calculateGradient(temperatureGrad,
769 elemCtx,
770 faceIdx,
771 temperatureCallback);
772
773 // scalar product of temperature gradient and scvf normal
774 const auto& face = elemCtx.stencil(/*timeIdx=*/0).interiorFace(faceIdx);
775
776 temperatureGradNormal_ = 0;
777 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
778 temperatureGradNormal_ += (face.normal()[dimIdx]*temperatureGrad[dimIdx]);
779
780 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
781 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
782 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
783
784 // arithmetic mean
785 thermalConductivity_ =
786 0.5 * (intQuantsInside.thermalConductivity() + intQuantsOutside.thermalConductivity());
787 Opm::Valgrind::CheckDefined(thermalConductivity_);
788 }
789
790 template <class Context, class FluidState>
791 void updateBoundary_(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fs)
792 {
793 const auto& stencil = context.stencil(timeIdx);
794 const auto& face = stencil.boundaryFace(bfIdx);
795
796 const auto& elemCtx = context.elementContext();
797 unsigned insideScvIdx = face.interiorIndex();
798 const auto& insideScv = elemCtx.stencil(timeIdx).subControlVolume(insideScvIdx);
799
800 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
801 const auto& fsInside = intQuantsInside.fluidState();
802
803 // distance between the center of the SCV and center of the boundary face
804 DimVector distVec = face.integrationPos();
805 distVec -= insideScv.geometry().center();
806
807 Scalar tmp = 0;
808 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
809 tmp += distVec[dimIdx] * face.normal()[dimIdx];
810 Scalar dist = tmp;
811
812 // if the following assertation triggers, the center of the
813 // center of the interior SCV was not inside the element!
814 assert(dist > 0);
815
816 // calculate the temperature gradient using two-point gradient
817 // appoximation
818 temperatureGradNormal_ =
819 (fs.temperature(/*phaseIdx=*/0) - fsInside.temperature(/*phaseIdx=*/0)) / dist;
820
821 // take the value for thermal conductivity from the interior finite volume
822 thermalConductivity_ = intQuantsInside.thermalConductivity();
823 }
824
825public:
829 const Evaluation& temperatureGradNormal() const
830 { return temperatureGradNormal_; }
831
836 const Evaluation& thermalConductivity() const
837 { return thermalConductivity_; }
838
839private:
840 Evaluation temperatureGradNormal_;
841 Evaluation thermalConductivity_;
842};
843
844} // namespace Opm
845
846#endif
void updateBoundary_(const Context &, unsigned, unsigned, const FluidState &)
Definition: energymodule.hh:716
Scalar temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition: energymodule.hh:726
Scalar thermalConductivity() const
The total thermal conductivity at the face .
Definition: energymodule.hh:735
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:710
const Evaluation & temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition: energymodule.hh:829
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition: energymodule.hh:762
const Evaluation & thermalConductivity() const
The total thermal conductivity at the face .
Definition: energymodule.hh:836
void updateBoundary_(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fs)
Definition: energymodule.hh:791
Provides the quantities required to calculate energy fluxes.
Definition: energymodule.hh:694
Evaluation thermalConductivity() const
Returns the total thermal conductivity of the solid matrix in the sub-control volume.
Definition: energymodule.hh:559
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:569
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:583
Evaluation solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition: energymodule.hh:550
const Evaluation & solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition: energymodule.hh:672
const Evaluation & thermalConductivity() const
Returns the total conductivity capacity of the solid matrix in the sub-control volume.
Definition: energymodule.hh:679
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:640
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:616
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:530
static Scalar eqWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a equation of the residual.
Definition: energymodule.hh:104
static void setPriVarTemperatures(PrimaryVariables &, const FluidState &)
Given a fluid state, set the temperature in the primary variables.
Definition: energymodule.hh:113
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:197
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:153
static void addToEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:138
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:210
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:174
static void setEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:131
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:185
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:163
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:81
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:89
static Scalar thermalConductionRate(const ExtensiveQuantities &)
Add the rate of the conductive energy flux to a rate vector.
Definition: energymodule.hh:145
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:122
static Scalar primaryVarWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a primary variable for calculating relative errors.
Definition: energymodule.hh:96
static void registerParameters()
Register all run-time parameters for the energy module.
Definition: energymodule.hh:73
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:432
static void addToEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Add the rate of the enthalpy flux to a rate vector.
Definition: energymodule.hh:310
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:373
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:356
static void setEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Set the rate of energy flux of a rate vector.
Definition: energymodule.hh:304
static void setPriVarTemperatures(PrimaryVariables &priVars, const FluidState &fs)
Given a fluid state, set the temperature in the primary variables.
Definition: energymodule.hh:339
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:468
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:324
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:266
static void registerParameters()
Register all run-time parameters for the energy module.
Definition: energymodule.hh:246
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:392
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:403
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:254
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:277
static Scalar eqWeight(const Model &, unsigned, unsigned eqIdx)
Returns the relative weight of a equation.
Definition: energymodule.hh:289
static Evaluation thermalConductionRate(const ExtensiveQuantities &extQuants)
Returns the conductive energy flux for a given flux integration point.
Definition: energymodule.hh:316
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:48
Callback class for temperature.
Definition: quantitycallbacks.hh:48
Declare the properties used by the infrastructure code of the finite volume discretizations.
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:242
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:489