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