blackoilenergymodules.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*/
28#ifndef EWOMS_BLACK_OIL_ENERGY_MODULE_HH
29#define EWOMS_BLACK_OIL_ENERGY_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/common/utility/gpuDecorators.hpp>
34
35#include <opm/material/common/Tabulated1DFunction.hpp>
36#include <opm/material/common/Valgrind.hpp>
37#include <opm/material/fluidstates/BlackOilFluidState.hpp>
38
43#include <opm/material/thermal/EnergyModuleType.hpp>
44
45#include <cassert>
46#include <cmath>
47#include <istream>
48#include <memory>
49#include <ostream>
50#include <stdexcept>
51#include <string>
52
53
54namespace Opm {
55
61template <class TypeTag, EnergyModules activeModule = getPropValue<TypeTag, Properties::EnergyModuleType>()>
63{
75
76 static constexpr unsigned temperatureIdx = Indices::temperatureIdx;
77 static constexpr unsigned contiEnergyEqIdx = Indices::contiEnergyEqIdx;
78
79 static constexpr unsigned enableFullyImplicitThermal = (activeModule == EnergyModules::FullyImplicitThermal);
80 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
81 static constexpr unsigned numPhases = FluidSystem::numPhases;
82
83public:
85
89 static void registerParameters()
90 {
91 if constexpr (enableFullyImplicitThermal) {
93 }
94 }
95
99 static void registerOutputModules(Model& model,
100 Simulator& simulator)
101 {
102 if constexpr (enableFullyImplicitThermal) {
103 model.addOutputModule(std::make_unique<VtkBlackOilEnergyModule<TypeTag>>(simulator));
104 }
105 }
106
107 static bool primaryVarApplies(unsigned pvIdx)
108 {
109 if constexpr (enableFullyImplicitThermal) {
110 return pvIdx == temperatureIdx;
111 }
112 else {
113 return false;
114 }
115 }
116
117 static std::string primaryVarName([[maybe_unused]] unsigned pvIdx)
118 {
119 assert(primaryVarApplies(pvIdx));
120
121 return "temperature";
122 }
123
124 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
125 {
126 assert(primaryVarApplies(pvIdx));
127
128 // TODO: it may be beneficial to chose this differently.
129 return static_cast<Scalar>(1.0);
130 }
131
132 static bool eqApplies(unsigned eqIdx)
133 {
134 if constexpr (enableFullyImplicitThermal) {
135 return eqIdx == contiEnergyEqIdx;
136 }
137 else {
138 return false;
139 }
140 }
141
142 static std::string eqName([[maybe_unused]] unsigned eqIdx)
143 {
144 assert(eqApplies(eqIdx));
145
146 return "conti^energy";
147 }
148
149 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
150 {
151 assert(eqApplies(eqIdx));
152
153 return 1.0;
154 }
155
156 // must be called after water storage is computed
157 template <class StorageType>
158 OPM_HOST_DEVICE static void addStorage(StorageType& storage,
159 const IntensiveQuantities& intQuants)
160 {
161 using LhsEval = typename StorageType::value_type;
162
163 if constexpr (enableFullyImplicitThermal) {
164 const FluidSystem& fsys = intQuants.getFluidSystem();
165
166 const auto& poro = decay<LhsEval>(intQuants.porosity());
167
168 // accumulate the internal energy of the fluids
169 const auto& fs = intQuants.fluidState();
170 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
171 if (!fsys.phaseIsActive(phaseIdx)) {
172 continue;
173 }
174
175 const auto& u = decay<LhsEval>(fs.internalEnergy(phaseIdx));
176 const auto& S = decay<LhsEval>(fs.saturation(phaseIdx));
177 const auto& rho = decay<LhsEval>(fs.density(phaseIdx));
178
179 storage[contiEnergyEqIdx] += poro*S*u*rho;
180 }
181
182 // add the internal energy of the rock
183 const Scalar rockFraction = intQuants.rockFraction();
184 const auto& uRock = decay<LhsEval>(intQuants.rockInternalEnergy());
185 storage[contiEnergyEqIdx] += rockFraction * uRock;
186 storage[contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
187 }
188 }
189
190 OPM_HOST_DEVICE static void computeFlux([[maybe_unused]] RateVector& flux,
191 [[maybe_unused]] const ElementContext& elemCtx,
192 [[maybe_unused]] unsigned scvfIdx,
193 [[maybe_unused]] unsigned timeIdx)
194 {
195 if constexpr (enableFullyImplicitThermal) {
196 flux[contiEnergyEqIdx] = 0.0;
197
198 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
199 const unsigned focusIdx = elemCtx.focusDofIndex();
200 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
201 if (!FluidSystem::phaseIsActive(phaseIdx)) {
202 continue;
203 }
204
205 const unsigned upIdx = extQuants.upstreamIndex(phaseIdx);
206 if (upIdx == focusIdx) {
207 addPhaseEnthalpyFlux_<Evaluation>(flux, phaseIdx, elemCtx, scvfIdx, timeIdx);
208 }
209 else {
210 addPhaseEnthalpyFlux_<Scalar>(flux, phaseIdx, elemCtx, scvfIdx, timeIdx);
211 }
212 }
213
214 // diffusive energy flux
215 flux[contiEnergyEqIdx] += extQuants.energyFlux();
216 flux[contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
217 }
218 }
219
220 template<class RateVectorT>
221 OPM_HOST_DEVICE static void addHeatFlux(RateVectorT& flux,
222 const Evaluation& heatFlux)
223 {
224 if constexpr (enableFullyImplicitThermal) {
225 // diffusive energy flux
226 flux[contiEnergyEqIdx] += heatFlux;
227 flux[contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
228 }
229 }
230
231 template <class UpEval, class RateVectorT, class Eval, class FluidState>
232 OPM_HOST_DEVICE static void addPhaseEnthalpyFluxes_(RateVectorT& flux,
233 unsigned phaseIdx,
234 const Eval& volumeFlux,
235 const FluidState& upFs)
236 {
237 flux[contiEnergyEqIdx] +=
238 decay<UpEval>(upFs.enthalpy(phaseIdx)) *
239 decay<UpEval>(upFs.density(phaseIdx)) *
240 volumeFlux;
241 }
242
243 template <class UpstreamEval>
244 static void addPhaseEnthalpyFlux_(RateVector& flux,
245 unsigned phaseIdx,
246 const ElementContext& elemCtx,
247 unsigned scvfIdx,
248 unsigned timeIdx)
249 {
250 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
251 const unsigned upIdx = extQuants.upstreamIndex(phaseIdx);
252 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
253 const auto& fs = up.fluidState();
254 const auto& volFlux = extQuants.volumeFlux(phaseIdx);
255 addPhaseEnthalpyFluxes_<UpstreamEval>(flux,
256 phaseIdx,
257 volFlux,
258 fs);
259 }
260
261 static void addToEnthalpyRate(RateVector& flux,
262 const Evaluation& hRate)
263 {
264 if constexpr (enableFullyImplicitThermal) {
265 flux[contiEnergyEqIdx] += hRate;
266 }
267 }
268
272 template <class FluidState>
273 static void assignPrimaryVars(PrimaryVariables& priVars,
274 const FluidState& fluidState)
275 {
276 if constexpr (enableFullyImplicitThermal) {
277 priVars[temperatureIdx] = getValue(fluidState.temperature(/*phaseIdx=*/0));
278 }
279 }
280
284 static void updatePrimaryVars(PrimaryVariables& newPv,
285 const PrimaryVariables& oldPv,
286 const EqVector& delta)
287 {
288 if constexpr (enableFullyImplicitThermal) {
289 // do a plain unchopped Newton update
290 newPv[temperatureIdx] = oldPv[temperatureIdx] - delta[temperatureIdx];
291 }
292 }
293
297 static Scalar computeUpdateError(const PrimaryVariables&,
298 const EqVector&)
299 {
300 // do not consider consider the cange of energy primary variables for
301 // convergence
302 // TODO: maybe this should be changed
303 return static_cast<Scalar>(0.0);
304 }
305
309 static Scalar computeResidualError(const EqVector& resid)
310 {
311 // do not weight the residual of energy when it comes to convergence
312 return std::abs(scalarValue(resid[contiEnergyEqIdx]));
313 }
314
315 template <class DofEntity>
316 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
317 {
318 if constexpr (enableFullyImplicitThermal) {
319 const unsigned dofIdx = model.dofMapper().index(dof);
320 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
321 outstream << priVars[temperatureIdx];
322 }
323 }
324
325 template <class DofEntity>
326 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
327 {
328 if constexpr (enableFullyImplicitThermal) {
329 const unsigned dofIdx = model.dofMapper().index(dof);
330 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
331 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
332
333 instream >> priVars0[temperatureIdx];
334
335 // set the primary variables for the beginning of the current time step.
336 priVars1 = priVars0[temperatureIdx];
337 }
338 }
339};
340
341template <class TypeTag, EnergyModules activeModule>
343
351template <class TypeTag>
352class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::FullyImplicitThermal>
353{
355
365
366 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
367 static constexpr int temperatureIdx = Indices::temperatureIdx;
368
369public:
373 BlackOilEnergyIntensiveQuantities(Evaluation rockInternalEnergy,
374 Evaluation totalThermalConductivity,
375 Scalar rockFraction)
376 : rockInternalEnergy_(rockInternalEnergy)
377 , totalThermalConductivity_(totalThermalConductivity)
378 , rockFraction_(rockFraction)
379 {
380 }
381
383
388 void updateTemperature_(const ElementContext& elemCtx,
389 unsigned dofIdx,
390 unsigned timeIdx)
391 {
392 auto& fs = asImp_().fluidState_;
393 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
394
395 // set temperature
396 fs.setTemperature(priVars.makeEvaluation(temperatureIdx, timeIdx, elemCtx.linearizationType()));
397 }
398
403 void updateTemperature_([[maybe_unused]] const Problem& problem,
404 const PrimaryVariables& priVars,
405 [[maybe_unused]] unsigned globalDofIdx,
406 const unsigned timeIdx,
407 const LinearizationType& lintype)
408 {
409 auto& fs = asImp_().fluidState_;
410 fs.setTemperature(priVars.makeEvaluation(temperatureIdx, timeIdx, lintype));
411 }
412
417 void updateEnergyQuantities_(const ElementContext& elemCtx,
418 unsigned dofIdx,
419 unsigned timeIdx)
420 {
421 updateEnergyQuantities_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
422 }
423
424 void updateEnergyQuantities_(const Problem& problem,
425 const unsigned globalSpaceIdx,
426 const unsigned timeIdx)
427 {
428 auto& fs = asImp_().fluidState_;
429
430 // compute the specific enthalpy of the fluids, the specific enthalpy of the rock
431 // and the thermal conductivity coefficients
432 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
433 if (!FluidSystem::phaseIsActive(phaseIdx)) {
434 continue;
435 }
436
437 const auto& h = FluidSystem::enthalpy(fs, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
438 fs.setEnthalpy(phaseIdx, h);
439 }
440
441 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
442 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fs);
443
444 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
445 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fs);
446
447 // Retrieve the rock fraction from the problem
448 // Usually 1 - porosity, but if pvmult is used to modify porosity
449 // we will apply the same multiplier to the rock fraction
450 // i.e. pvmult*(1 - porosity) and thus interpret multpv as a volume
451 // multiplier. This is to avoid negative rock volume for pvmult*porosity > 1
452 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
453 }
454
455 OPM_HOST_DEVICE const Evaluation& rockInternalEnergy() const
456 { return rockInternalEnergy_; }
457
458 OPM_HOST_DEVICE const Evaluation& totalThermalConductivity() const
459 { return totalThermalConductivity_; }
460
461 OPM_HOST_DEVICE Scalar rockFraction() const
462 { return rockFraction_; }
463
464protected:
465 Implementation& asImp_()
466 { return *static_cast<Implementation*>(this); }
467
471};
472
473template <class TypeTag>
474class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::ConstantTemperature>
475{
482
484
485public:
486
487 void updateTemperature_(const ElementContext& elemCtx,
488 unsigned dofIdx,
489 unsigned timeIdx)
490 {
491 updateTemperature_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
492 }
493
494 template<class Problem>
495 void updateTemperature_(const Problem& problem,
496 [[maybe_unused]] const PrimaryVariables& priVars,
497 unsigned globalDofIdx,
498 unsigned timeIdx,
499 [[maybe_unused]] const LinearizationType& lintype
500 )
501 {
502 updateTemperature_(problem, globalDofIdx, timeIdx);
503 }
504
505 void updateTemperature_(const Problem& problem, unsigned globalDofIdx, unsigned timeIdx)
506 {
507 auto& fs = asImp_().fluidState_;
508 const Scalar T = problem.temperature(globalDofIdx, timeIdx);
509 fs.setTemperature(T);
510 }
511
512 void updateEnergyQuantities_(const ElementContext&,
513 unsigned,
514 unsigned,
515 const typename FluidSystem::template ParameterCache<Evaluation>&)
516 {}
517
518 const Evaluation& rockInternalEnergy() const
519 {
520 throw std::logic_error("Requested the rock internal energy, which is "
521 "unavailable because energy is not conserved");
522 }
523
524 const Evaluation& totalThermalConductivity() const
525 {
526 throw std::logic_error("Requested the total thermal conductivity, which is "
527 "unavailable because energy is not conserved");
528 }
529
530protected:
531 Implementation& asImp_()
532 { return *static_cast<Implementation*>(this); }
533};
534
535template <class TypeTag>
536class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::SequentialImplicitThermal>
537{
548
549public:
550
551 void updateTemperature_(const Problem& problem, unsigned globalDofIdx, unsigned timeIdx)
552 {
553 // update the temperature for output (without derivatives)
554 auto& fs = asImp_().fluidState_;
555 fs.setTemperature(problem.temperature(globalDofIdx, timeIdx));
556 }
557
558 void updateTemperature_(const ElementContext& elemCtx,
559 unsigned dofIdx,
560 unsigned timeIdx)
561 {
562 updateTemperature_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
563 }
564
565 template<class Problem>
566 void updateTemperature_(const Problem& problem,
567 [[maybe_unused]] const PrimaryVariables& priVars,
568 unsigned globalDofIdx,
569 unsigned timeIdx,
570 [[maybe_unused]] const LinearizationType& lintype)
571 {
572 updateTemperature_(problem, globalDofIdx, timeIdx);
573 }
574
579 void updateEnergyQuantities_([[maybe_unused]] const ElementContext& elemCtx,
580 [[maybe_unused]] unsigned dofIdx,
581 [[maybe_unused]] unsigned timeIdx)
582 {
583 }
584
585 void updateEnergyQuantities_([[maybe_unused]] const Problem& problem,
586 [[maybe_unused]] const unsigned globalSpaceIdx,
587 [[maybe_unused]] const unsigned timeIdx)
588 {
589 }
590
591 const Evaluation& rockInternalEnergy() const
592 {
593 throw std::logic_error("Requested the rock internal energy, which is "
594 "unavailable because energy is not conserved");
595 }
596
597 const Evaluation& totalThermalConductivity() const
598 {
599 throw std::logic_error("Requested the total thermal conductivity, which is "
600 "unavailable because energy is not conserved");
601 }
602
603protected:
604 Implementation& asImp_()
605 { return *static_cast<Implementation*>(this); }
606
607};
608
609template <class TypeTag>
610class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::NoTemperature>
611{
618
619public:
620 void updateTemperature_([[maybe_unused]] const ElementContext& elemCtx,
621 [[maybe_unused]] unsigned dofIdx,
622 [[maybe_unused]] unsigned timeIdx)
623 {
624 }
625
626 template<class Problem>
627 void updateTemperature_([[maybe_unused]] const Problem& problem,
628 [[maybe_unused]] const PrimaryVariables& priVars,
629 [[maybe_unused]] unsigned globalDofIdx,
630 [[maybe_unused]] unsigned timeIdx,
631 [[maybe_unused]] const LinearizationType& lintype)
632 {
633 }
634
635 void updateEnergyQuantities_(const ElementContext&,
636 unsigned,
637 unsigned,
638 const typename FluidSystem::template ParameterCache<Evaluation>&)
639 {}
640
641 const Evaluation& rockInternalEnergy() const
642 {
643 throw std::logic_error("Requested the rock internal energy, which is "
644 "unavailable because energy is not conserved");
645 }
646
647 const Evaluation& totalThermalConductivity() const
648 {
649 throw std::logic_error("Requested the total thermal conductivity, which is "
650 "unavailable because energy is not conserved");
651 }
652
653protected:
654 Implementation& asImp_()
655 { return *static_cast<Implementation*>(this); }
656};
657
658template <class TypeTag, EnergyModules activeModule>
660
668template <class TypeTag>
669class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::FullyImplicitThermal>
670{
672
676
677public:
678 template<class Evaluation, class FluidState, class IntensiveQuantities>
679 OPM_HOST_DEVICE static void updateEnergy(Evaluation& energyFlux,
680 const unsigned& focusDofIndex,
681 const unsigned& inIdx,
682 const unsigned& exIdx,
683 const IntensiveQuantities& inIq,
684 const IntensiveQuantities& exIq,
685 const FluidState& inFs,
686 const FluidState& exFs,
687 const Scalar& inAlpha,
688 const Scalar& outAlpha,
689 const Scalar& faceArea)
690 {
691 Evaluation deltaT;
692 if (focusDofIndex == inIdx) {
693 deltaT = decay<Scalar>(exFs.temperature(/*phaseIdx=*/0)) -
694 inFs.temperature(/*phaseIdx=*/0);
695 }
696 else if (focusDofIndex == exIdx) {
697 deltaT = exFs.temperature(/*phaseIdx=*/0) -
698 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
699 }
700 else {
701 deltaT = decay<Scalar>(exFs.temperature(/*phaseIdx=*/0)) -
702 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
703 }
704
705 Evaluation inLambda;
706 if (focusDofIndex == inIdx) {
707 inLambda = inIq.totalThermalConductivity();
708 }
709 else {
710 inLambda = decay<Scalar>(inIq.totalThermalConductivity());
711 }
712
713 Evaluation exLambda;
714 if (focusDofIndex == exIdx) {
715 exLambda = exIq.totalThermalConductivity();
716 }
717 else {
718 exLambda = decay<Scalar>(exIq.totalThermalConductivity());
719 }
720
721 Evaluation H;
722 const Evaluation& inH = inLambda*inAlpha;
723 const Evaluation& exH = exLambda*outAlpha;
724 if (inH > 0 && exH > 0) {
725 // compute the "thermal transmissibility". In contrast to the normal
726 // transmissibility this cannot be done as a preprocessing step because the
727 // average thermal conductivity is analogous to the permeability but
728 // depends on the solution.
729 H = 1.0 / (1.0 / inH + 1.0 / exH);
730 }
731 else {
732 H = 0.0;
733 }
734
735 energyFlux = deltaT * (-H / faceArea);
736 }
737
738 void updateEnergy(const ElementContext& elemCtx,
739 unsigned scvfIdx,
740 unsigned timeIdx)
741 {
742 const auto& stencil = elemCtx.stencil(timeIdx);
743 const auto& scvf = stencil.interiorFace(scvfIdx);
744
745 const Scalar faceArea = scvf.area();
746 const unsigned inIdx = scvf.interiorIndex();
747 const unsigned exIdx = scvf.exteriorIndex();
748 const auto& inIq = elemCtx.intensiveQuantities(inIdx, timeIdx);
749 const auto& exIq = elemCtx.intensiveQuantities(exIdx, timeIdx);
750 const auto& inFs = inIq.fluidState();
751 const auto& exFs = exIq.fluidState();
752 const Scalar inAlpha = elemCtx.problem().thermalHalfTransmissibilityIn(elemCtx, scvfIdx, timeIdx);
753 const Scalar outAlpha = elemCtx.problem().thermalHalfTransmissibilityOut(elemCtx, scvfIdx, timeIdx);
754 updateEnergy(energyFlux_,
755 elemCtx.focusDofIndex(),
756 inIdx,
757 exIdx,
758 inIq,
759 exIq,
760 inFs,
761 exFs,
762 inAlpha,
763 outAlpha,
764 faceArea);
765 }
766
767 template <class Context, class BoundaryFluidState>
768 void updateEnergyBoundary(const Context& ctx,
769 unsigned scvfIdx,
770 unsigned timeIdx,
771 const BoundaryFluidState& boundaryFs)
772 {
773 const auto& stencil = ctx.stencil(timeIdx);
774 const auto& scvf = stencil.boundaryFace(scvfIdx);
775
776 const unsigned inIdx = scvf.interiorIndex();
777 const auto& inIq = ctx.intensiveQuantities(inIdx, timeIdx);
778 const auto& focusDofIdx = ctx.focusDofIndex();
779 const Scalar alpha = ctx.problem().thermalHalfTransmissibilityBoundary(ctx, scvfIdx);
780 updateEnergyBoundary(energyFlux_, inIq, focusDofIdx, inIdx, alpha, boundaryFs);
781 }
782
783 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
784 OPM_HOST_DEVICE static void updateEnergyBoundary(Evaluation& energyFlux,
785 const IntensiveQuantities& inIq,
786 unsigned focusDofIndex,
787 unsigned inIdx,
788 Scalar alpha,
789 const BoundaryFluidState& boundaryFs)
790 {
791 const auto& inFs = inIq.fluidState();
792 Evaluation deltaT;
793 if (focusDofIndex == inIdx) {
794 deltaT = boundaryFs.temperature(/*phaseIdx=*/0) -
795 inFs.temperature(/*phaseIdx=*/0);
796 }
797 else {
798 deltaT = decay<Scalar>(boundaryFs.temperature(/*phaseIdx=*/0)) -
799 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
800 }
801
802 Evaluation lambda;
803 if (focusDofIndex == inIdx) {
804 lambda = inIq.totalThermalConductivity();
805 }
806 else {
807 lambda = decay<Scalar>(inIq.totalThermalConductivity());
808 }
809
810 if (lambda > 0.0) {
811 // compute the "thermal transmissibility". In contrast to the normal
812 // transmissibility this cannot be done as a preprocessing step because the
813 // average thermal conductivity is analogous to the permeability but depends
814 // on the solution.
815 energyFlux = deltaT * lambda * -alpha;
816 }
817 else {
818 energyFlux = 0.0;
819 }
820 }
821
822 const Evaluation& energyFlux() const
823 { return energyFlux_; }
824
825private:
826 Implementation& asImp_()
827 { return *static_cast<Implementation*>(this); }
828
829 Evaluation energyFlux_;
830};
831
832template <class TypeTag>
833class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::ConstantTemperature>
834{
838
839public:
840 template<class Evaluation, class FluidState, class IntensiveQuantities>
841 static void updateEnergy(Evaluation& /*energyFlux*/,
842 const unsigned& /*focusDofIndex*/,
843 const unsigned& /*inIdx*/,
844 const unsigned& /*exIdx*/,
845 const IntensiveQuantities& /*inIq*/,
846 const IntensiveQuantities& /*exIq*/,
847 const FluidState& /*inFs*/,
848 const FluidState& /*exFs*/,
849 const Scalar& /*inAlpha*/,
850 const Scalar& /*outAlpha*/,
851 const Scalar& /*faceArea*/)
852 {}
853
854 void updateEnergy(const ElementContext&,
855 unsigned,
856 unsigned)
857 {}
858
859 template <class Context, class BoundaryFluidState>
860 void updateEnergyBoundary(const Context&,
861 unsigned,
862 unsigned,
863 const BoundaryFluidState&)
864 {}
865
866 template <class BoundaryFluidState,class IntensiveQuantities>
867 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
868 const IntensiveQuantities& /*inIq*/,
869 unsigned /*focusDofIndex*/,
870 unsigned /*inIdx*/,
871 unsigned /*timeIdx*/,
872 Scalar /*alpha*/,
873 const BoundaryFluidState& /*boundaryFs*/)
874 {}
875
876 const Evaluation& energyFlux() const
877 { throw std::logic_error("Requested the energy flux, but energy is not conserved"); }
878};
879
880template <class TypeTag>
881class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::SequentialImplicitThermal>
882{
886
887public:
888 template<class Evaluation, class FluidState, class IntensiveQuantities>
889 static void updateEnergy(Evaluation& energyFlux,
890 const unsigned& focusDofIndex,
891 const unsigned& inIdx,
892 const unsigned& exIdx,
893 const IntensiveQuantities& inIq,
894 const IntensiveQuantities& exIq,
895 const FluidState& inFs,
896 const FluidState& exFs,
897 const Scalar& inAlpha,
898 const Scalar& outAlpha,
899 const Scalar& faceArea)
900 {
901 // Uses the same caculations as the FullyImplicitThermal approach
903 focusDofIndex,
904 inIdx,
905 exIdx,
906 inIq,
907 exIq,
908 inFs,
909 exFs,
910 inAlpha,
911 outAlpha,
912 faceArea);
913 }
914
915 void updateEnergy(const ElementContext&,
916 unsigned,
917 unsigned)
918 { } // Old interface still used output code for fluxes. But energy flux is not used. i.e. do nothing
919
920 template <class Context, class BoundaryFluidState>
921 void updateEnergyBoundary(const Context&,
922 unsigned,
923 unsigned,
924 const BoundaryFluidState&)
925 { }
926
927 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
928 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
929 const IntensiveQuantities& /*inIq*/,
930 unsigned /*focusDofIndex*/,
931 unsigned /*inIdx*/,
932 unsigned /*timeIdx*/,
933 Scalar /*alpha*/,
934 const BoundaryFluidState& /*boundaryFs*/)
935 { }
936
937 const Evaluation& energyFlux() const
938 { throw std::logic_error("Requested the energy flux, but energy is not conserved"); }
939};
940
941template <class TypeTag>
942class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::NoTemperature>
943{
947
948public:
949 template<class Evaluation, class FluidState, class IntensiveQuantities>
950 static void updateEnergy(Evaluation& /*energyFlux*/,
951 const unsigned& /*focusDofIndex*/,
952 const unsigned& /*inIdx*/,
953 const unsigned& /*exIdx*/,
954 const IntensiveQuantities& /*inIq*/,
955 const IntensiveQuantities& /*exIq*/,
956 const FluidState& /*inFs*/,
957 const FluidState& /*exFs*/,
958 const Scalar& /*inAlpha*/,
959 const Scalar& /*outAlpha*/,
960 const Scalar& /*faceArea*/)
961 {}
962
963 void updateEnergy(const ElementContext&,
964 unsigned,
965 unsigned)
966 {}
967
968 template <class Context, class BoundaryFluidState>
969 void updateEnergyBoundary(const Context&,
970 unsigned,
971 unsigned,
972 const BoundaryFluidState&)
973 {}
974
975 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
976 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
977 const IntensiveQuantities& /*inIq*/,
978 unsigned /*focusDofIndex*/,
979 unsigned /*inIdx*/,
980 unsigned /*timeIdx*/,
981 Scalar /*alpha*/,
982 const BoundaryFluidState& /*boundaryFs*/)
983 {}
984
985 const Evaluation& energyFlux() const
986 { throw std::logic_error("Requested the energy flux, but energy is not conserved"); }
987};
988
989} // namespace Opm
990
991#endif
Declares the properties required by the black oil model.
static void updateEnergy(Evaluation &, const unsigned &, const unsigned &, const unsigned &, const IntensiveQuantities &, const IntensiveQuantities &, const FluidState &, const FluidState &, const Scalar &, const Scalar &, const Scalar &)
Definition: blackoilenergymodules.hh:841
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:860
const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:876
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:854
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:867
static OPM_HOST_DEVICE void updateEnergy(Evaluation &energyFlux, const unsigned &focusDofIndex, const unsigned &inIdx, const unsigned &exIdx, const IntensiveQuantities &inIq, const IntensiveQuantities &exIq, const FluidState &inFs, const FluidState &exFs, const Scalar &inAlpha, const Scalar &outAlpha, const Scalar &faceArea)
Definition: blackoilenergymodules.hh:679
const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:822
void updateEnergyBoundary(const Context &ctx, unsigned scvfIdx, unsigned timeIdx, const BoundaryFluidState &boundaryFs)
Definition: blackoilenergymodules.hh:768
void updateEnergy(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:738
static OPM_HOST_DEVICE void updateEnergyBoundary(Evaluation &energyFlux, const IntensiveQuantities &inIq, unsigned focusDofIndex, unsigned inIdx, Scalar alpha, const BoundaryFluidState &boundaryFs)
Definition: blackoilenergymodules.hh:784
const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:985
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:976
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:969
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:963
static void updateEnergy(Evaluation &, const unsigned &, const unsigned &, const unsigned &, const IntensiveQuantities &, const IntensiveQuantities &, const FluidState &, const FluidState &, const Scalar &, const Scalar &, const Scalar &)
Definition: blackoilenergymodules.hh:950
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:928
static void updateEnergy(Evaluation &energyFlux, const unsigned &focusDofIndex, const unsigned &inIdx, const unsigned &exIdx, const IntensiveQuantities &inIq, const IntensiveQuantities &exIq, const FluidState &inFs, const FluidState &exFs, const Scalar &inAlpha, const Scalar &outAlpha, const Scalar &faceArea)
Definition: blackoilenergymodules.hh:889
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:915
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:921
const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:937
Provides the energy specific extensive quantities to the generic black-oil module's extensive quantit...
Definition: blackoilenergymodules.hh:659
Implementation & asImp_()
Definition: blackoilenergymodules.hh:531
const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:518
void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:487
void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:505
const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:524
void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:495
void updateEnergyQuantities_(const ElementContext &, unsigned, unsigned, const typename FluidSystem::template ParameterCache< Evaluation > &)
Definition: blackoilenergymodules.hh:512
void updateEnergyQuantities_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities needed to handle energy conservation.
Definition: blackoilenergymodules.hh:417
OPM_HOST_DEVICE Scalar rockFraction() const
Definition: blackoilenergymodules.hh:461
OPM_HOST_DEVICE const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:458
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:403
Evaluation totalThermalConductivity_
Definition: blackoilenergymodules.hh:469
BlackOilEnergyIntensiveQuantities(Evaluation rockInternalEnergy, Evaluation totalThermalConductivity, Scalar rockFraction)
Construct the energy intensive quantities for the fully implicit thermal module.
Definition: blackoilenergymodules.hh:373
void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the temperature of the intensive quantity's fluid state.
Definition: blackoilenergymodules.hh:388
Implementation & asImp_()
Definition: blackoilenergymodules.hh:465
OPM_HOST_DEVICE const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:455
void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: blackoilenergymodules.hh:424
void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:627
const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:641
void updateEnergyQuantities_(const ElementContext &, unsigned, unsigned, const typename FluidSystem::template ParameterCache< Evaluation > &)
Definition: blackoilenergymodules.hh:635
const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:647
Implementation & asImp_()
Definition: blackoilenergymodules.hh:654
void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:620
void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: blackoilenergymodules.hh:585
void updateEnergyQuantities_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities needed to handle energy conservation.
Definition: blackoilenergymodules.hh:579
void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:551
void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:566
const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:591
const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:597
void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:558
Provides the volumetric quantities required for the equations needed by the energys extension of the ...
Definition: blackoilenergymodules.hh:342
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:63
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:158
static std::string eqName(unsigned eqIdx)
Definition: blackoilenergymodules.hh:142
static bool eqApplies(unsigned eqIdx)
Definition: blackoilenergymodules.hh:132
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilenergymodules.hh:117
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilenergymodules.hh:316
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition: blackoilenergymodules.hh:309
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilenergymodules.hh:124
static OPM_HOST_DEVICE void addPhaseEnthalpyFluxes_(RateVectorT &flux, unsigned phaseIdx, const Eval &volumeFlux, const FluidState &upFs)
Definition: blackoilenergymodules.hh:232
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilenergymodules.hh:326
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilenergymodules.hh:149
static void registerOutputModules(Model &model, Simulator &simulator)
Register all energy specific VTK and ECL output modules.
Definition: blackoilenergymodules.hh:99
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilenergymodules.hh:107
static OPM_HOST_DEVICE void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:190
GetPropType< TypeTag, Properties::ExtensiveQuantities > ExtensiveQuantities
Definition: blackoilenergymodules.hh:84
static void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the energy specific primary variables to a PrimaryVariables object.
Definition: blackoilenergymodules.hh:273
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:284
static OPM_HOST_DEVICE void addHeatFlux(RateVectorT &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:221
static void registerParameters()
Register all run-time parameters for the black-oil energy module.
Definition: blackoilenergymodules.hh:89
static void addPhaseEnthalpyFlux_(RateVector &flux, unsigned phaseIdx, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:244
static void addToEnthalpyRate(RateVector &flux, const Evaluation &hRate)
Definition: blackoilenergymodules.hh:261
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilenergymodules.hh:297
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:84
Definition: blackoilbioeffectsmodules.hh:45
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