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>
115 template <
class Flu
idState>
124 template <
class RateVector,
class Flu
idState>
155 template <
class LhsEval>
157 const IntensiveQuantities&,
165 template <
class LhsEval,
class Scv>
167 const IntensiveQuantities&,
176 template <
class LhsEval>
178 const IntensiveQuantities&)
187 template <
class Context>
199 template <
class Context>
212 template <
class Context>
223template <
class TypeTag>
237 enum { numPhases = FluidSystem::numPhases };
238 enum { energyEqIdx = Indices::energyEqIdx };
239 enum { temperatureIdx = Indices::temperatureIdx };
241 using Toolbox = Opm::MathToolbox<Evaluation>;
257 if (pvIdx == temperatureIdx) {
258 return "temperature";
268 static std::string
eqName(
unsigned eqIdx)
270 if (eqIdx == energyEqIdx) {
282 if (pvIdx != temperatureIdx) {
287 return std::max(
static_cast<Scalar
>(1.0) / 1000,
288 1.0 / model.solution(0)[globalDofIdx][temperatureIdx]);
298 if (eqIdx != energyEqIdx) {
304 return static_cast<Scalar
>(1.0) / (4.184e3 * 30.0);
311 { rateVec[energyEqIdx] = rate; }
317 { rateVec[energyEqIdx] += rate; }
323 {
return -extQuants.temperatureGradNormal() * extQuants.thermalConductivity(); }
329 template <
class RateVector,
class Flu
idState>
331 const FluidState& fluidState,
333 const Evaluation& volume)
335 rateVec[energyEqIdx] =
337 fluidState.density(phaseIdx) *
338 fluidState.enthalpy(phaseIdx);
344 template <
class Flu
idState>
346 const FluidState& fs)
348 priVars[temperatureIdx] = Toolbox::value(fs.temperature(0));
350 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
351 assert(std::abs(Toolbox::value(fs.temperature(0))
352 - Toolbox::value(fs.temperature(phaseIdx))) < 1e-30);
361 template <
class LhsEval>
363 const IntensiveQuantities& intQuants,
366 const auto& fs = intQuants.fluidState();
367 storage[energyEqIdx] +=
368 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
369 Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx)) *
370 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
371 Toolbox::template decay<LhsEval>(intQuants.porosity());
378 template <
class Scv,
class LhsEval>
380 const IntensiveQuantities& intQuants,
384 const auto& fs = intQuants.fractureFluidState();
385 storage[energyEqIdx] +=
386 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
387 Toolbox::template decay<LhsEval>(fs.internalEnergy(phaseIdx)) *
388 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
389 Toolbox::template decay<LhsEval>(intQuants.fracturePorosity()) *
390 Toolbox::template decay<LhsEval>(intQuants.fractureVolume()) / scv.volume();
397 template <
class LhsEval>
399 const IntensiveQuantities& intQuants)
400 { storage[energyEqIdx] += Opm::decay<LhsEval>(intQuants.solidInternalEnergy()); }
408 template <
class Context>
410 const Context& context,
414 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
417 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
418 if (!context.model().phaseIsConsidered(phaseIdx)) {
423 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
424 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
427 extQuants.volumeFlux(phaseIdx) *
428 up.fluidState().enthalpy(phaseIdx) *
429 up.fluidState().density(phaseIdx);
438 template <
class Context>
440 const Context& context,
444 const auto& scvf = context.stencil(timeIdx).interiorFace(spaceIdx);
445 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
449 flux[energyEqIdx] *= 1 - extQuants.fractureWidth() / (2 * scvf.area());
452 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
453 if (!context.model().phaseIsConsidered(phaseIdx)) {
458 const unsigned upIdx =
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
459 const IntensiveQuantities& up = context.intensiveQuantities(upIdx, timeIdx);
462 extQuants.fractureVolumeFlux(phaseIdx) *
463 up.fluidState().enthalpy(phaseIdx) *
464 up.fluidState().density(phaseIdx);
474 template <
class Context>
476 const Context& context,
480 const auto& extQuants = context.extensiveQuantities(spaceIdx, timeIdx);
483 flux[energyEqIdx] += -extQuants.temperatureGradNormal() * extQuants.thermalConductivity();
493template <
unsigned PVOffset,
bool enableEnergy>
499template <
unsigned PVOffset>
503 enum { temperatureIdx = -1000 };
506 enum { energyEqIdx = -1000 };
515template <
unsigned PVOffset>
519 enum { temperatureIdx = PVOffset };
522 enum { energyEqIdx = PVOffset };
534template <
class TypeTag,
bool enableEnergy>
540template <
class TypeTag>
548 using Toolbox = Opm::MathToolbox<Evaluation>;
557 throw std::logic_error(
"solidInternalEnergy() does not make sense for isothermal models");
566 throw std::logic_error(
"thermalConductivity() does not make sense for isothermal models");
573 template <
class Flu
idState,
class Context>
575 const Context& context,
579 const Scalar T = context.problem().temperature(context, spaceIdx, timeIdx);
580 fluidState.setTemperature(Toolbox::createConstant(T));
587 template <
class Flu
idState>
589 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>&,
590 const ElementContext&,
599template <
class TypeTag>
610 enum { numPhases = FluidSystem::numPhases };
611 enum { temperatureIdx = Indices::temperatureIdx };
613 using Toolbox = Opm::MathToolbox<Evaluation>;
619 template <
class Flu
idState,
class Context>
621 const Context& context,
625 const auto& priVars = context.primaryVars(spaceIdx, timeIdx);
626 if constexpr (std::is_same_v<Evaluation, Scalar>) {
627 fluidState.setTemperature(Toolbox::createConstant(priVars[temperatureIdx]));
632 fluidState.setTemperature(Toolbox::createVariable(priVars[temperatureIdx], temperatureIdx));
635 fluidState.setTemperature(Toolbox::createConstant(priVars[temperatureIdx]));
644 template <
class Flu
idState>
646 typename FluidSystem::template ParameterCache<typename FluidState::Scalar>& paramCache,
647 const ElementContext& elemCtx,
652 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
653 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
657 fs.setEnthalpy(phaseIdx,
658 FluidSystem::enthalpy(fs, paramCache, phaseIdx));
662 const auto& problem = elemCtx.problem();
663 const auto& solidEnergyParams = problem.solidEnergyLawParams(elemCtx, dofIdx, timeIdx);
664 const auto& thermalCondParams = problem.thermalConductionLawParams(elemCtx, dofIdx, timeIdx);
666 solidInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyParams, fs);
667 thermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalCondParams, fs);
669 Opm::Valgrind::CheckDefined(solidInternalEnergy_);
670 Opm::Valgrind::CheckDefined(thermalConductivity_);
679 {
return solidInternalEnergy_; }
686 {
return thermalConductivity_; }
689 Evaluation solidInternalEnergy_;
690 Evaluation thermalConductivity_;
699template <
class TypeTag,
bool enableEnergy>
705template <
class TypeTag>
721 template <
class Context,
class Flu
idState>
722 void updateBoundary_(
const Context&,
734 throw std::logic_error(
"Calling temperatureGradNormal() does not make sense "
735 "for isothermal models");
743 throw std::logic_error(
"Calling thermalConductivity() does not make sense for "
744 "isothermal models");
751template <
class TypeTag>
759 enum { dimWorld = GridView::dimensionworld };
760 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
761 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
768 void update_(
const ElementContext& elemCtx,
unsigned faceIdx,
unsigned timeIdx)
770 const auto& gradCalc = elemCtx.gradientCalculator();
773 EvalDimVector temperatureGrad;
774 gradCalc.calculateGradient(temperatureGrad,
777 temperatureCallback);
780 const auto& face = elemCtx.stencil(0).interiorFace(faceIdx);
782 temperatureGradNormal_ = 0;
783 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
784 temperatureGradNormal_ += (face.normal()[dimIdx]*temperatureGrad[dimIdx]);
787 const auto& extQuants = elemCtx.extensiveQuantities(faceIdx, timeIdx);
788 const auto& intQuantsInside = elemCtx.intensiveQuantities(extQuants.interiorIndex(), timeIdx);
789 const auto& intQuantsOutside = elemCtx.intensiveQuantities(extQuants.exteriorIndex(), timeIdx);
792 thermalConductivity_ = 0.5 * (intQuantsInside.thermalConductivity() +
793 intQuantsOutside.thermalConductivity());
794 Opm::Valgrind::CheckDefined(thermalConductivity_);
797 template <
class Context,
class Flu
idState>
798 void updateBoundary_(
const Context& context,
unsigned bfIdx,
unsigned timeIdx,
const FluidState& fs)
800 const auto& stencil = context.stencil(timeIdx);
801 const auto& face = stencil.boundaryFace(bfIdx);
803 const auto& elemCtx = context.elementContext();
804 const unsigned insideScvIdx = face.interiorIndex();
805 const auto& insideScv = elemCtx.stencil(timeIdx).subControlVolume(insideScvIdx);
807 const auto& intQuantsInside = elemCtx.intensiveQuantities(insideScvIdx, timeIdx);
808 const auto& fsInside = intQuantsInside.fluidState();
811 DimVector distVec = face.integrationPos();
812 distVec -= insideScv.geometry().center();
815 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
816 dist += distVec[dimIdx] * face.normal()[dimIdx];
825 temperatureGradNormal_ =
826 (fs.temperature(0) - fsInside.temperature(0)) / dist;
829 thermalConductivity_ = intQuantsInside.thermalConductivity();
837 {
return temperatureGradNormal_; }
844 {
return thermalConductivity_; }
847 Evaluation temperatureGradNormal_;
848 Evaluation thermalConductivity_;
Scalar temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition energymodule.hh:732
Scalar thermalConductivity() const
The total thermal conductivity at the face .
Definition energymodule.hh:741
void update_(const ElementContext &, unsigned, unsigned)
Update the quantities required to calculate energy fluxes.
Definition energymodule.hh:716
const Evaluation & temperatureGradNormal() const
The temperature gradient times the face normal [K m^2 / m].
Definition energymodule.hh:836
void update_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Update the quantities required to calculate energy fluxes.
Definition energymodule.hh:768
const Evaluation & thermalConductivity() const
The total thermal conductivity at the face .
Definition energymodule.hh:843
Provides the quantities required to calculate energy fluxes.
Definition energymodule.hh:700
Evaluation thermalConductivity() const
Returns the total thermal conductivity of the solid matrix in the sub-control volume.
Definition energymodule.hh:564
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:574
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:588
Evaluation solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition energymodule.hh:555
const Evaluation & solidInternalEnergy() const
Returns the volumetric internal energy of the solid matrix in the sub-control volume.
Definition energymodule.hh:678
const Evaluation & thermalConductivity() const
Returns the total conductivity capacity of the solid matrix in the sub-control volume.
Definition energymodule.hh:685
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:645
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:620
Provides the volumetric quantities required for the energy equation.
Definition energymodule.hh:535
static Scalar eqWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a equation of the residual.
Definition energymodule.hh:107
static void setPriVarTemperatures(PrimaryVariables &, const FluidState &)
Given a fluid state, set the temperature in the primary variables.
Definition energymodule.hh:116
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:200
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:156
static void addToEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:141
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:213
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:177
static void setEnthalpyRate(RateVector &, const Evaluation &)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:134
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:188
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:166
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:84
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:92
static Scalar thermalConductionRate(const ExtensiveQuantities &)
Add the rate of the conductive energy flux to a rate vector.
Definition energymodule.hh:148
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:125
static Scalar primaryVarWeight(const Model &, unsigned, unsigned)
Returns the relative weight of a primary variable for calculating relative errors.
Definition energymodule.hh:99
static void registerParameters()
Register all run-time parameters for the energy module.
Definition energymodule.hh:76
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:439
static void addToEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Add the rate of the enthalpy flux to a rate vector.
Definition energymodule.hh:316
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:379
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:362
static void setEnthalpyRate(RateVector &rateVec, const Evaluation &rate)
Set the rate of energy flux of a rate vector.
Definition energymodule.hh:310
static void setPriVarTemperatures(PrimaryVariables &priVars, const FluidState &fs)
Given a fluid state, set the temperature in the primary variables.
Definition energymodule.hh:345
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:475
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:330
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:247
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:398
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:409
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:255
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:280
static Scalar eqWeight(const Model &, unsigned, unsigned eqIdx)
Returns the relative weight of a equation.
Definition energymodule.hh:294
static Evaluation thermalConductionRate(const ExtensiveQuantities &extQuants)
Returns the conductive energy flux for a given flux integration point.
Definition energymodule.hh:322
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
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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:494