blackoillocalresidual.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_RESIDUAL_HH
29#define EWOMS_BLACK_OIL_LOCAL_RESIDUAL_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/fluidstates/BlackOilFluidState.hpp>
33
36
38
39#include <cassert>
40#include <limits>
41
42namespace Opm {
43
49template <class TypeTag>
50class BlackOilLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
51{
60
61 enum { conti0EqIdx = Indices::conti0EqIdx };
62 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
63 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
64 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
65
66 enum { gasPhaseIdx = FluidSystem::gasPhaseIdx };
67 enum { oilPhaseIdx = FluidSystem::oilPhaseIdx };
68 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
69
70 enum { gasCompIdx = FluidSystem::gasCompIdx };
71 enum { oilCompIdx = FluidSystem::oilCompIdx };
72 enum { waterCompIdx = FluidSystem::waterCompIdx };
73 static constexpr unsigned compositionSwitchIdx = Indices::compositionSwitchIdx;
74
75 static constexpr bool waterEnabled = Indices::waterEnabled;
76 static constexpr bool gasEnabled = Indices::gasEnabled;
77 static constexpr bool oilEnabled = Indices::oilEnabled;
78 static constexpr bool compositionSwitchEnabled =
79 compositionSwitchIdx != std::numeric_limits<unsigned>::max();
80
81 static constexpr bool blackoilConserveSurfaceVolume =
82 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
83 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
84 static constexpr bool enableBrine = getPropValue<TypeTag, Properties::EnableBrine>();
85 static constexpr bool enableConvectiveMixing = getPropValue<TypeTag, Properties::EnableConvectiveMixing>();
86 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
87 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
88 static constexpr bool enableFoam = getPropValue<TypeTag, Properties::EnableFoam>();
89 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
90 static constexpr bool enableFullyImplicitThermal = energyModuleType == EnergyModules::FullyImplicitThermal;
91 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
92 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
93
94 using Toolbox = MathToolbox<Evaluation>;
95
96 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
97 using BrineModule = BlackOilBrineModule<TypeTag, enableBrine>;
98 using ConvectiveMixingModule = BlackOilConvectiveMixingModule<TypeTag, enableConvectiveMixing>;
101 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
103 using PolymerModule = BlackOilPolymerModule<TypeTag, enablePolymer>;
104 using SolventModule = BlackOilSolventModule<TypeTag, enableSolvent>;
105
106public:
110 template <class LhsEval>
111 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
112 const ElementContext& elemCtx,
113 unsigned dofIdx,
114 unsigned timeIdx) const
115 {
116 // retrieve the intensive quantities for the SCV at the specified point in time
117 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
118 const auto& fs = intQuants.fluidState();
119
120 storage = 0.0;
121
122 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
123 if (!FluidSystem::phaseIsActive(phaseIdx)) {
124 if (Indices::numPhases == 3) { // add trivial equation for the pseudo phase
125 const unsigned activeCompIdx =
126 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
127 if (timeIdx == 0) {
128 storage[conti0EqIdx + activeCompIdx] = variable<LhsEval>(0.0, conti0EqIdx + activeCompIdx);
129 }
130 else {
131 storage[conti0EqIdx + activeCompIdx] = 0.0;
132 }
133 }
134 continue;
135 }
136
137 const unsigned activeCompIdx =
138 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
139 const LhsEval surfaceVolume =
140 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
141 Toolbox::template decay<LhsEval>(fs.invB(phaseIdx)) *
142 Toolbox::template decay<LhsEval>(intQuants.porosity());
143
144 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
145
146 // account for dissolved gas
147 if (phaseIdx == oilPhaseIdx && FluidSystem::enableDissolvedGas()) {
148 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
149 storage[conti0EqIdx + activeGasCompIdx] +=
150 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs()) *
151 surfaceVolume;
152 }
153
154 // account for dissolved gas in water phase
155 if (phaseIdx == waterPhaseIdx && FluidSystem::enableDissolvedGasInWater()) {
156 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
157 storage[conti0EqIdx + activeGasCompIdx] +=
158 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw()) *
159 surfaceVolume;
160 }
161
162 // account for vaporized oil
163 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedOil()) {
164 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
165 storage[conti0EqIdx + activeOilCompIdx] +=
166 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv()) *
167 surfaceVolume;
168 }
169
170 // account for vaporized water
171 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedWater()) {
172 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
173 storage[conti0EqIdx + activeWaterCompIdx] +=
174 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw()) *
175 surfaceVolume;
176 }
177 }
178
179 adaptMassConservationQuantities_(storage, intQuants.pvtRegionIndex());
180
181 // deal with solvents (if present)
182 if constexpr (enableSolvent) {
183 SolventModule::addStorage(storage, intQuants);
184 }
185
186 // deal with zFracton (if present)
187 if constexpr (enableExtbo) {
188 ExtboModule::addStorage(storage, intQuants);
189 }
190
191 // deal with polymer (if present)
192 if constexpr (enablePolymer) {
193 PolymerModule::addStorage(storage, intQuants);
194 }
195
196 // deal with energy (if present)
197 if constexpr (enableFullyImplicitThermal) {
198 EnergyModule::addStorage(storage, intQuants);
199 }
200
201 // deal with foam (if present)
202 if constexpr (enableFoam) {
203 FoamModule::addStorage(storage, intQuants);
204 }
205
206 // deal with salt (if present)
207 if constexpr (enableBrine) {
208 BrineModule::addStorage(storage, intQuants);
209 }
210
211 // deal with bioeffects (if present)
212 if constexpr (enableBioeffects) {
213 BioeffectsModule::addStorage(storage, intQuants);
214 }
215 }
216
220 void computeFlux(RateVector& flux,
221 const ElementContext& elemCtx,
222 unsigned scvfIdx,
223 unsigned timeIdx) const
224 {
225 assert(timeIdx == 0);
226
227 flux = 0.0;
228
229 const ExtensiveQuantities& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
230 unsigned focusDofIdx = elemCtx.focusDofIndex();
231 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
232 if (!FluidSystem::phaseIsActive(phaseIdx)) {
233 continue;
234 }
235
236 const unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
237 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
238 const unsigned pvtRegionIdx = up.pvtRegionIndex();
239 if (upIdx == focusDofIdx) {
240 evalPhaseFluxes_<Evaluation>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
241 }
242 else {
243 evalPhaseFluxes_<Scalar>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
244 }
245 }
246
247 // deal with solvents (if present)
248 if constexpr (enableSolvent) {
249 SolventModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
250 }
251
252 // deal with zFracton (if present)
253 if constexpr (enableExtbo) {
254 ExtboModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
255 }
256
257 // deal with polymer (if present)
258 if constexpr (enablePolymer) {
259 PolymerModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
260 }
261
262 // deal with energy (if present)
263 if constexpr (enableFullyImplicitThermal) {
264 EnergyModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
265 }
266
267 // deal with foam (if present)
268 if constexpr (enableFoam) {
269 FoamModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
270 }
271
272 // deal with salt (if present)
273 if constexpr (enableBrine) {
274 BrineModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
275 }
276
277 // deal with bioeffects (if present)
278 if constexpr (enableBioeffects) {
279 BioeffectsModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
280 }
281
282 // deal with diffusion (if present)
283 if constexpr (enableDiffusion) {
284 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
285 }
286
287 // deal with convective mixing (if enabled)
288 if constexpr (enableConvectiveMixing) {
289 ConvectiveMixingModule::addConvectiveMixingFlux(flux, elemCtx, scvfIdx, timeIdx);
290 }
291 }
292
296 void computeSource(RateVector& source,
297 const ElementContext& elemCtx,
298 unsigned dofIdx,
299 unsigned timeIdx) const
300 {
301 // retrieve the source term intrinsic to the problem
302 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
303
304 // deal with MICP (if present)
305 if constexpr (enableBioeffects) {
306 BioeffectsModule::addSource(source, elemCtx, dofIdx, timeIdx);
307 }
308
309 // scale the source term of the energy equation
310 if constexpr (enableFullyImplicitThermal) {
311 source[Indices::contiEnergyEqIdx] *=
312 getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
313 }
314 }
315
320 template <class UpEval, class FluidState>
321 static void evalPhaseFluxes_(RateVector& flux,
322 unsigned phaseIdx,
323 unsigned pvtRegionIdx,
324 const ExtensiveQuantities& extQuants,
325 const FluidState& upFs)
326 {
327 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
328 const auto& surfaceVolumeFlux = invB*extQuants.volumeFlux(phaseIdx);
329 const unsigned activeCompIdx =
330 FluidSystem::canonicalToActiveCompIdx(FluidSystem::solventComponentIndex(phaseIdx));
331
332 if constexpr (blackoilConserveSurfaceVolume) {
333 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
334 }
335 else {
336 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux *
337 FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
338 }
339
340 if (phaseIdx == oilPhaseIdx) {
341 // dissolved gas (in the oil phase).
342 if (FluidSystem::enableDissolvedGas()) {
343 const auto& Rs = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
344
345 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
346 if constexpr (blackoilConserveSurfaceVolume) {
347 flux[conti0EqIdx + activeGasCompIdx] += Rs * surfaceVolumeFlux;
348 }
349 else {
350 flux[conti0EqIdx + activeGasCompIdx] +=
351 Rs * surfaceVolumeFlux *
352 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
353 }
354 }
355 } else if (phaseIdx == waterPhaseIdx) {
356 // dissolved gas (in the water phase).
357 if (FluidSystem::enableDissolvedGasInWater()) {
358 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
359
360 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
361 if constexpr (blackoilConserveSurfaceVolume) {
362 flux[conti0EqIdx + activeGasCompIdx] += Rsw * surfaceVolumeFlux;
363 }
364 else {
365 flux[conti0EqIdx + activeGasCompIdx] +=
366 Rsw * surfaceVolumeFlux *
367 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
368 }
369 }
370 }
371 else if (phaseIdx == gasPhaseIdx) {
372 // vaporized oil (in the gas phase).
373 if (FluidSystem::enableVaporizedOil()) {
374 const auto& Rv = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
375
376 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
377 if constexpr (blackoilConserveSurfaceVolume) {
378 flux[conti0EqIdx + activeOilCompIdx] += Rv * surfaceVolumeFlux;
379 }
380 else {
381 flux[conti0EqIdx + activeOilCompIdx] +=
382 Rv * surfaceVolumeFlux *
383 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
384 }
385 }
386 // vaporized water (in the gas phase).
387 if (FluidSystem::enableVaporizedWater()) {
388 const auto& Rvw = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
389
390 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
391 if constexpr (blackoilConserveSurfaceVolume) {
392 flux[conti0EqIdx + activeWaterCompIdx] += Rvw * surfaceVolumeFlux;
393 }
394 else {
395 flux[conti0EqIdx + activeWaterCompIdx] +=
396 Rvw * surfaceVolumeFlux *
397 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
398 }
399 }
400 }
401 }
402
414 template <class Scalar>
415 static void adaptMassConservationQuantities_(Dune::FieldVector<Scalar, numEq>& container,
416 unsigned pvtRegionIdx)
417 {
418 if constexpr (!blackoilConserveSurfaceVolume) {
419 // convert "surface volume" to mass. this is complicated a bit by the fact that
420 // not all phases are necessarily enabled. (we here assume that if a fluid phase
421 // is disabled, its respective "main" component is not considered as well.)
422
423 if constexpr (waterEnabled) {
424 const unsigned activeWaterCompIdx = FluidSystem::canonicalToActiveCompIdx(waterCompIdx);
425 container[conti0EqIdx + activeWaterCompIdx] *=
426 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
427 }
428
429 if constexpr (gasEnabled) {
430 const unsigned activeGasCompIdx = FluidSystem::canonicalToActiveCompIdx(gasCompIdx);
431 container[conti0EqIdx + activeGasCompIdx] *=
432 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
433 }
434
435 if constexpr (oilEnabled) {
436 const unsigned activeOilCompIdx = FluidSystem::canonicalToActiveCompIdx(oilCompIdx);
437 container[conti0EqIdx + activeOilCompIdx] *=
438 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
439 }
440 }
441 }
442};
443
444} // namespace Opm
445
446#endif
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: blackoilmodules.hpp:63
Definition: blackoilfoammodules.hh:57
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:51
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: blackoillocalresidual.hh:220
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidual.hh:296
static void evalPhaseFluxes_(RateVector &flux, unsigned phaseIdx, unsigned pvtRegionIdx, const ExtensiveQuantities &extQuants, const FluidState &upFs)
Helper function to calculate the flux of mass in terms of conservation quantities via specific fluid ...
Definition: blackoillocalresidual.hh:321
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: blackoillocalresidual.hh:415
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: blackoillocalresidual.hh:111
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233