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
384
389 OPM_HOST_DEVICE void updateTemperature_(const ElementContext& elemCtx,
390 unsigned dofIdx,
391 unsigned timeIdx)
392 {
393 auto& fs = asImp_().fluidState_;
394 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
395
396 // set temperature
397 fs.setTemperature(priVars.makeEvaluation(temperatureIdx, timeIdx, elemCtx.linearizationType()));
398 }
399
404 OPM_HOST_DEVICE void updateTemperature_([[maybe_unused]] const Problem& problem,
405 const PrimaryVariables& priVars,
406 [[maybe_unused]] unsigned globalDofIdx,
407 const unsigned timeIdx,
408 const LinearizationType& lintype)
409 {
410 auto& fs = asImp_().fluidState_;
411 fs.setTemperature(priVars.makeEvaluation(temperatureIdx, timeIdx, lintype));
412 }
413
418 OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext& elemCtx,
419 unsigned dofIdx,
420 unsigned timeIdx)
421 {
422 updateEnergyQuantities_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
423 }
424
425 OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem& problem,
426 const unsigned globalSpaceIdx,
427 const unsigned timeIdx)
428 {
429 auto& fs = asImp_().fluidState_;
430
431 // compute the specific enthalpy of the fluids, the specific enthalpy of the rock
432 // and the thermal conductivity coefficients
433 for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
434 if (!asImp_().getFluidSystem().phaseIsActive(phaseIdx)) {
435 continue;
436 }
437
438 const auto& h = asImp_().getFluidSystem().enthalpy(fs, phaseIdx, problem.pvtRegionIndex(globalSpaceIdx));
439 fs.setEnthalpy(phaseIdx, h);
440 }
441
442 const auto& solidEnergyLawParams = problem.solidEnergyLawParams(globalSpaceIdx, timeIdx);
443 rockInternalEnergy_ = SolidEnergyLaw::solidInternalEnergy(solidEnergyLawParams, fs);
444
445 const auto& thermalConductionLawParams = problem.thermalConductionLawParams(globalSpaceIdx, timeIdx);
446 totalThermalConductivity_ = ThermalConductionLaw::thermalConductivity(thermalConductionLawParams, fs);
447
448 // Retrieve the rock fraction from the problem
449 // Usually 1 - porosity, but if pvmult is used to modify porosity
450 // we will apply the same multiplier to the rock fraction
451 // i.e. pvmult*(1 - porosity) and thus interpret multpv as a volume
452 // multiplier. This is to avoid negative rock volume for pvmult*porosity > 1
453 rockFraction_ = problem.rockFraction(globalSpaceIdx, timeIdx);
454 }
455
456 OPM_HOST_DEVICE const Evaluation& rockInternalEnergy() const
457 { return rockInternalEnergy_; }
458
459 OPM_HOST_DEVICE const Evaluation& totalThermalConductivity() const
460 { return totalThermalConductivity_; }
461
462 OPM_HOST_DEVICE Scalar rockFraction() const
463 { return rockFraction_; }
464
465protected:
466 OPM_HOST_DEVICE Implementation& asImp_()
467 { return *static_cast<Implementation*>(this); }
468
472};
473
474template <class TypeTag>
475class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::ConstantTemperature>
476{
483
485
486public:
487
488 OPM_HOST_DEVICE void updateTemperature_(const ElementContext& elemCtx,
489 unsigned dofIdx,
490 unsigned timeIdx)
491 {
492 updateTemperature_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
493 }
494
495 template<class Problem>
496 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
497 [[maybe_unused]] const PrimaryVariables& priVars,
498 unsigned globalDofIdx,
499 unsigned timeIdx,
500 [[maybe_unused]] const LinearizationType& lintype)
501 {
502 updateTemperature_(problem, globalDofIdx, timeIdx);
503 }
504
505 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
506 unsigned globalDofIdx,
507 unsigned timeIdx)
508 {
509 auto& fs = asImp_().fluidState_;
510 const Scalar T = problem.temperature(globalDofIdx, timeIdx);
511 fs.setTemperature(T);
512 }
513
514 OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext&,
515 unsigned,
516 unsigned,
517 const typename FluidSystem::template ParameterCache<Evaluation>&)
518 {}
519
520 OPM_HOST_DEVICE const Evaluation& rockInternalEnergy() const
521 {
522 OPM_THROW(std::logic_error,
523 "Requested the rock internal energy, which is "
524 "unavailable because energy is not conserved");
525 }
526
527 OPM_HOST_DEVICE const Evaluation& totalThermalConductivity() const
528 {
529 OPM_THROW(std::logic_error,
530 "Requested the total thermal conductivity, which is "
531 "unavailable because energy is not conserved");
532 }
533
534protected:
535 OPM_HOST_DEVICE Implementation& asImp_()
536 { return *static_cast<Implementation*>(this); }
537};
538
539template <class TypeTag>
540class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::SequentialImplicitThermal>
541{
552
553public:
554
555 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
556 unsigned globalDofIdx,
557 unsigned timeIdx)
558 {
559 // update the temperature for output (without derivatives)
560 auto& fs = asImp_().fluidState_;
561 fs.setTemperature(problem.temperature(globalDofIdx, timeIdx));
562 }
563
564 OPM_HOST_DEVICE void updateTemperature_(const ElementContext& elemCtx,
565 unsigned dofIdx,
566 unsigned timeIdx)
567 {
568 updateTemperature_(elemCtx.problem(), elemCtx.globalSpaceIndex(dofIdx, timeIdx), timeIdx);
569 }
570
571 template<class Problem>
572 OPM_HOST_DEVICE void updateTemperature_(const Problem& problem,
573 [[maybe_unused]] const PrimaryVariables& priVars,
574 unsigned globalDofIdx,
575 unsigned timeIdx,
576 [[maybe_unused]] const LinearizationType& lintype)
577 {
578 updateTemperature_(problem, globalDofIdx, timeIdx);
579 }
580
585 OPM_HOST_DEVICE void updateEnergyQuantities_([[maybe_unused]] const ElementContext& elemCtx,
586 [[maybe_unused]] unsigned dofIdx,
587 [[maybe_unused]] unsigned timeIdx)
588 {
589 }
590
591 OPM_HOST_DEVICE void updateEnergyQuantities_([[maybe_unused]] const Problem& problem,
592 [[maybe_unused]] const unsigned globalSpaceIdx,
593 [[maybe_unused]] const unsigned timeIdx)
594 {
595 }
596
597 OPM_HOST_DEVICE const Evaluation& rockInternalEnergy() const
598 {
599 OPM_THROW(std::logic_error,
600 "Requested the rock internal energy, which is "
601 "unavailable because energy is not conserved");
602 }
603
604 OPM_HOST_DEVICE const Evaluation& totalThermalConductivity() const
605 {
606 OPM_THROW(std::logic_error,
607 "Requested the total thermal conductivity, which is "
608 "unavailable because energy is not conserved");
609 }
610
611protected:
612 OPM_HOST_DEVICE Implementation& asImp_()
613 { return *static_cast<Implementation*>(this); }
614
615};
616
617template <class TypeTag>
618class BlackOilEnergyIntensiveQuantities<TypeTag, EnergyModules::NoTemperature>
619{
626
627public:
628 OPM_HOST_DEVICE void updateTemperature_([[maybe_unused]] const ElementContext& elemCtx,
629 [[maybe_unused]] unsigned dofIdx,
630 [[maybe_unused]] unsigned timeIdx)
631 {
632 }
633
634 template<class Problem>
635 OPM_HOST_DEVICE void updateTemperature_([[maybe_unused]] const Problem& problem,
636 [[maybe_unused]] const PrimaryVariables& priVars,
637 [[maybe_unused]] unsigned globalDofIdx,
638 [[maybe_unused]] unsigned timeIdx,
639 [[maybe_unused]] const LinearizationType& lintype)
640 {
641 }
642
643 OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext&,
644 unsigned,
645 unsigned,
646 const typename FluidSystem::template ParameterCache<Evaluation>&)
647 {}
648
649 OPM_HOST_DEVICE const Evaluation& rockInternalEnergy() const
650 {
651 OPM_THROW(std::logic_error,
652 "Requested the rock internal energy, which is "
653 "unavailable because energy is not conserved");
654 }
655
656 OPM_HOST_DEVICE const Evaluation& totalThermalConductivity() const
657 {
658 OPM_THROW(std::logic_error,
659 "Requested the total thermal conductivity, which is "
660 "unavailable because energy is not conserved");
661 }
662
663protected:
664 OPM_HOST_DEVICE Implementation& asImp_()
665 { return *static_cast<Implementation*>(this); }
666};
667
668template <class TypeTag, EnergyModules activeModule>
670
678template <class TypeTag>
679class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::FullyImplicitThermal>
680{
682
686
687public:
688 template<class Evaluation, class FluidState, class IntensiveQuantities>
689 OPM_HOST_DEVICE static void updateEnergy(Evaluation& energyFlux,
690 const unsigned& focusDofIndex,
691 const unsigned& inIdx,
692 const unsigned& exIdx,
693 const IntensiveQuantities& inIq,
694 const IntensiveQuantities& exIq,
695 const FluidState& inFs,
696 const FluidState& exFs,
697 const Scalar& inAlpha,
698 const Scalar& outAlpha,
699 const Scalar& faceArea)
700 {
701 Evaluation deltaT;
702 if (focusDofIndex == inIdx) {
703 deltaT = decay<Scalar>(exFs.temperature(/*phaseIdx=*/0)) -
704 inFs.temperature(/*phaseIdx=*/0);
705 }
706 else if (focusDofIndex == exIdx) {
707 deltaT = exFs.temperature(/*phaseIdx=*/0) -
708 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
709 }
710 else {
711 deltaT = decay<Scalar>(exFs.temperature(/*phaseIdx=*/0)) -
712 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
713 }
714
715 Evaluation inLambda;
716 if (focusDofIndex == inIdx) {
717 inLambda = inIq.totalThermalConductivity();
718 }
719 else {
720 inLambda = decay<Scalar>(inIq.totalThermalConductivity());
721 }
722
723 Evaluation exLambda;
724 if (focusDofIndex == exIdx) {
725 exLambda = exIq.totalThermalConductivity();
726 }
727 else {
728 exLambda = decay<Scalar>(exIq.totalThermalConductivity());
729 }
730
731 Evaluation H;
732 const Evaluation& inH = inLambda*inAlpha;
733 const Evaluation& exH = exLambda*outAlpha;
734 if (inH > 0 && exH > 0) {
735 // compute the "thermal transmissibility". In contrast to the normal
736 // transmissibility this cannot be done as a preprocessing step because the
737 // average thermal conductivity is analogous to the permeability but
738 // depends on the solution.
739 H = 1.0 / (1.0 / inH + 1.0 / exH);
740 }
741 else {
742 H = 0.0;
743 }
744
745 energyFlux = deltaT * (-H / faceArea);
746 }
747
748 void updateEnergy(const ElementContext& elemCtx,
749 unsigned scvfIdx,
750 unsigned timeIdx)
751 {
752 const auto& stencil = elemCtx.stencil(timeIdx);
753 const auto& scvf = stencil.interiorFace(scvfIdx);
754
755 const Scalar faceArea = scvf.area();
756 const unsigned inIdx = scvf.interiorIndex();
757 const unsigned exIdx = scvf.exteriorIndex();
758 const auto& inIq = elemCtx.intensiveQuantities(inIdx, timeIdx);
759 const auto& exIq = elemCtx.intensiveQuantities(exIdx, timeIdx);
760 const auto& inFs = inIq.fluidState();
761 const auto& exFs = exIq.fluidState();
762 const Scalar inAlpha = elemCtx.problem().thermalHalfTransmissibilityIn(elemCtx, scvfIdx, timeIdx);
763 const Scalar outAlpha = elemCtx.problem().thermalHalfTransmissibilityOut(elemCtx, scvfIdx, timeIdx);
764 updateEnergy(energyFlux_,
765 elemCtx.focusDofIndex(),
766 inIdx,
767 exIdx,
768 inIq,
769 exIq,
770 inFs,
771 exFs,
772 inAlpha,
773 outAlpha,
774 faceArea);
775 }
776
777 template <class Context, class BoundaryFluidState>
778 void updateEnergyBoundary(const Context& ctx,
779 unsigned scvfIdx,
780 unsigned timeIdx,
781 const BoundaryFluidState& boundaryFs)
782 {
783 const auto& stencil = ctx.stencil(timeIdx);
784 const auto& scvf = stencil.boundaryFace(scvfIdx);
785
786 const unsigned inIdx = scvf.interiorIndex();
787 const auto& inIq = ctx.intensiveQuantities(inIdx, timeIdx);
788 const auto& focusDofIdx = ctx.focusDofIndex();
789 const Scalar alpha = ctx.problem().thermalHalfTransmissibilityBoundary(ctx, scvfIdx);
790 updateEnergyBoundary(energyFlux_, inIq, focusDofIdx, inIdx, alpha, boundaryFs);
791 }
792
793 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
794 OPM_HOST_DEVICE static void updateEnergyBoundary(Evaluation& energyFlux,
795 const IntensiveQuantities& inIq,
796 unsigned focusDofIndex,
797 unsigned inIdx,
798 Scalar alpha,
799 const BoundaryFluidState& boundaryFs)
800 {
801 const auto& inFs = inIq.fluidState();
802 Evaluation deltaT;
803 if (focusDofIndex == inIdx) {
804 deltaT = boundaryFs.temperature(/*phaseIdx=*/0) -
805 inFs.temperature(/*phaseIdx=*/0);
806 }
807 else {
808 deltaT = decay<Scalar>(boundaryFs.temperature(/*phaseIdx=*/0)) -
809 decay<Scalar>(inFs.temperature(/*phaseIdx=*/0));
810 }
811
812 Evaluation lambda;
813 if (focusDofIndex == inIdx) {
814 lambda = inIq.totalThermalConductivity();
815 }
816 else {
817 lambda = decay<Scalar>(inIq.totalThermalConductivity());
818 }
819
820 if (lambda > 0.0) {
821 // compute the "thermal transmissibility". In contrast to the normal
822 // transmissibility this cannot be done as a preprocessing step because the
823 // average thermal conductivity is analogous to the permeability but depends
824 // on the solution.
825 energyFlux = deltaT * lambda * -alpha;
826 }
827 else {
828 energyFlux = 0.0;
829 }
830 }
831
832 const Evaluation& energyFlux() const
833 { return energyFlux_; }
834
835private:
836 Implementation& asImp_()
837 { return *static_cast<Implementation*>(this); }
838
839 Evaluation energyFlux_;
840};
841
842template <class TypeTag>
843class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::ConstantTemperature>
844{
848
849public:
850 template<class Evaluation, class FluidState, class IntensiveQuantities>
851 static void updateEnergy(Evaluation& /*energyFlux*/,
852 const unsigned& /*focusDofIndex*/,
853 const unsigned& /*inIdx*/,
854 const unsigned& /*exIdx*/,
855 const IntensiveQuantities& /*inIq*/,
856 const IntensiveQuantities& /*exIq*/,
857 const FluidState& /*inFs*/,
858 const FluidState& /*exFs*/,
859 const Scalar& /*inAlpha*/,
860 const Scalar& /*outAlpha*/,
861 const Scalar& /*faceArea*/)
862 {}
863
864 void updateEnergy(const ElementContext&,
865 unsigned,
866 unsigned)
867 {}
868
869 template <class Context, class BoundaryFluidState>
870 void updateEnergyBoundary(const Context&,
871 unsigned,
872 unsigned,
873 const BoundaryFluidState&)
874 {}
875
876 template <class BoundaryFluidState,class IntensiveQuantities>
877 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
878 const IntensiveQuantities& /*inIq*/,
879 unsigned /*focusDofIndex*/,
880 unsigned /*inIdx*/,
881 unsigned /*timeIdx*/,
882 Scalar /*alpha*/,
883 const BoundaryFluidState& /*boundaryFs*/)
884 {}
885
886 OPM_HOST_DEVICE const Evaluation& energyFlux() const
887 { OPM_THROW(std::logic_error, "Requested the energy flux, but energy is not conserved"); }
888};
889
890template <class TypeTag>
891class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::SequentialImplicitThermal>
892{
896
897public:
898 template<class Evaluation, class FluidState, class IntensiveQuantities>
899 static void updateEnergy(Evaluation& energyFlux,
900 const unsigned& focusDofIndex,
901 const unsigned& inIdx,
902 const unsigned& exIdx,
903 const IntensiveQuantities& inIq,
904 const IntensiveQuantities& exIq,
905 const FluidState& inFs,
906 const FluidState& exFs,
907 const Scalar& inAlpha,
908 const Scalar& outAlpha,
909 const Scalar& faceArea)
910 {
911 // Uses the same caculations as the FullyImplicitThermal approach
913 focusDofIndex,
914 inIdx,
915 exIdx,
916 inIq,
917 exIq,
918 inFs,
919 exFs,
920 inAlpha,
921 outAlpha,
922 faceArea);
923 }
924
925 void updateEnergy(const ElementContext&,
926 unsigned,
927 unsigned)
928 { } // Old interface still used output code for fluxes. But energy flux is not used. i.e. do nothing
929
930 template <class Context, class BoundaryFluidState>
931 void updateEnergyBoundary(const Context&,
932 unsigned,
933 unsigned,
934 const BoundaryFluidState&)
935 { }
936
937 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
938 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
939 const IntensiveQuantities& /*inIq*/,
940 unsigned /*focusDofIndex*/,
941 unsigned /*inIdx*/,
942 unsigned /*timeIdx*/,
943 Scalar /*alpha*/,
944 const BoundaryFluidState& /*boundaryFs*/)
945 { }
946
947 OPM_HOST_DEVICE const Evaluation& energyFlux() const
948 { OPM_THROW(std::logic_error, "Requested the energy flux, but energy is not conserved"); }
949};
950
951template <class TypeTag>
952class BlackOilEnergyExtensiveQuantities<TypeTag, EnergyModules::NoTemperature>
953{
957
958public:
959 template<class Evaluation, class FluidState, class IntensiveQuantities>
960 static void updateEnergy(Evaluation& /*energyFlux*/,
961 const unsigned& /*focusDofIndex*/,
962 const unsigned& /*inIdx*/,
963 const unsigned& /*exIdx*/,
964 const IntensiveQuantities& /*inIq*/,
965 const IntensiveQuantities& /*exIq*/,
966 const FluidState& /*inFs*/,
967 const FluidState& /*exFs*/,
968 const Scalar& /*inAlpha*/,
969 const Scalar& /*outAlpha*/,
970 const Scalar& /*faceArea*/)
971 {}
972
973 void updateEnergy(const ElementContext&,
974 unsigned,
975 unsigned)
976 {}
977
978 template <class Context, class BoundaryFluidState>
979 void updateEnergyBoundary(const Context&,
980 unsigned,
981 unsigned,
982 const BoundaryFluidState&)
983 {}
984
985 template <class Evaluation, class BoundaryFluidState, class IntensiveQuantities>
986 static void updateEnergyBoundary(Evaluation& /*heatFlux*/,
987 const IntensiveQuantities& /*inIq*/,
988 unsigned /*focusDofIndex*/,
989 unsigned /*inIdx*/,
990 unsigned /*timeIdx*/,
991 Scalar /*alpha*/,
992 const BoundaryFluidState& /*boundaryFs*/)
993 {}
994
995 OPM_HOST_DEVICE const Evaluation& energyFlux() const
996 { OPM_THROW(std::logic_error, "Requested the energy flux, but energy is not conserved"); }
997};
998
999} // namespace Opm
1000
1001#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:851
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:870
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:864
OPM_HOST_DEVICE const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:886
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:877
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:689
const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:832
void updateEnergyBoundary(const Context &ctx, unsigned scvfIdx, unsigned timeIdx, const BoundaryFluidState &boundaryFs)
Definition: blackoilenergymodules.hh:778
void updateEnergy(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:748
static OPM_HOST_DEVICE void updateEnergyBoundary(Evaluation &energyFlux, const IntensiveQuantities &inIq, unsigned focusDofIndex, unsigned inIdx, Scalar alpha, const BoundaryFluidState &boundaryFs)
Definition: blackoilenergymodules.hh:794
OPM_HOST_DEVICE const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:995
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:986
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:979
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:973
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:960
static void updateEnergyBoundary(Evaluation &, const IntensiveQuantities &, unsigned, unsigned, unsigned, Scalar, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:938
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:899
void updateEnergy(const ElementContext &, unsigned, unsigned)
Definition: blackoilenergymodules.hh:925
OPM_HOST_DEVICE const Evaluation & energyFlux() const
Definition: blackoilenergymodules.hh:947
void updateEnergyBoundary(const Context &, unsigned, unsigned, const BoundaryFluidState &)
Definition: blackoilenergymodules.hh:931
Provides the energy specific extensive quantities to the generic black-oil module's extensive quantit...
Definition: blackoilenergymodules.hh:669
OPM_HOST_DEVICE const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:527
OPM_HOST_DEVICE Implementation & asImp_()
Definition: blackoilenergymodules.hh:535
OPM_HOST_DEVICE const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:520
OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext &, unsigned, unsigned, const typename FluidSystem::template ParameterCache< Evaluation > &)
Definition: blackoilenergymodules.hh:514
OPM_HOST_DEVICE void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:488
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:496
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:505
OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: blackoilenergymodules.hh:425
OPM_HOST_DEVICE Scalar rockFraction() const
Definition: blackoilenergymodules.hh:462
OPM_HOST_DEVICE const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:459
OPM_HOST_DEVICE Implementation & asImp_()
Definition: blackoilenergymodules.hh:466
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:404
OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities needed to handle energy conservation.
Definition: blackoilenergymodules.hh:418
Evaluation totalThermalConductivity_
Definition: blackoilenergymodules.hh:470
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:389
OPM_HOST_DEVICE const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:456
OPM_HOST_DEVICE Implementation & asImp_()
Definition: blackoilenergymodules.hh:664
OPM_HOST_DEVICE const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:649
OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext &, unsigned, unsigned, const typename FluidSystem::template ParameterCache< Evaluation > &)
Definition: blackoilenergymodules.hh:643
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:635
OPM_HOST_DEVICE void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:628
OPM_HOST_DEVICE const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:656
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, const PrimaryVariables &priVars, unsigned globalDofIdx, unsigned timeIdx, const LinearizationType &lintype)
Definition: blackoilenergymodules.hh:572
OPM_HOST_DEVICE const Evaluation & totalThermalConductivity() const
Definition: blackoilenergymodules.hh:604
OPM_HOST_DEVICE void updateEnergyQuantities_(const Problem &problem, const unsigned globalSpaceIdx, const unsigned timeIdx)
Definition: blackoilenergymodules.hh:591
OPM_HOST_DEVICE const Evaluation & rockInternalEnergy() const
Definition: blackoilenergymodules.hh:597
OPM_HOST_DEVICE void updateTemperature_(const Problem &problem, unsigned globalDofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:555
OPM_HOST_DEVICE Implementation & asImp_()
Definition: blackoilenergymodules.hh:612
OPM_HOST_DEVICE void updateTemperature_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:564
OPM_HOST_DEVICE void updateEnergyQuantities_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Compute the intensive quantities needed to handle energy conservation.
Definition: blackoilenergymodules.hh:585
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