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