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