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