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