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