29#ifndef EWOMS_ENERGY_MODULE_HH
30#define EWOMS_ENERGY_MODULE_HH
32#include <dune/common/fvector.hh>
34#include <opm/material/common/Valgrind.hpp>
53template <
class TypeTag,
bool enableEnergy>
59template <
class TypeTag>
73 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
118 template <
class Flu
idState>
127 template <
class RateVector,
class Flu
idState>
158 template <
class LhsEval>
160 const IntensiveQuantities&,
168 template <
class LhsEval,
class Scv>
170 const IntensiveQuantities&,
179 template <
class LhsEval>
181 const IntensiveQuantities&)
190 template <
class Context>
202 template <
class Context>
215 template <
class Context>
226template <
class TypeTag>
241 enum { numPhases = FluidSystem::numPhases };
242 enum { energyEqIdx = Indices::energyEqIdx };
243 enum { temperatureIdx = Indices::temperatureIdx };
245 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
262 if (
pvIdx == temperatureIdx) {
263 return "temperature";
275 if (
eqIdx == energyEqIdx) {
287 if (
pvIdx != temperatureIdx) {
292 return std::max(
static_cast<Scalar
>(1.0) / 1000,
293 1.0 / model.solution(0)[globalDofIdx][temperatureIdx]);
303 if (
eqIdx != energyEqIdx) {
309 return static_cast<Scalar
>(1.0) / (4.184e3 * 30.0);
316 {
rateVec[energyEqIdx] = rate; }
322 {
rateVec[energyEqIdx] += rate; }
334 template <
class RateVector,
class Flu
idState>
336 const FluidState& fluidState,
338 const Evaluation& volume)
349 template <
class Flu
idState>
351 const FluidState& fs)
353 priVars[temperatureIdx] = Toolbox::value(fs.temperature(0));
356 assert(std::abs(Toolbox::value(fs.temperature(0))
357 - Toolbox::value(fs.temperature(
phaseIdx))) < 1
e-30);
366 template <
class LhsEval>
368 const IntensiveQuantities& intQuants,
371 const auto& fs = intQuants.fluidState();
383 template <
class Scv,
class LhsEval>
385 const IntensiveQuantities& intQuants,
389 const auto& fs = intQuants.fractureFluidState();
402 template <
class LhsEval>
404 const IntensiveQuantities& intQuants)
413 template <
class Context>
415 const Context& context,
423 if (!context.model().phaseIsConsidered(
phaseIdx)) {
429 const IntensiveQuantities&
up = context.intensiveQuantities(
upIdx,
timeIdx);
443 template <
class Context>
445 const Context& context,
458 if (!context.model().phaseIsConsidered(
phaseIdx)) {
464 const IntensiveQuantities&
up = context.intensiveQuantities(
upIdx,
timeIdx);
479 template <
class Context>
481 const Context& context,
498template <
unsigned PVOffset,
bool enableEnergy>
504template <
unsigned PVOffset>
508 enum { temperatureIdx = -1000 };
511 enum { energyEqIdx = -1000 };
520template <
unsigned PVOffset>
539template <
class TypeTag,
bool enableEnergy>
545template <
class TypeTag>
562 throw std::logic_error(
"solidInternalEnergy() does not make sense for isothermal models");
571 throw std::logic_error(
"thermalConductivity() does not make sense for isothermal models");
578 template <
class Flu
idState,
class Context>
580 const Context& context,
584 const Scalar T = context.problem().temperature(context,
spaceIdx,
timeIdx);
585 fluidState.setTemperature(Toolbox::createConstant(T));
592 template <
class Flu
idState>
595 const ElementContext&,
604template <
class TypeTag>
615 enum { numPhases = FluidSystem::numPhases };
616 enum { energyEqIdx = Indices::energyEqIdx };
617 enum { temperatureIdx = Indices::temperatureIdx };
625 template <
class Flu
idState,
class Context>
627 const Context& context,
632 if constexpr (std::is_same_v<Evaluation, Scalar>) {
633 fluidState.setTemperature(Toolbox::createConstant(priVars[temperatureIdx]));
638 fluidState.setTemperature(Toolbox::createVariable(priVars[temperatureIdx], temperatureIdx));
641 fluidState.setTemperature(Toolbox::createConstant(priVars[temperatureIdx]));
650 template <
class Flu
idState>
653 const ElementContext& elemCtx,
659 if (!elemCtx.model().phaseIsConsidered(
phaseIdx)) {
668 const auto& problem = elemCtx.problem();
669 const auto& solidEnergyParams = problem.solidEnergyLawParams(elemCtx, dofIdx,
timeIdx);
672 solidInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyParams, fs);
673 thermalConductivity_ = ThermalConductionLaw::thermalConductivity(
thermalCondParams, fs);
675 Opm::Valgrind::CheckDefined(solidInternalEnergy_);
676 Opm::Valgrind::CheckDefined(thermalConductivity_);
685 {
return solidInternalEnergy_; }
692 {
return thermalConductivity_; }
695 Evaluation solidInternalEnergy_;
696 Evaluation thermalConductivity_;
705template <
class TypeTag,
bool enableEnergy>
711template <
class TypeTag>
727 template <
class Context,
class Flu
idState>
728 void updateBoundary_(
const Context&,
740 throw std::logic_error(
"Calling temperatureGradNormal() does not make sense "
741 "for isothermal models");
749 throw std::logic_error(
"Calling thermalConductivity() does not make sense for "
750 "isothermal models");
757template <
class TypeTag>
765 enum { dimWorld = GridView::dimensionworld };
766 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
767 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
776 const auto&
gradCalc = elemCtx.gradientCalculator();
786 const auto&
face = elemCtx.stencil(0).interiorFace(faceIdx);
788 temperatureGradNormal_ = 0;
800 Opm::Valgrind::CheckDefined(thermalConductivity_);
803 template <
class Context,
class Flu
idState>
804 void updateBoundary_(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fs)
806 const auto& stencil = context.stencil(
timeIdx);
807 const auto&
face = stencil.boundaryFace(
bfIdx);
809 const auto& elemCtx = context.elementContext();
831 temperatureGradNormal_ =
843 {
return temperatureGradNormal_; }
850 {
return thermalConductivity_; }
853 Evaluation temperatureGradNormal_;
854 Evaluation thermalConductivity_;
Scalar temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition energymodule.hh:738
Scalar thermalConductivity() const
The total thermal conductivity at the face .
Definition energymodule.hh:747
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate energy fluxes.
Definition energymodule.hh:722
const Evaluation & temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition energymodule.hh:842
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition energymodule.hh:774
const Evaluation & thermalConductivity() const
The total thermal conductivity at the face .
Definition energymodule.hh:849
Provides the quantities required to calculate energy fluxes.
Definition energymodule.hh:706
Evaluation thermalConductivity() const
Returns the total thermal conductivity of the solid matrix in the sub-control volume.
Definition energymodule.hh:569
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:579
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:593
Evaluation solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition energymodule.hh:560
const Evaluation & solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition energymodule.hh:684
const Evaluation & thermalConductivity() const
Returns the total conductivity capacity of the solid matrix in the sub-control volume.
Definition energymodule.hh:691
void update_(FluidState &fs, typename FluidSystem::template ParameterCache< typename FluidState::Scalar > ¶mCache, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition energymodule.hh:651
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:626
Provides the volumetric quantities required for the energy equation.
Definition energymodule.hh:540
static Scalar eqWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a equation of the residual.
Definition energymodule.hh:110
static void setPriVarTemperatures(PrimaryVariables &, const FluidState &)
Given a fluid state, set the temperature in the primary variables.
Definition energymodule.hh:119
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:203
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:159
static void addToEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:144
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:216
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:180
static void setEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:137
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:191
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:169
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:87
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:95
static Scalar thermalConductionRate(const ExtensiveQuantities &)
Add the rate of the conductive energy flux to a rate vector.
Definition energymodule.hh:151
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:128
static Scalar primaryVarWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a primary variable for calculating relative errors.
Definition energymodule.hh:102
static void registerParameters()
Register all run-time parameters for the energy module.
Definition energymodule.hh:79
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:444
static void addToEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:321
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:384
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:367
static void setEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Set the rate of energy flux of a rate vector.
Definition energymodule.hh:315
static void setPriVarTemperatures(PrimaryVariables &priVars, const FluidState &fs)
Given a fluid state, set the temperature in the primary variables.
Definition energymodule.hh:350
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:480
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:335
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:273
static void registerParameters()
Register all run-time parameters for the energy module.
Definition energymodule.hh:252
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:403
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:414
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:260
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:285
static Scalar eqWeight(const Model &, unsigned, unsigned eqIdx)
Returns the relative weight of a equation.
Definition energymodule.hh:299
static Evaluation thermalConductionRate(const ExtensiveQuantities &extQuants)
Returns the conductive energy flux for a given flux integration point.
Definition energymodule.hh:327
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:54
Callback class for temperature.
Definition quantitycallbacks.hh:49
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.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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:233
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:499