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 "blackoilproperties.hh"
40#include <opm/material/fluidstates/BlackOilFluidState.hpp>
41
42namespace Opm {
48template <class TypeTag>
49class BlackOilLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
50{
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 const bool waterEnabled = Indices::waterEnabled;
76 static const bool gasEnabled = Indices::gasEnabled;
77 static const bool oilEnabled = Indices::oilEnabled;
78 static const bool compositionSwitchEnabled = (compositionSwitchIdx >= 0);
79
80 static constexpr bool blackoilConserveSurfaceVolume = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
81 static constexpr bool enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>();
82 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
83
84 using Toolbox = MathToolbox<Evaluation>;
93
94public:
98 template <class LhsEval>
99 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
100 const ElementContext& elemCtx,
101 unsigned dofIdx,
102 unsigned timeIdx) const
103 {
104 // retrieve the intensive quantities for the SCV at the specified point in time
105 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
106 const auto& fs = intQuants.fluidState();
107
108 storage = 0.0;
109
110 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
111 if (!FluidSystem::phaseIsActive(phaseIdx)) {
112 if (Indices::numPhases == 3) { // add trivial equation for the pseudo phase
113 unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
114 if (timeIdx == 0)
115 storage[conti0EqIdx + activeCompIdx] = variable<LhsEval>(0.0, conti0EqIdx + activeCompIdx);
116 else
117 storage[conti0EqIdx + activeCompIdx] = 0.0;
118 }
119 continue;
120 }
121
122 unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
123 LhsEval surfaceVolume =
124 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
125 * Toolbox::template decay<LhsEval>(fs.invB(phaseIdx))
126 * Toolbox::template decay<LhsEval>(intQuants.porosity());
127
128 storage[conti0EqIdx + activeCompIdx] += surfaceVolume;
129
130 // account for dissolved gas
131 if (phaseIdx == oilPhaseIdx && FluidSystem::enableDissolvedGas()) {
132 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(gasCompIdx);
133 storage[conti0EqIdx + activeGasCompIdx] +=
134 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rs())
135 * surfaceVolume;
136 }
137
138 // account for dissolved gas in water phase
139 if (phaseIdx == waterPhaseIdx && FluidSystem::enableDissolvedGasInWater()) {
140 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(gasCompIdx);
141 storage[conti0EqIdx + activeGasCompIdx] +=
142 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rsw())
143 * surfaceVolume;
144 }
145
146 // account for vaporized oil
147 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedOil()) {
148 unsigned activeOilCompIdx = Indices::canonicalToActiveComponentIndex(oilCompIdx);
149 storage[conti0EqIdx + activeOilCompIdx] +=
150 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rv())
151 * surfaceVolume;
152 }
153
154 // account for vaporized water
155 if (phaseIdx == gasPhaseIdx && FluidSystem::enableVaporizedWater()) {
156 unsigned activeWaterCompIdx = Indices::canonicalToActiveComponentIndex(waterCompIdx);
157 storage[conti0EqIdx + activeWaterCompIdx] +=
158 Toolbox::template decay<LhsEval>(intQuants.fluidState().Rvw())
159 * surfaceVolume;
160 }
161 }
162
163 adaptMassConservationQuantities_(storage, intQuants.pvtRegionIndex());
164
165 // deal with solvents (if present)
166 SolventModule::addStorage(storage, intQuants);
167
168 // deal with zFracton (if present)
169 ExtboModule::addStorage(storage, intQuants);
170
171 // deal with polymer (if present)
172 PolymerModule::addStorage(storage, intQuants);
173
174 // deal with energy (if present)
175 EnergyModule::addStorage(storage, intQuants);
176
177 // deal with foam (if present)
178 FoamModule::addStorage(storage, intQuants);
179
180 // deal with salt (if present)
181 BrineModule::addStorage(storage, intQuants);
182
183 // deal with micp (if present)
184 MICPModule::addStorage(storage, intQuants);
185 }
186
190 void computeFlux(RateVector& flux,
191 const ElementContext& elemCtx,
192 unsigned scvfIdx,
193 unsigned timeIdx) const
194 {
195 assert(timeIdx == 0);
196
197 flux = 0.0;
198
199 const ExtensiveQuantities& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
200 unsigned focusDofIdx = elemCtx.focusDofIndex();
201 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
202 if (!FluidSystem::phaseIsActive(phaseIdx))
203 continue;
204
205 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
206 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
207 unsigned pvtRegionIdx = up.pvtRegionIndex();
208 if (upIdx == focusDofIdx)
209 evalPhaseFluxes_<Evaluation>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
210 else
211 evalPhaseFluxes_<Scalar>(flux, phaseIdx, pvtRegionIdx, extQuants, up.fluidState());
212 }
213
214 // deal with solvents (if present)
215 SolventModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
216
217 // deal with zFracton (if present)
218 ExtboModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
219
220 // deal with polymer (if present)
221 PolymerModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
222
223 // deal with energy (if present)
224 EnergyModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
225
226 // deal with foam (if present)
227 FoamModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
228
229 // deal with salt (if present)
230 BrineModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
231
232 // deal with micp (if present)
233 MICPModule::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
234
235 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
236 }
237
241 void computeSource(RateVector& source,
242 const ElementContext& elemCtx,
243 unsigned dofIdx,
244 unsigned timeIdx) const
245 {
246 // retrieve the source term intrinsic to the problem
247 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
248
249 // deal with MICP (if present)
250 MICPModule::addSource(source, elemCtx, dofIdx, timeIdx);
251
252 // scale the source term of the energy equation
253 if (enableEnergy)
254 source[Indices::contiEnergyEqIdx] *= getPropValue<TypeTag, Properties::BlackOilEnergyScalingFactor>();
255 }
256
261 template <class UpEval, class FluidState>
262 static void evalPhaseFluxes_(RateVector& flux,
263 unsigned phaseIdx,
264 unsigned pvtRegionIdx,
265 const ExtensiveQuantities& extQuants,
266 const FluidState& upFs)
267 {
268 const auto& invB = getInvB_<FluidSystem, FluidState, UpEval>(upFs, phaseIdx, pvtRegionIdx);
269 const auto& surfaceVolumeFlux = invB*extQuants.volumeFlux(phaseIdx);
270 unsigned activeCompIdx = Indices::canonicalToActiveComponentIndex(FluidSystem::solventComponentIndex(phaseIdx));
271
272 if (blackoilConserveSurfaceVolume)
273 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux;
274 else
275 flux[conti0EqIdx + activeCompIdx] += surfaceVolumeFlux*FluidSystem::referenceDensity(phaseIdx, pvtRegionIdx);
276
277 if (phaseIdx == oilPhaseIdx) {
278 // dissolved gas (in the oil phase).
279 if (FluidSystem::enableDissolvedGas()) {
280 const auto& Rs = BlackOil::getRs_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
281
282 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(gasCompIdx);
283 if (blackoilConserveSurfaceVolume)
284 flux[conti0EqIdx + activeGasCompIdx] += Rs*surfaceVolumeFlux;
285 else
286 flux[conti0EqIdx + activeGasCompIdx] += Rs*surfaceVolumeFlux*FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
287 }
288 } else if (phaseIdx == waterPhaseIdx) {
289 // dissolved gas (in the water phase).
290 if (FluidSystem::enableDissolvedGasInWater()) {
291 const auto& Rsw = BlackOil::getRsw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
292
293 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(gasCompIdx);
294 if (blackoilConserveSurfaceVolume)
295 flux[conti0EqIdx + activeGasCompIdx] += Rsw*surfaceVolumeFlux;
296 else
297 flux[conti0EqIdx + activeGasCompIdx] += Rsw*surfaceVolumeFlux*FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
298 }
299 }
300 else if (phaseIdx == gasPhaseIdx) {
301 // vaporized oil (in the gas phase).
302 if (FluidSystem::enableVaporizedOil()) {
303 const auto& Rv = BlackOil::getRv_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
304
305 unsigned activeOilCompIdx = Indices::canonicalToActiveComponentIndex(oilCompIdx);
306 if (blackoilConserveSurfaceVolume)
307 flux[conti0EqIdx + activeOilCompIdx] += Rv*surfaceVolumeFlux;
308 else
309 flux[conti0EqIdx + activeOilCompIdx] += Rv*surfaceVolumeFlux*FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
310 }
311 // vaporized water (in the gas phase).
312 if (FluidSystem::enableVaporizedWater()) {
313 const auto& Rvw = BlackOil::getRvw_<FluidSystem, FluidState, UpEval>(upFs, pvtRegionIdx);
314
315 unsigned activeWaterCompIdx = Indices::canonicalToActiveComponentIndex(waterCompIdx);
316 if (blackoilConserveSurfaceVolume)
317 flux[conti0EqIdx + activeWaterCompIdx] += Rvw*surfaceVolumeFlux;
318 else
319 flux[conti0EqIdx + activeWaterCompIdx] += Rvw*surfaceVolumeFlux*FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
320 }
321 }
322 }
323
335 template <class Scalar>
336 static void adaptMassConservationQuantities_(Dune::FieldVector<Scalar, numEq>& container, unsigned pvtRegionIdx)
337 {
338 if (blackoilConserveSurfaceVolume)
339 return;
340
341 // convert "surface volume" to mass. this is complicated a bit by the fact that
342 // not all phases are necessarily enabled. (we here assume that if a fluid phase
343 // is disabled, its respective "main" component is not considered as well.)
344
345 if (waterEnabled) {
346 unsigned activeWaterCompIdx = Indices::canonicalToActiveComponentIndex(waterCompIdx);
347 container[conti0EqIdx + activeWaterCompIdx] *=
348 FluidSystem::referenceDensity(waterPhaseIdx, pvtRegionIdx);
349 }
350
351 if (gasEnabled) {
352 unsigned activeGasCompIdx = Indices::canonicalToActiveComponentIndex(gasCompIdx);
353 container[conti0EqIdx + activeGasCompIdx] *=
354 FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx);
355 }
356
357 if (oilEnabled) {
358 unsigned activeOilCompIdx = Indices::canonicalToActiveComponentIndex(oilCompIdx);
359 container[conti0EqIdx + activeOilCompIdx] *=
360 FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx);
361 }
362 }
363};
364
365} // namespace Opm
366
367#endif
Contains the classes required to extend the black-oil model by brine.
Classes required for molecular diffusion.
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 computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:260
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
Contains the high level supplements required to extend the black oil model by energy.
Definition: blackoilenergymodules.hh:52
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilenergymodules.hh:168
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilenergymodules.hh:141
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
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilextbomodules.hh:354
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
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilfoammodules.hh:286
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:50
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:190
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidual.hh:241
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:262
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:336
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:99
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:56
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilmicpmodules.hh:265
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
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:512
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
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:462
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