28#ifndef EWOMS_BLACK_OIL_ENERGY_MODULE_HH
29#define EWOMS_BLACK_OIL_ENERGY_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/common/Tabulated1DFunction.hpp>
34#include <opm/material/common/Valgrind.hpp>
71 static constexpr unsigned temperatureIdx = Indices::temperatureIdx;
72 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
76 static constexpr unsigned numPhases = FluidSystem::numPhases;
86 if constexpr (enableEnergy) {
97 if constexpr (enableEnergy) {
102 static bool primaryVarApplies(
unsigned pvIdx)
104 if constexpr (enableEnergy) {
105 return pvIdx == temperatureIdx;
116 return "temperature";
124 return static_cast<Scalar
>(1.0);
127 static bool eqApplies(
unsigned eqIdx)
129 if constexpr (enableEnergy) {
130 return eqIdx == contiEnergyEqIdx;
141 return "conti^energy";
152 template <
class LhsEval>
153 static void addStorage(Dune::FieldVector<LhsEval, numEq>&
storage,
154 const IntensiveQuantities& intQuants)
156 if constexpr (enableEnergy) {
160 const auto& fs = intQuants.fluidState();
162 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
170 storage[contiEnergyEqIdx] += poro*
S*
u*rho;
174 const Scalar rockFraction = intQuants.rockFraction();
186 if constexpr (enableEnergy) {
187 flux[contiEnergyEqIdx] = 0.0;
190 const unsigned focusIdx = elemCtx.focusDofIndex();
192 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
211 static void addHeatFlux(RateVector&
flux,
214 if constexpr (enableEnergy) {
221 template <
class UpEval,
class Eval,
class Flu
idState>
222 static void addPhaseEnthalpyFluxes_(RateVector&
flux,
224 const Eval& volumeFlux,
225 const FluidState&
upFs)
227 flux[contiEnergyEqIdx] +=
233 template <
class UpstreamEval>
234 static void addPhaseEnthalpyFlux_(RateVector&
flux,
236 const ElementContext& elemCtx,
243 const auto& fs =
up.fluidState();
251 static void addToEnthalpyRate(RateVector&
flux,
252 const Evaluation&
hRate)
254 if constexpr (enableEnergy) {
262 template <
class Flu
idState>
264 const FluidState& fluidState)
266 if constexpr (enableEnergy) {
267 priVars[temperatureIdx] =
getValue(fluidState.temperature(0));
275 const PrimaryVariables&
oldPv,
276 const EqVector&
delta)
278 if constexpr (enableEnergy) {
293 return static_cast<Scalar
>(0.0);
305 template <
class DofEntity>
308 if constexpr (enableEnergy) {
309 const unsigned dofIdx = model.dofMapper().index(
dof);
310 const PrimaryVariables& priVars = model.solution(0)[dofIdx];
315 template <
class DofEntity>
318 if constexpr (enableEnergy) {
319 const unsigned dofIdx = model.dofMapper().index(
dof);
320 PrimaryVariables&
priVars0 = model.solution(0)[dofIdx];
321 PrimaryVariables&
priVars1 = model.solution(1)[dofIdx];
331template <
class TypeTag,
bool enableEnergyV>
341template <
class TypeTag>
357 static constexpr int temperatureIdx = Indices::temperatureIdx;
368 auto& fs = asImp_().fluidState_;
369 const auto& priVars = elemCtx.primaryVars(dofIdx,
timeIdx);
372 fs.setTemperature(priVars.makeEvaluation(temperatureIdx,
timeIdx, elemCtx.linearizationType()));
380 const PrimaryVariables& priVars,
385 auto& fs = asImp_().fluidState_;
386 fs.setTemperature(priVars.makeEvaluation(temperatureIdx,
timeIdx,
lintype));
397 updateEnergyQuantities_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx,
timeIdx),
timeIdx);
400 void updateEnergyQuantities_(
const Problem& problem,
404 auto& fs = asImp_().fluidState_;
409 if (!FluidSystem::phaseIsActive(
phaseIdx)) {
418 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fs);
420 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(
globalSpaceIdx,
timeIdx);
421 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fs);
431 const Evaluation& rockInternalEnergy()
const
432 {
return rockInternalEnergy_; }
434 const Evaluation& totalThermalConductivity()
const
435 {
return totalThermalConductivity_; }
437 Scalar rockFraction()
const
438 {
return rockFraction_; }
441 Implementation& asImp_()
442 {
return *
static_cast<Implementation*
>(
this); }
444 Evaluation rockInternalEnergy_;
445 Evaluation totalThermalConductivity_;
446 Scalar rockFraction_;
449template <
class TypeTag>
463 void updateTemperature_([[
maybe_unused]]
const ElementContext& elemCtx,
467 if constexpr (enableTemperature) {
470 auto& fs = asImp_().fluidState_;
471 const Scalar T = elemCtx.problem().temperature(elemCtx, dofIdx,
timeIdx);
472 fs.setTemperature(T);
476 template<
class Problem>
477 void updateTemperature_([[
maybe_unused]]
const Problem& problem,
484 if constexpr (enableTemperature) {
485 auto& fs = asImp_().fluidState_;
488 const Scalar T = problem.temperature(globalDofIdx,
timeIdx);
489 fs.setTemperature(T);
493 void updateEnergyQuantities_(
const ElementContext&,
499 const Evaluation& rockInternalEnergy()
const
501 throw std::logic_error(
"Requested the rock internal energy, which is "
502 "unavailable because energy is not conserved");
505 const Evaluation& totalThermalConductivity()
const
507 throw std::logic_error(
"Requested the total thermal conductivity, which is "
508 "unavailable because energy is not conserved");
512 Implementation& asImp_()
513 {
return *
static_cast<Implementation*
>(
this); }
516template <
class TypeTag,
bool enableEnergyV>
526template <
class TypeTag>
537 template<
class Flu
idState>
538 static void updateEnergy(Evaluation& energyFlux,
539 const unsigned& focusDofIndex,
540 const unsigned&
inIdx,
541 const unsigned&
exIdx,
542 const IntensiveQuantities&
inIq,
543 const IntensiveQuantities&
exIq,
544 const FluidState&
inFs,
545 const FluidState&
exFs,
546 const Scalar& inAlpha,
547 const Scalar& outAlpha,
548 const Scalar& faceArea)
551 if (focusDofIndex ==
inIdx) {
555 else if (focusDofIndex ==
exIdx) {
565 if (focusDofIndex ==
inIdx) {
573 if (focusDofIndex ==
exIdx) {
588 H = 1.0 / (1.0 /
inH + 1.0 /
exH);
594 energyFlux =
deltaT * (-
H / faceArea);
597 void updateEnergy(
const ElementContext& elemCtx,
601 const auto& stencil = elemCtx.stencil(
timeIdx);
604 const Scalar faceArea =
scvf.area();
605 const unsigned inIdx =
scvf.interiorIndex();
606 const unsigned exIdx =
scvf.exteriorIndex();
609 const auto&
inFs =
inIq.fluidState();
610 const auto&
exFs =
exIq.fluidState();
611 const Scalar inAlpha = elemCtx.problem().thermalHalfTransmissibilityIn(elemCtx,
scvfIdx,
timeIdx);
612 const Scalar outAlpha = elemCtx.problem().thermalHalfTransmissibilityOut(elemCtx,
scvfIdx,
timeIdx);
613 updateEnergy(energyFlux_,
614 elemCtx.focusDofIndex(),
626 template <
class Context,
class BoundaryFlu
idState>
627 void updateEnergyBoundary(
const Context&
ctx,
635 const unsigned inIdx =
scvf.interiorIndex();
638 const Scalar alpha =
ctx.problem().thermalHalfTransmissibilityBoundary(
ctx,
scvfIdx);
642 template <
class BoundaryFlu
idState>
643 static void updateEnergyBoundary(Evaluation& energyFlux,
644 const IntensiveQuantities&
inIq,
645 unsigned focusDofIndex,
650 const auto&
inFs =
inIq.fluidState();
652 if (focusDofIndex ==
inIdx) {
662 if (focusDofIndex ==
inIdx) {
681 const Evaluation& energyFlux()
const
682 {
return energyFlux_; }
685 Implementation& asImp_()
686 {
return *
static_cast<Implementation*
>(
this); }
688 Evaluation energyFlux_;
691template <
class TypeTag>
700 template<
class Flu
idState>
701 static void updateEnergy(Evaluation& ,
705 const IntensiveQuantities& ,
706 const IntensiveQuantities& ,
714 void updateEnergy(
const ElementContext&,
719 template <
class Context,
class BoundaryFlu
idState>
720 void updateEnergyBoundary(
const Context&,
726 template <
class BoundaryFlu
idState>
727 static void updateEnergyBoundary(Evaluation& ,
728 const IntensiveQuantities& ,
736 const Evaluation& energyFlux()
const
737 {
throw std::logic_error(
"Requested the energy flux, but energy is not conserved"); }
Declares the properties required by the black oil model.
Provides the energy specific extensive quantities to the generic black-oil module's extensive quantit...
Definition blackoilenergymodules.hh:517
void updateEnergyQuantities_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities needed to handle energy conservation.
Definition blackoilenergymodules.hh:393
void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the temperature of the intensive quantity's fluid state.
Definition blackoilenergymodules.hh:364
void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, const unsigned timeIdx, const LinearizationType &lintype)
Update the temperature of the intensive quantity's fluid state.
Definition blackoilenergymodules.hh:379
Provides the volumetric quantities required for the equations needed by the energys extension of the ...
Definition blackoilenergymodules.hh:332
Contains the high level supplements required to extend the black oil model by energy.
Definition blackoilenergymodules.hh:58
static void registerParameters()
Register all run-time parameters for the black-oil energy module.
Definition blackoilenergymodules.hh:84
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition blackoilenergymodules.hh:299
static void updatePrimaryVars(PrimaryVariables &newPv, const PrimaryVariables &oldPv, const EqVector &delta)
Do a Newton-Raphson update the primary variables of the energys.
Definition blackoilenergymodules.hh:274
static void registerOutputModules(Model &model, Simulator &simulator)
Register all energy specific VTK and ECL output modules.
Definition blackoilenergymodules.hh:94
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition blackoilenergymodules.hh:287
static void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the energy specific primary variables to a PrimaryVariables object.
Definition blackoilenergymodules.hh:263
VTK output module for the black oil model's energy related quantities.
Definition vtkblackoilenergymodule.hpp:54
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkblackoilenergymodule.hpp:87
The common code for the linearizers of non-linear systems of equations.
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...
Definition linearizationtype.hh:34
VTK output module for the black oil model's energy related quantities.