blackoillocalresidualtpfa.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_LOCAL_TPFA_RESIDUAL_HH
29#define EWOMS_BLACK_OIL_LOCAL_TPFA_RESIDUAL_HH
30
31#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
32#include <opm/input/eclipse/Schedule/BCProp.hpp>
33
34#include <opm/material/common/ConditionalStorage.hpp>
35#include <opm/material/common/MathToolbox.hpp>
36#include <opm/material/fluidstates/BlackOilFluidState.hpp>
37
50
51#include <opm/common/ErrorMacros.hpp>
52#include <opm/common/utility/gpuDecorators.hpp>
53
54#include <array>
55#include <cassert>
56#include <stdexcept>
57#include <string>
58
59#include <opm/common/utility/gpuistl_if_available.hpp>
60
61namespace Opm
62{
69template <class TypeTag>
70class BlackOilLocalResidualTPFA : public GetPropType<TypeTag, Properties::DiscLocalResidual>
71{
82 using FluidState = typename IntensiveQuantities::FluidState;
83
84 enum { conti0EqIdx = Indices::conti0EqIdx };
85 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
86 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
87 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
88
89 enum { dimWorld = GridView::dimensionworld };
90 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
91 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
92 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
93
94 enum { gasCompIdx = FluidSystem::gasCompIdx };
95 enum { oilCompIdx = FluidSystem::oilCompIdx };
96 enum { waterCompIdx = FluidSystem::waterCompIdx };
97 enum { compositionSwitchIdx = Indices::compositionSwitchIdx };
98
99 static constexpr bool waterEnabled = Indices::waterEnabled;
100 static constexpr bool gasEnabled = Indices::gasEnabled;
101 static constexpr bool oilEnabled = Indices::oilEnabled;
102 static constexpr bool compositionSwitchEnabled = compositionSwitchIdx >= 0;
103
104 static constexpr bool blackoilConserveSurfaceVolume
105 = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
106
107 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
108 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
109 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
110 static constexpr bool enableFullyImplicitThermal
111 = (getPropValue<TypeTag, Properties::EnergyModuleType>()
112 == EnergyModules::FullyImplicitThermal);
113 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
114 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
115 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
116 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
117 static constexpr bool enableConvectiveMixing
118 = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
119 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
120 static constexpr bool enableSaltPrecipitation
121 = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
122 static constexpr bool enableMICP = Indices::enableMICP;
123 static constexpr bool runAssemblyOnGpu = getPropValue<TypeTag, Properties::RunAssemblyOnGpu>();
124
134
137
138 using Toolbox = MathToolbox<Evaluation>;
139
140public:
142 double trans;
143 double faceArea;
144 double thpres;
145 double dZg;
146 FaceDir::DirEnum faceDir;
147 double Vin;
148 double Vex;
149 ConditionalStorage<enableFullyImplicitThermal, double> inAlpha;
150 ConditionalStorage<enableFullyImplicitThermal, double> outAlpha;
151 ConditionalStorage<enableDiffusion, double> diffusivity;
152 ConditionalStorage<enableDispersion, double> dispersivity;
153 };
154
158 template <class LhsEval>
159 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
160 const ElementContext& elemCtx,
161 unsigned dofIdx,
162 unsigned timeIdx) const
163 {
164 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
165 computeStorage<LhsEval>(storage, intQuants);
166 }
167
168 template <class LhsEval, class StorageType, class IntensiveQuantitiesType = IntensiveQuantities>
169 OPM_HOST_DEVICE static void computeStorage(StorageType& storage,
170 const IntensiveQuantitiesType& intQuants)
171 {
172 OPM_TIMEBLOCK_LOCAL(computeStorage, Subsystem::Assembly);
173 // retrieve the intensive quantities for the SCV at the specified point in time
174 const auto& fs = intQuants.fluidState();
175 storage = 0.0;
176
177 const FluidSystem& fsys = intQuants.getFluidSystem();
178
179 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
180 if (!fsys.phaseIsActive(phaseIdx)) {
181 continue;
182 }
183 unsigned activeCompIdx
184 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
185 LhsEval surfaceVolume = Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
186 * Toolbox::template decay<LhsEval>(fs.invB(phaseIdx))
187 * Toolbox::template decay<LhsEval>(intQuants.porosity());
188
189 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
190
191 // account for dissolved gas
192 if (phaseIdx == oilPhaseIdx && fsys.enableDissolvedGas()) {
193 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
194 storage[conti0EqIdx + activeGasCompIdx]
195 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs())
196 * surfaceVolume;
197 }
198
199 // account for dissolved gas in water
200 if (phaseIdx == waterPhaseIdx && fsys.enableDissolvedGasInWater()) {
201 unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
202 storage[conti0EqIdx + activeGasCompIdx]
203 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw())
204 * surfaceVolume;
205 }
206
207 // account for vaporized oil
208 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedOil()) {
209 unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
210 storage[conti0EqIdx + activeOilCompIdx]
211 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv())
212 * surfaceVolume;
213 }
214
215 // account for vaporized water
216 if (phaseIdx == gasPhaseIdx && fsys.enableVaporizedWater()) {
217 unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
218 storage[conti0EqIdx + activeWaterCompIdx]
219 += Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw())
220 * surfaceVolume;
221 }
222 }
223
224 adaptMassConservationQuantities_(storage, intQuants.pvtRegionIndex(), fsys);
225
226 // deal with solvents (if present)
227 SolventModule::addStorage(storage, intQuants);
228
229 // deal with zFracton (if present)
230 ExtboModule::addStorage(storage, intQuants);
231
232 // deal with polymer (if present)
233 PolymerModule::addStorage(storage, intQuants);
234
235 // deal with energy (if present)
236 EnergyModule::addStorage(storage, intQuants);
237
238 // deal with foam (if present)
239 FoamModule::addStorage(storage, intQuants);
240
241 // deal with salt (if present)
242 BrineModule::addStorage(storage, intQuants);
243
244 // deal with bioeffects (if present)
245 BioeffectsModule::addStorage(storage, intQuants);
246 }
247
253 template <class ModuleParamsT,
254 class RateVectorT,
255 class IntensiveQuantitiesT,
256 class ResidualNBInfoT>
257 OPM_HOST_DEVICE static void computeFlux(RateVectorT& flux,
258 RateVectorT& darcy,
259 const unsigned globalIndexIn,
260 const unsigned globalIndexEx,
261 const IntensiveQuantitiesT& intQuantsIn,
262 const IntensiveQuantitiesT& intQuantsEx,
263 const ResidualNBInfoT& nbInfo,
264 const ModuleParamsT& moduleParams)
265 {
266 OPM_TIMEBLOCK_LOCAL(computeFlux, Subsystem::Assembly);
267 flux = 0.0;
268 darcy = 0.0;
269
270 calculateFluxes_(flux,
271 darcy,
272 intQuantsIn,
273 intQuantsEx,
274 globalIndexIn,
275 globalIndexEx,
276 nbInfo,
277 moduleParams);
278 }
279
280 // This function demonstrates compatibility with the ElementContext-based interface.
281 // Actually using it will lead to double work since the element context already contains
282 // fluxes through its stored ExtensiveQuantities.
283 static void
284 computeFlux(RateVector& flux, const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
285 {
286 OPM_TIMEBLOCK_LOCAL(computeFlux, Subsystem::Assembly);
287 assert(timeIdx == 0);
288
289 flux = 0.0;
290 RateVector darcy = 0.0;
291 // need for dary flux calculation
292 const auto& problem = elemCtx.problem();
293 const auto& stencil = elemCtx.stencil(timeIdx);
294 const auto& scvf = stencil.interiorFace(scvfIdx);
295
296 unsigned interiorDofIdx = scvf.interiorIndex();
297 unsigned exteriorDofIdx = scvf.exteriorIndex();
298 assert(interiorDofIdx != exteriorDofIdx);
299
300 // unsigned I = stencil.globalSpaceIndex(interiorDofIdx);
301 // unsigned J = stencil.globalSpaceIndex(exteriorDofIdx);
302 Scalar Vin = elemCtx.dofVolume(interiorDofIdx, /*timeIdx=*/0);
303 Scalar Vex = elemCtx.dofVolume(exteriorDofIdx, /*timeIdx=*/0);
304 const auto& globalIndexIn = stencil.globalSpaceIndex(interiorDofIdx);
305 const auto& globalIndexEx = stencil.globalSpaceIndex(exteriorDofIdx);
306 Scalar trans = problem.transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
307 Scalar faceArea = scvf.area();
308 const auto faceDir = faceDirFromDirId(scvf.dirId());
309 Scalar thpres = problem.thresholdPressure(globalIndexIn, globalIndexEx);
310
311 // estimate the gravity correction: for performance reasons we use a simplified
312 // approach for this flux module that assumes that gravity is constant and always
313 // acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
314 const Scalar g = problem.gravity()[dimWorld - 1];
315 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
316 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
317
318 // this is quite hacky because the dune grid interface does not provide a
319 // cellCenterDepth() method (so we ask the problem to provide it). The "good"
320 // solution would be to take the Z coordinate of the element centroids, but since
321 // ECL seems to like to be inconsistent on that front, it needs to be done like
322 // here...
323 const Scalar zIn = problem.dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
324 const Scalar zEx = problem.dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
325 // the distances from the DOF's depths. (i.e., the additional depth of the
326 // exterior DOF)
327 const Scalar distZ = zIn - zEx;
328 // for thermal harmonic mean of half trans
329 const Scalar inAlpha = problem.thermalHalfTransmissibility(globalIndexIn, globalIndexEx);
330 const Scalar outAlpha = problem.thermalHalfTransmissibility(globalIndexEx, globalIndexIn);
331 const Scalar diffusivity = problem.diffusivity(globalIndexEx, globalIndexIn);
332 const Scalar dispersivity = problem.dispersivity(globalIndexEx, globalIndexIn);
333
334 const ResidualNBInfo res_nbinfo {trans,
335 faceArea,
336 thpres,
337 distZ * g,
338 faceDir,
339 Vin,
340 Vex,
341 inAlpha,
342 outAlpha,
343 diffusivity,
344 dispersivity};
345
346 calculateFluxes_(flux,
347 darcy,
348 intQuantsIn,
349 intQuantsEx,
350 globalIndexIn,
351 globalIndexEx,
352 res_nbinfo,
353 problem.moduleParams());
354 }
355
356 template <class RateVectorT,
357 class IntensiveQuantitiesT,
358 class ResidualNBInfoT,
359 class ModuleParamsT>
360 OPM_HOST_DEVICE static void calculateFluxes_(RateVectorT& flux,
361 RateVectorT& darcy,
362 const IntensiveQuantitiesT& intQuantsIn,
363 const IntensiveQuantitiesT& intQuantsEx,
364 const unsigned& globalIndexIn,
365 const unsigned& globalIndexEx,
366 const ResidualNBInfoT& nbInfo,
367 const ModuleParamsT& moduleParams)
368 {
369 OPM_TIMEBLOCK_LOCAL(calculateFluxes, Subsystem::Assembly);
370 const Scalar Vin = nbInfo.Vin;
371 const Scalar Vex = nbInfo.Vex;
372 const Scalar distZg = nbInfo.dZg;
373 const Scalar thpres = nbInfo.thpres;
374 const Scalar trans = nbInfo.trans;
375 const Scalar faceArea = nbInfo.faceArea;
376 FaceDir::DirEnum facedir = nbInfo.faceDir;
377
378 const FluidSystem& fsys = intQuantsIn.getFluidSystem();
379
380 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
381 if (!fsys.phaseIsActive(phaseIdx)) {
382 continue;
383 }
384 // darcy flux calculation
385 short dnIdx;
386 //
387 short upIdx;
388 // fake intices should only be used to get upwind anc compatibility with old functions
389 short interiorDofIdx = 0; // NB
390 short exteriorDofIdx = 1; // NB
391 Evaluation pressureDifference;
392 ExtensiveQuantities::calculatePhasePressureDiff_(upIdx,
393 dnIdx,
394 pressureDifference,
395 intQuantsIn,
396 intQuantsEx,
397 phaseIdx, // input
398 interiorDofIdx, // input
399 exteriorDofIdx, // input
400 Vin,
401 Vex,
402 globalIndexIn,
403 globalIndexEx,
404 distZg,
405 thpres,
406 moduleParams);
407
408 const IntensiveQuantities& up = (upIdx == interiorDofIdx) ? intQuantsIn : intQuantsEx;
409 unsigned globalUpIndex = (upIdx == interiorDofIdx) ? globalIndexIn : globalIndexEx;
410 // Use arithmetic average (more accurate with harmonic, but that requires recomputing
411 // the transmissbility)
412 Evaluation transMult = (intQuantsIn.rockCompTransMultiplier()
413 + Toolbox::value(intQuantsEx.rockCompTransMultiplier()))
414 / 2;
415 if constexpr (enableBioeffects || enableSaltPrecipitation) {
416 transMult
417 *= (intQuantsIn.permFactor() + Toolbox::value(intQuantsEx.permFactor())) / 2;
418 }
419
420 Evaluation darcyFlux;
421 if (globalUpIndex == globalIndexIn) {
422 darcyFlux = pressureDifference * up.mobility(phaseIdx, facedir) * transMult
423 * (-trans / faceArea);
424 } else {
425 darcyFlux = pressureDifference
426 * (Toolbox::value(up.mobility(phaseIdx, facedir)) * transMult
427 * (-trans / faceArea));
428 }
429
430 unsigned activeCompIdx
431 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
432 // NB! For the FLORES fluxes without derivatives
433 darcy[conti0EqIdx + activeCompIdx] = darcyFlux.value() * faceArea;
434
435 unsigned pvtRegionIdx = up.pvtRegionIndex();
436 // if (upIdx == globalFocusDofIdx){
437 if (globalUpIndex == globalIndexIn) {
438 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
439 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
440 const auto& surfaceVolumeFlux = invB * darcyFlux;
441
442 evalPhaseFluxes_<Evaluation>(
443 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
444 if constexpr (enableFullyImplicitThermal) {
445 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
446 flux, phaseIdx, darcyFlux, up.fluidState());
447 }
448 if constexpr (enableBioeffects) {
449 BioeffectsModule::template addBioeffectsFluxes_<Evaluation>(
450 flux, phaseIdx, darcyFlux, up);
451 }
452 if constexpr (enableBrine) {
453 BrineModule::template addBrineFluxes_<Evaluation, FluidState>(
454 flux, phaseIdx, darcyFlux, up.fluidState());
455 }
456 } else {
457 const auto& invB = getInvB_<FluidSystem, FluidState, Scalar>(
458 up.fluidState(), phaseIdx, pvtRegionIdx, fsys);
459 const auto& surfaceVolumeFlux = invB * darcyFlux;
460 evalPhaseFluxes_<Scalar>(
461 flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, up.fluidState());
462 if constexpr (enableFullyImplicitThermal) {
463 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
464 flux, phaseIdx, darcyFlux, up.fluidState());
465 }
466 if constexpr (enableBioeffects) {
467 BioeffectsModule::template addBioeffectsFluxes_<Scalar>(
468 flux, phaseIdx, darcyFlux, up);
469 }
470 if constexpr (enableBrine) {
471 BrineModule::template addBrineFluxes_<Scalar, FluidState>(
472 flux, phaseIdx, darcyFlux, up.fluidState());
473 }
474 }
475 }
476
477 // deal with solvents (if present)
478 static_assert(
479 !enableSolvent,
480 "Relevant computeFlux() method must be implemented for this module before enabling.");
481 // SolventModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
482
483 // deal with zFracton (if present)
484 static_assert(
485 !enableExtbo,
486 "Relevant computeFlux() method must be implemented for this module before enabling.");
487 // ExtboModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
488
489 // deal with polymer (if present)
490 static_assert(
491 !enablePolymer,
492 "Relevant computeFlux() method must be implemented for this module before enabling.");
493 // PolymerModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
494
495 // deal with convective mixing
496 if constexpr (enableConvectiveMixing) {
497 ConvectiveMixingModule::addConvectiveMixingFlux(
498 flux,
499 intQuantsIn,
500 intQuantsEx,
501 globalIndexIn,
502 globalIndexEx,
503 nbInfo.dZg,
504 nbInfo.trans,
505 nbInfo.faceArea,
506 moduleParams.convectiveMixingModuleParam);
507 }
508
509 // deal with energy (if present)
510 if constexpr (enableFullyImplicitThermal) {
511 const Scalar inAlpha = nbInfo.inAlpha;
512 const Scalar outAlpha = nbInfo.outAlpha;
513 Evaluation heatFlux;
514
515 short interiorDofIdx = 0; // NB
516 short exteriorDofIdx = 1; // NB
517
518 EnergyModule::ExtensiveQuantities::updateEnergy(heatFlux,
519 interiorDofIdx, // focusDofIndex,
520 interiorDofIdx,
521 exteriorDofIdx,
522 intQuantsIn,
523 intQuantsEx,
524 intQuantsIn.fluidState(),
525 intQuantsEx.fluidState(),
526 inAlpha,
527 outAlpha,
528 faceArea);
529 EnergyModule::addHeatFlux(flux, heatFlux);
530 }
531 // NB need to be tha last energy call since it does scaling
532 // EnergyModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
533
534 // deal with foam (if present)
535 static_assert(
536 !enableFoam,
537 "Relevant computeFlux() method must be implemented for this module before enabling.");
538 // FoamModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
539
540 // deal with diffusion (if present). opm-models expects per area flux (added in the
541 // tmpdiffusivity).
542 if constexpr (enableDiffusion) {
543 typename DiffusionModule::ExtensiveQuantities::EvaluationArray
544 effectiveDiffusionCoefficient;
545 DiffusionModule::ExtensiveQuantities::update(
546 effectiveDiffusionCoefficient, intQuantsIn, intQuantsEx);
547 const Scalar diffusivity = nbInfo.diffusivity;
548 const Scalar tmpdiffusivity = diffusivity / faceArea;
549 DiffusionModule::addDiffusiveFlux(
550 flux, intQuantsIn, intQuantsEx, tmpdiffusivity, effectiveDiffusionCoefficient);
551 }
552
553 // deal with dispersion (if present). opm-models expects per area flux (added in the
554 // tmpdispersivity).
555 if constexpr (enableDispersion) {
556 typename DispersionModule::ExtensiveQuantities::ScalarArray normVelocityAvg;
557 DispersionModule::ExtensiveQuantities::update(
558 normVelocityAvg, intQuantsIn, intQuantsEx);
559 const Scalar dispersivity = nbInfo.dispersivity;
560 const Scalar tmpdispersivity = dispersivity / faceArea;
561 DispersionModule::addDispersiveFlux(
562 flux, intQuantsIn, intQuantsEx, tmpdispersivity, normVelocityAvg);
563 }
564
565 // apply the scaling for the urea equation in MICP
566 if constexpr (enableMICP) {
568 }
569 }
570
571 template <class BoundaryConditionData, class RateVectorLocal, class LocalProblem>
572 OPM_HOST_DEVICE static void computeBoundaryFlux(RateVectorLocal& bdyFlux,
573 const LocalProblem& problem,
574 const BoundaryConditionData& bdyInfo,
575 const IntensiveQuantities& insideIntQuants,
576 unsigned globalSpaceIdx)
577 {
578#if OPM_IS_INSIDE_HOST_FUNCTION
579 switch (bdyInfo.type) {
580 case BCType::NONE:
581 bdyFlux = 0.0;
582 break;
583 case BCType::RATE:
584 computeBoundaryFluxRate(bdyFlux, bdyInfo);
585 break;
586 case BCType::FREE:
587 case BCType::DIRICHLET:
588 computeBoundaryFluxFree(problem, bdyFlux, bdyInfo, insideIntQuants, globalSpaceIdx);
589 break;
590 case BCType::THERMAL:
591 computeBoundaryThermal(problem, bdyFlux, bdyInfo, insideIntQuants, globalSpaceIdx);
592 break;
593 default:
594 throw std::logic_error("Unknown boundary condition type "
595 + std::to_string(static_cast<int>(bdyInfo.type))
596 + " in computeBoundaryFlux().");
597 }
598#else // TODO: support all boundary conditions on GPU as well to unify this code
599 switch (bdyInfo.type) {
600 case BCType::NONE:
601 bdyFlux = 0.0;
602 break;
603 case BCType::THERMAL:
604 computeBoundaryThermal(problem, bdyFlux, bdyInfo, insideIntQuants, globalSpaceIdx);
605 break;
606 default:
607 OPM_THROW(std::logic_error,
608 "Boundary condition type " + std::to_string(static_cast<int>(bdyInfo.type))
609 + " is not supported for GPU fluid systems in computeBoundaryFlux().");
610 }
611#endif
612 }
613
614 template <class BoundaryConditionData>
615 static void computeBoundaryFluxRate(RateVector& bdyFlux, const BoundaryConditionData& bdyInfo)
616 {
617 bdyFlux.setMassRate(bdyInfo.massRate, bdyInfo.pvtRegionIdx);
618 }
619
620 template <class BoundaryConditionData>
621 static void computeBoundaryFluxFree(const Problem& problem,
622 RateVector& bdyFlux,
623 const BoundaryConditionData& bdyInfo,
624 const IntensiveQuantities& insideIntQuants,
625 unsigned globalSpaceIdx)
626 {
627 OPM_TIMEBLOCK_LOCAL(computeBoundaryFluxFree, Subsystem::Assembly);
628 std::array<short, numPhases> upIdx;
629 std::array<short, numPhases> dnIdx;
630 std::array<Evaluation, numPhases> volumeFlux;
631 std::array<Evaluation, numPhases> pressureDifference;
632
633 ExtensiveQuantities::calculateBoundaryGradients_(problem,
634 globalSpaceIdx,
635 insideIntQuants,
636 bdyInfo.boundaryFaceIndex,
637 bdyInfo.faceArea,
638 bdyInfo.faceZCoord,
639 bdyInfo.exFluidState,
640 upIdx,
641 dnIdx,
642 volumeFlux,
643 pressureDifference);
644
646 // advective fluxes of all components in all phases
648 bdyFlux = 0.0;
649 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
650 if (!FluidSystem::phaseIsActive(phaseIdx)) {
651 continue;
652 }
653 const auto& pBoundary = bdyInfo.exFluidState.pressure(phaseIdx);
654 const Evaluation& pInside = insideIntQuants.fluidState().pressure(phaseIdx);
655 const unsigned pvtRegionIdx = insideIntQuants.pvtRegionIndex();
656
657 RateVector tmp(0.0);
658 const auto& darcyFlux = volumeFlux[phaseIdx];
659 // mass conservation
660 if (pBoundary < pInside) {
661 // outflux
662 const auto& invB = getInvB_<FluidSystem, FluidState, Evaluation>(
663 insideIntQuants.fluidState(), phaseIdx, pvtRegionIdx);
664 Evaluation surfaceVolumeFlux = invB * darcyFlux;
665 evalPhaseFluxes_<Evaluation>(tmp,
666 phaseIdx,
667 insideIntQuants.pvtRegionIndex(),
668 surfaceVolumeFlux,
669 insideIntQuants.fluidState());
670 if constexpr (enableFullyImplicitThermal) {
671 EnergyModule::template addPhaseEnthalpyFluxes_<Evaluation>(
672 tmp, phaseIdx, darcyFlux, insideIntQuants.fluidState());
673 }
674 } else if (pBoundary > pInside) {
675 // influx
676 using ScalarFluidState = decltype(bdyInfo.exFluidState);
677 const auto& invB = getInvB_<FluidSystem, ScalarFluidState, Scalar>(
678 bdyInfo.exFluidState, phaseIdx, pvtRegionIdx);
679 Evaluation surfaceVolumeFlux = invB * darcyFlux;
680 evalPhaseFluxes_<Scalar>(tmp,
681 phaseIdx,
682 insideIntQuants.pvtRegionIndex(),
683 surfaceVolumeFlux,
684 bdyInfo.exFluidState);
685 if constexpr (enableFullyImplicitThermal) {
686 EnergyModule::template addPhaseEnthalpyFluxes_<Scalar>(
687 tmp, phaseIdx, darcyFlux, bdyInfo.exFluidState);
688 }
689 }
690
691 for (unsigned i = 0; i < tmp.size(); ++i) {
692 bdyFlux[i] += tmp[i];
693 }
694 }
695
696 // conductive heat flux from boundary
697 if constexpr (enableFullyImplicitThermal) {
698 Evaluation heatFlux;
699 // avoid overload of functions with same number of elements in eclproblem
700 Scalar alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
701 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
702 unsigned inIdx = 0; // dummy
703 // always calculated with derivatives of this cell
704 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
705 insideIntQuants,
706 /*focusDofIndex*/ inIdx,
707 inIdx,
708 alpha,
709 bdyInfo.exFluidState);
710 EnergyModule::addHeatFlux(bdyFlux, heatFlux);
711 }
712
713 static_assert(
714 !enableSolvent,
715 "Relevant treatment of boundary conditions must be implemented before enabling.");
716 static_assert(
717 !enablePolymer,
718 "Relevant treatment of boundary conditions must be implemented before enabling.");
719
720 const FluidSystem& fsys = insideIntQuants.getFluidSystem();
721
722 // make sure that the right mass conservation quantities are used
723 adaptMassConservationQuantities_(bdyFlux, insideIntQuants.pvtRegionIndex(), fsys);
724
725#ifndef NDEBUG
726 for (unsigned i = 0; i < numEq; ++i) {
727 Valgrind::CheckDefined(bdyFlux[i]);
728 }
729 Valgrind::CheckDefined(bdyFlux);
730#endif
731 }
732
733 template <class ProblemLocal, class BoundaryConditionData, class RateVectorLocal>
734 OPM_HOST_DEVICE static void computeBoundaryThermal(const ProblemLocal& problem,
735 RateVectorLocal& bdyFlux,
736 const BoundaryConditionData& bdyInfo,
737 const IntensiveQuantities& insideIntQuants,
738 [[maybe_unused]] unsigned globalSpaceIdx)
739 {
740 OPM_TIMEBLOCK_LOCAL(computeBoundaryThermal, Subsystem::Assembly);
741 // only heat is allowed to flow through this boundary
742 bdyFlux = 0.0;
743
744 // conductive heat flux from boundary
745 if constexpr (enableFullyImplicitThermal) {
746 Evaluation heatFlux;
747 // avoid overload of functions with same numeber of elements in eclproblem
748
749 Scalar alpha;
750 if constexpr (runAssemblyOnGpu) {
751 // This path is currently only intended for the SimplifiedBlackoilModel for GPUs
752 // which currently does not aim to reproduce the full problem object on the GPU.
753 alpha = problem.getAlpha(globalSpaceIdx, bdyInfo.boundaryFaceIndex);
754 } else {
755 alpha = problem.eclTransmissibilities().thermalHalfTransBoundary(
756 globalSpaceIdx, bdyInfo.boundaryFaceIndex);
757 }
758
759 unsigned inIdx = 0; // dummy
760 // always calculated with derivatives of this cell
761 EnergyModule::ExtensiveQuantities::updateEnergyBoundary(heatFlux,
762 insideIntQuants,
763 /*focusDofIndex*/ inIdx,
764 inIdx,
765 alpha,
766 bdyInfo.exFluidState);
767 EnergyModule::addHeatFlux(bdyFlux, heatFlux);
768 }
769
770#ifndef NDEBUG
771 for (unsigned i = 0; i < numEq; ++i) {
772 Valgrind::CheckDefined(bdyFlux[i]);
773 }
774 Valgrind::CheckDefined(bdyFlux);
775#endif
776 }
777
778 static void computeSource(RateVector& source,
779 const Problem& problem,
780 const IntensiveQuantities& insideIntQuants,
781 unsigned globalSpaceIdex,
782 unsigned timeIdx)
783 {
784 OPM_TIMEBLOCK_LOCAL(computeSource, Subsystem::Assembly);
785 // retrieve the source term intrinsic to the problem
786 problem.source(source, globalSpaceIdex, timeIdx);
787
788 // deal with bioeffects (if present)
789 BioeffectsModule::addSource(source, problem, insideIntQuants, globalSpaceIdex);
790
791 // scale the source term of the energy equation
792 if constexpr (enableFullyImplicitThermal) {
793 source[Indices::contiEnergyEqIdx]
794 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
795 }
796 }
797
798 static void computeSourceDense(RateVector& source,
799 const Problem& problem,
800 const IntensiveQuantities& insideIntQuants,
801 unsigned globalSpaceIdex,
802 unsigned timeIdx)
803 {
804 source = 0.0;
805 problem.addToSourceDense(source, globalSpaceIdex, timeIdx);
806
807 // deal with bioeffects (if present)
808 BioeffectsModule::addSource(source, problem, insideIntQuants, globalSpaceIdex);
809
810 // scale the source term of the energy equation
811 if constexpr (enableFullyImplicitThermal) {
812 source[Indices::contiEnergyEqIdx]
813 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
814 }
815 }
816
820 void computeSource(RateVector& source,
821 const ElementContext& elemCtx,
822 unsigned dofIdx,
823 unsigned timeIdx) const
824 {
825 OPM_TIMEBLOCK_LOCAL(computeSource, Subsystem::Assembly);
826 // retrieve the source term intrinsic to the problem
827 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
828
829 // deal with bioeffects (if present)
830 BioeffectsModule::addSource(source, elemCtx, dofIdx, timeIdx);
831
832 // scale the source term of the energy equation
833 if constexpr (enableFullyImplicitThermal) {
834 source[Indices::contiEnergyEqIdx]
835 *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
836 }
837 }
838
839 template <class UpEval, class FluidState>
840 static void evalPhaseFluxes_(RateVector& flux,
841 unsigned phaseIdx,
842 unsigned pvtRegionIdx,
843 const ExtensiveQuantities& extQuants,
844 const FluidState& upFs)
845 {
846 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
847 const auto& surfaceVolumeFlux = invB * extQuants.volumeFlux(phaseIdx);
848 evalPhaseFluxes_<UpEval>(flux, phaseIdx, pvtRegionIdx, surfaceVolumeFlux, upFs);
849 }
850
855 template <class UpEval, class Eval, class FluidState, class RateVectorT = RateVector>
856 OPM_HOST_DEVICE static void evalPhaseFluxes_(RateVectorT& flux,
857 unsigned phaseIdx,
858 unsigned pvtRegionIdx,
859 const Eval& surfaceVolumeFlux,
860 const FluidState& upFs)
861 {
862 const FluidSystem& fsys = upFs.fluidSystem();
863
864 unsigned activeCompIdx
865 = fsys.canonicalToActiveCompIdx(fsys.solventComponentIndex(phaseIdx));
866
867 if constexpr (blackoilConserveSurfaceVolume) {
868 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
869 } else {
870 flux[conti0EqIdx + activeCompIdx]
871 += surfaceVolumeFlux * fsys.referenceDensity(phaseIdx, pvtRegionIdx);
872 }
873
874 if (phaseIdx == oilPhaseIdx) {
875 // dissolved gas (in the oil phase).
876 if (fsys.enableDissolvedGas()) {
877 const auto& Rs
878 = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
879
880 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
881 if constexpr (blackoilConserveSurfaceVolume) {
882 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
883 } else {
884 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux
885 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
886 }
887 }
888 } else if (phaseIdx == waterPhaseIdx) {
889 // dissolved gas (in the water phase).
890 if (fsys.enableDissolvedGasInWater()) {
891 const auto& Rsw
892 = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
893
894 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
895 if constexpr (blackoilConserveSurfaceVolume) {
896 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
897 } else {
898 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux
899 * fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
900 }
901 }
902 } else if (phaseIdx == gasPhaseIdx) {
903 // vaporized oil (in the gas phase).
904 if (fsys.enableVaporizedOil()) {
905 const auto& Rv
906 = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
907
908 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
909 if constexpr (blackoilConserveSurfaceVolume) {
910 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
911 } else {
912 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux
913 * fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
914 }
915 }
916 // vaporized water (in the gas phase).
917 if (fsys.enableVaporizedWater()) {
918 const auto& Rvw
919 = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
920
921 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
922 if constexpr (blackoilConserveSurfaceVolume) {
923 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
924 } else {
925 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux
926 * fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
927 }
928 }
929 }
930 }
931
940 template <class Scalar>
941 static void adaptMassConservationQuantities_(Dune::FieldVector<Scalar, numEq>& container,
942 unsigned pvtRegionIdx)
943 {
944 adaptMassConservationQuantities_(container, pvtRegionIdx, FluidSystem {});
945 }
946
961 template <class ScalarVector, class FsysType>
962 OPM_HOST_DEVICE static void adaptMassConservationQuantities_(ScalarVector& container,
963 unsigned pvtRegionIdx,
964 const FsysType& fsys)
965 {
966 if constexpr (!blackoilConserveSurfaceVolume) {
967 // convert "surface volume" to mass. this is complicated a bit by the fact that
968 // not all phases are necessarily enabled. (we here assume that if a fluid phase
969 // is disabled, its respective "main" component is not considered as well.)
970
971 if constexpr (waterEnabled) {
972 const unsigned activeWaterCompIdx = fsys.canonicalToActiveCompIdx(waterCompIdx);
973 container[conti0EqIdx + activeWaterCompIdx]
974 *= fsys.referenceDensity(waterPhaseIdx, pvtRegionIdx);
975 }
976
977 if constexpr (gasEnabled) {
978 const unsigned activeGasCompIdx = fsys.canonicalToActiveCompIdx(gasCompIdx);
979 container[conti0EqIdx + activeGasCompIdx]
980 *= fsys.referenceDensity(gasPhaseIdx, pvtRegionIdx);
981 }
982
983 if constexpr (oilEnabled) {
984 const unsigned activeOilCompIdx = fsys.canonicalToActiveCompIdx(oilCompIdx);
985 container[conti0EqIdx + activeOilCompIdx]
986 *= fsys.referenceDensity(oilPhaseIdx, pvtRegionIdx);
987 }
988 }
989 }
990
991 // NNC does not have a direction
992 static FaceDir::DirEnum faceDirFromDirId(const int dirId)
993 {
994 return dirId < 0 ? FaceDir::DirEnum::Unknown : FaceDir::FromIntersectionIndex(dirId);
995 }
996};
997
998} // namespace Opm
999
1000#endif
Contains the classes required to extend the black-oil model by bioeffects.
Contains the classes required to extend the black-oil model by brine.
Classes required for dynamic convective mixing.
Classes required for molecular diffusion.
Classes required for mechanical dispersion.
Contains the classes required to extend the black-oil model by energy.
Contains the classes required to extend the black-oil model by solvent component. For details,...
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by polymer.
Declares the properties required by the black oil model.
Contains the classes required to extend the black-oil model by solvents.
Contains the high level supplements required to extend the black oil model by bioeffects.
Definition: blackoilbioeffectsmodules.hh:95
static void applyScaling(RateVector &flux)
Definition: blackoilbioeffectsmodules.hh:238
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbioeffectsmodules.hh:179
static void addSource(RateVector &source, const Problem &problem, const IntensiveQuantities &intQuants, unsigned globalSpaceIdex)
Definition: blackoilbioeffectsmodules.hh:279
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:58
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:165
Definition: blackoilconvectivemixingmodule.hh:99
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: blackoildiffusionmodule.hh:54
Provides the auxiliary methods required for consideration of the dispersion equation.
Definition: blackoildispersionmodule.hh:56
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 OPM_HOST_DEVICE void addHeatFlux(RateVectorT &flux, const Evaluation &heatFlux)
Definition: blackoilenergymodules.hh:221
Contains the high level supplements required to extend the black oil model.
Definition: blackoilextbomodules.hh:64
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilextbomodules.hh:158
Contains the high level supplements required to extend the black oil model to include the effects of ...
Definition: blackoilfoammodules.hh:59
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilfoammodules.hh:165
Calculates the local residual of the black oil model.
Definition: blackoillocalresidualtpfa.hh:71
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidualtpfa.hh:820
static void computeBoundaryFluxFree(const Problem &problem, RateVector &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:621
static OPM_HOST_DEVICE void computeStorage(StorageType &storage, const IntensiveQuantitiesType &intQuants)
Definition: blackoillocalresidualtpfa.hh:169
static void computeBoundaryFluxRate(RateVector &bdyFlux, const BoundaryConditionData &bdyInfo)
Definition: blackoillocalresidualtpfa.hh:615
static FaceDir::DirEnum faceDirFromDirId(const int dirId)
Definition: blackoillocalresidualtpfa.hh:992
static OPM_HOST_DEVICE void computeBoundaryThermal(const ProblemLocal &problem, RateVectorLocal &bdyFlux, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:734
static OPM_HOST_DEVICE void calculateFluxes_(RateVectorT &flux, RateVectorT &darcy, const IntensiveQuantitiesT &intQuantsIn, const IntensiveQuantitiesT &intQuantsEx, const unsigned &globalIndexIn, const unsigned &globalIndexEx, const ResidualNBInfoT &nbInfo, const ModuleParamsT &moduleParams)
Definition: blackoillocalresidualtpfa.hh:360
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Definition: blackoillocalresidualtpfa.hh:840
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:284
static OPM_HOST_DEVICE void evalPhaseFluxes_(RateVectorT &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const Eval &surfaceVolumeFlux, const FluidState &upFs)
Helper function to calculate the flux of mass in terms of conservation quantities via specific fluid ...
Definition: blackoillocalresidualtpfa.hh:856
static void adaptMassConservationQuantities_(Dune::FieldVector< Scalar, numEq > &container, unsigned pvtRegionIdx)
Helper function to convert the mass-related parts of a Dune::FieldVector that stores conservation qua...
Definition: blackoillocalresidualtpfa.hh:941
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: blackoillocalresidualtpfa.hh:159
static OPM_HOST_DEVICE void computeFlux(RateVectorT &flux, RateVectorT &darcy, const unsigned globalIndexIn, const unsigned globalIndexEx, const IntensiveQuantitiesT &intQuantsIn, const IntensiveQuantitiesT &intQuantsEx, const ResidualNBInfoT &nbInfo, const ModuleParamsT &moduleParams)
Definition: blackoillocalresidualtpfa.hh:257
static OPM_HOST_DEVICE void computeBoundaryFlux(RateVectorLocal &bdyFlux, const LocalProblem &problem, const BoundaryConditionData &bdyInfo, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdx)
Definition: blackoillocalresidualtpfa.hh:572
static void computeSource(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:778
static void computeSourceDense(RateVector &source, const Problem &problem, const IntensiveQuantities &insideIntQuants, unsigned globalSpaceIdex, unsigned timeIdx)
Definition: blackoillocalresidualtpfa.hh:798
static OPM_HOST_DEVICE void adaptMassConservationQuantities_(ScalarVector &container, unsigned pvtRegionIdx, const FsysType &fsys)
Helper function to convert the mass-related parts of a vector that stores conservation quantities in ...
Definition: blackoillocalresidualtpfa.hh:962
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:65
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilpolymermodules.hh:237
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:185
@ NONE
Definition: DeferredLogger.hpp:46
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: blackoillocalresidualtpfa.hh:141
FaceDir::DirEnum faceDir
Definition: blackoillocalresidualtpfa.hh:146
double faceArea
Definition: blackoillocalresidualtpfa.hh:143
ConditionalStorage< enableFullyImplicitThermal, double > outAlpha
Definition: blackoillocalresidualtpfa.hh:150
ConditionalStorage< enableFullyImplicitThermal, double > inAlpha
Definition: blackoillocalresidualtpfa.hh:149
ConditionalStorage< enableDiffusion, double > diffusivity
Definition: blackoillocalresidualtpfa.hh:151
double dZg
Definition: blackoillocalresidualtpfa.hh:145
double Vin
Definition: blackoillocalresidualtpfa.hh:147
ConditionalStorage< enableDispersion, double > dispersivity
Definition: blackoillocalresidualtpfa.hh:152
double thpres
Definition: blackoillocalresidualtpfa.hh:144
double trans
Definition: blackoillocalresidualtpfa.hh:142
double Vex
Definition: blackoillocalresidualtpfa.hh:148