blackoilbrinemodules.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_BRINE_MODULE_HH
29#define EWOMS_BLACK_OIL_BRINE_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/common/utility/gpuDecorators.hpp>
34
37
40
42
43#include <cassert>
44#include <istream>
45#include <ostream>
46#include <stdexcept>
47#include <string>
48
49namespace Opm {
50
56template <class TypeTag, bool enableBrineV = getPropValue<TypeTag, Properties::EnableBrine>()>
58{
71
72 using Toolbox = MathToolbox<Evaluation>;
73
74 using TabulatedFunction = typename BlackOilBrineParams<Scalar>::TabulatedFunction;
75
76 static constexpr unsigned saltConcentrationIdx = Indices::saltConcentrationIdx;
77 static constexpr unsigned contiBrineEqIdx = Indices::contiBrineEqIdx;
78 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
79 static constexpr bool gasEnabled = Indices::gasEnabled;
80 static constexpr bool oilEnabled = Indices::oilEnabled;
81 static constexpr bool enableBrine = enableBrineV;
82 static constexpr bool enableSaltPrecipitation =
83 getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
84
85 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
86 static constexpr unsigned numPhases = FluidSystem::numPhases;
87
88public:
91 {
92 params_ = params;
93 }
94
98 static void registerParameters()
99 {}
100
101 static bool primaryVarApplies(unsigned pvIdx)
102 {
103 if constexpr (enableBrine) {
104 return pvIdx == saltConcentrationIdx;
105 }
106 else {
107 return false;
108 }
109 }
110
114 template <class FluidState>
115 static void assignPrimaryVars(PrimaryVariables& priVars,
116 const FluidState& fluidState)
117 {
118 if constexpr (enableBrine) {
119 priVars[saltConcentrationIdx] = fluidState.saltConcentration();
120 }
121 }
122
123 static std::string primaryVarName([[maybe_unused]] unsigned pvIdx)
124 {
125 assert(primaryVarApplies(pvIdx));
126
127 return "saltConcentration";
128 }
129
130 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
131 {
132 assert(primaryVarApplies(pvIdx));
133
134 // TODO: it may be beneficial to choose this differently.
135 return static_cast<Scalar>(1.0);
136 }
137
138 static bool eqApplies(unsigned eqIdx)
139 {
140 if constexpr (enableBrine) {
141 return eqIdx == contiBrineEqIdx;
142 }
143 else {
144 return false;
145 }
146 }
147
148 static std::string eqName([[maybe_unused]] unsigned eqIdx)
149 {
150 assert(eqApplies(eqIdx));
151
152 return "conti^brine";
153 }
154
155 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
156 {
157 assert(eqApplies(eqIdx));
158
159 // TODO: it may be beneficial to choose this differently.
160 return static_cast<Scalar>(1.0);
161 }
162
163 // must be called after water storage is computed
164 template <class StorageType>
165 OPM_HOST_DEVICE static void addStorage(StorageType& storage,
166 const IntensiveQuantities& intQuants)
167 {
168 using LhsEval = typename StorageType::value_type;
169
170 if constexpr (enableBrine) {
171 const auto& fs = intQuants.fluidState();
172
173 // avoid singular matrix if no water is present.
174 const LhsEval surfaceVolumeWater =
175 max(Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx)) *
176 Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx)) *
177 Toolbox::template decay<LhsEval>(intQuants.porosity()),
178 1e-10);
179
180 // Brine in water phase
181 const LhsEval massBrine = surfaceVolumeWater *
182 Toolbox::template decay<LhsEval>(fs.saltConcentration());
183
184 if constexpr (enableSaltPrecipitation) {
185 const double saltDensity = intQuants.saltDensity(); // Solid salt density kg/m3
186 const LhsEval solidSalt =
187 Toolbox::template decay<LhsEval>(intQuants.porosity()) /
188 (1.0 - Toolbox::template decay<LhsEval>(fs.saltSaturation()) + 1.e-8) *
189 saltDensity *
190 Toolbox::template decay<LhsEval>(fs.saltSaturation());
191
192 storage[contiBrineEqIdx] += massBrine + solidSalt;
193 }
194 else {
195 storage[contiBrineEqIdx] += massBrine;
196 }
197 }
198 }
199
200 static void computeFlux([[maybe_unused]] RateVector& flux,
201 [[maybe_unused]] const ElementContext& elemCtx,
202 [[maybe_unused]] unsigned scvfIdx,
203 [[maybe_unused]] unsigned timeIdx)
204 {
205 if constexpr (enableBrine) {
206 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
207 unsigned focusIdx = elemCtx.focusDofIndex();
208 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
209 flux[contiBrineEqIdx] = 0.0;
210 if (upIdx == focusIdx)
211 addBrineFluxes_<Evaluation>(flux, elemCtx, scvfIdx, timeIdx);
212 else
213 addBrineFluxes_<Scalar>(flux, elemCtx, scvfIdx, timeIdx);
214 }
215 }
216
217 template <class UpstreamEval>
218 static void addBrineFluxes_(RateVector& flux,
219 const ElementContext& elemCtx,
220 unsigned scvfIdx,
221 unsigned timeIdx)
222 {
223 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
224 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
225 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
226 const auto& upFs = up.fluidState();
227 const auto& volFlux = extQuants.volumeFlux(waterPhaseIdx);
228 addBrineFluxes_<UpstreamEval>(flux, waterPhaseIdx, volFlux, upFs);
229 }
230
231 template <class UpEval, class FluidState>
232 static void addBrineFluxes_(RateVector& flux,
233 unsigned phaseIdx,
234 const Evaluation& volFlux,
235 const FluidState& upFs)
236 {
237 if constexpr (enableBrine) {
238 if (phaseIdx == waterPhaseIdx) {
239 flux[contiBrineEqIdx] =
240 decay<UpEval>(upFs.saltConcentration())
241 * decay<UpEval>(upFs.invB(waterPhaseIdx))
242 * volFlux;
243 }
244 }
245 }
246
250 static Scalar computeUpdateError(const PrimaryVariables&,
251 const EqVector&)
252 {
253 // do not consider the change of Brine primary variables for
254 // convergence
255 // TODO: maybe this should be changed
256 return static_cast<Scalar>(0.0);
257 }
258
259 template <class DofEntity>
260 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
261 {
262 if constexpr (enableBrine) {
263 const unsigned dofIdx = model.dofMapper().index(dof);
264 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
265 outstream << priVars[saltConcentrationIdx];
266 }
267 }
268
269 template <class DofEntity>
270 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
271 {
272 if constexpr (enableBrine) {
273 const unsigned dofIdx = model.dofMapper().index(dof);
274 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
275 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
276
277 instream >> priVars0[saltConcentrationIdx];
278
279 // set the primary variables for the beginning of the current time step.
280 priVars1[saltConcentrationIdx] = priVars0[saltConcentrationIdx];
281 }
282 }
283
284 static Scalar referencePressure(const ElementContext& elemCtx,
285 unsigned scvIdx,
286 unsigned timeIdx)
287 {
288 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
289 return params_.referencePressure_[pvtnumRegionIdx];
290 }
291
292 static const TabulatedFunction& bdensityTable(const ElementContext& elemCtx,
293 unsigned scvIdx,
294 unsigned timeIdx)
295 {
296 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
297 return params_.bdensityTable_[pvtnumRegionIdx];
298 }
299
300 static const TabulatedFunction& pcfactTable(unsigned satnumRegionIdx)
301 { return params_.pcfactTable_[satnumRegionIdx]; }
302
303 static const TabulatedFunction& permfactTable(const ElementContext& elemCtx,
304 unsigned scvIdx,
305 unsigned timeIdx)
306 {
307 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
308 return params_.permfactTable_[satnumRegionIdx];
309 }
310
311 static const TabulatedFunction& permfactTable(unsigned satnumRegionIdx)
312 { return params_.permfactTable_[satnumRegionIdx]; }
313
314 static Scalar saltsolTable(const ElementContext& elemCtx,
315 unsigned scvIdx,
316 unsigned timeIdx)
317 {
318 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
319 return params_.saltsolTable_[pvtnumRegionIdx];
320 }
321
322 static Scalar saltsolTable(const unsigned pvtnumRegionIdx)
323 {
324 return params_.saltsolTable_[pvtnumRegionIdx];
325 }
326
327 static Scalar saltdenTable(const ElementContext& elemCtx,
328 unsigned scvIdx,
329 unsigned timeIdx)
330 {
331 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
332 return params_.saltdenTable_[pvtnumRegionIdx];
333 }
334
335 static Scalar saltdenTable(const unsigned pvtnumRegionIdx)
336 {
337 return params_.saltdenTable_[pvtnumRegionIdx];
338 }
339
340 static bool hasBDensityTables()
341 { return !params_.bdensityTable_.empty(); }
342
343 static bool hasSaltsolTables()
344 { return !params_.saltsolTable_.empty(); }
345
346 static bool hasPcfactTables()
347 {
348 if constexpr (enableSaltPrecipitation) {
349 return !params_.pcfactTable_.empty();
350 }
351 else {
352 return false;
353 }
354 }
355
356 static Scalar saltSol(unsigned regionIdx)
357 { return params_.saltsolTable_[regionIdx]; }
358
359private:
360 static BlackOilBrineParams<Scalar> params_;
361};
362
363template <class TypeTag, bool enableBrineV>
364BlackOilBrineParams<typename BlackOilBrineModule<TypeTag, enableBrineV>::Scalar>
365BlackOilBrineModule<TypeTag, enableBrineV>::params_;
366
367template <class TypeTag, bool enableBrineV>
369
377template <class TypeTag>
378class BlackOilBrineIntensiveQuantities<TypeTag, /*enableBrineV=*/true>
379{
381
389
391
392 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
393 static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
394 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
395 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
396 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
397 static constexpr bool enableBrine = true;
398 static constexpr bool enableSaltPrecipitation =
399 getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
400 static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
401
402public:
408 void updateSaltConcentration_(const ElementContext& elemCtx,
409 unsigned dofIdx,
410 unsigned timeIdx)
411 {
412 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
413 const LinearizationType lintype = elemCtx.linearizationType();
414 updateSaltConcentration_(priVars, timeIdx, lintype);
415 }
416
417 void updateSaltConcentration_(const PrimaryVariables& priVars,
418 const unsigned timeIdx,
419 const LinearizationType lintype)
420 {
421 const unsigned pvtnumRegionIdx = priVars.pvtRegionIndex();
422 auto& fs = asImp_().fluidState_;
423
424 if constexpr (enableSaltPrecipitation) {
425 saltSolubility_ = BrineModule::saltsolTable(pvtnumRegionIdx);
426 saltDensity_ = BrineModule::saltdenTable(pvtnumRegionIdx);
427
428 if (priVars.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
429 saltSaturation_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
430 fs.setSaltConcentration(saltSolubility_);
431 }
432 else {
433 saltConcentration_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
434 fs.setSaltConcentration(saltConcentration_);
435 saltSaturation_ = 0.0;
436 }
437 fs.setSaltSaturation(saltSaturation_);
438 }
439 else {
440 saltConcentration_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
441 fs.setSaltConcentration(saltConcentration_);
442 }
443 }
444
445 void saltPropertiesUpdate_([[maybe_unused]] const ElementContext& elemCtx,
446 [[maybe_unused]] unsigned dofIdx,
447 [[maybe_unused]] unsigned timeIdx)
448 {
449 if constexpr (enableSaltPrecipitation) {
450 const Evaluation porosityFactor = min(1.0 - asImp_().fluidState_.saltSaturation(), 1.0); //phi/phi_0
451
452 const auto& permfactTable = BrineModule::permfactTable(elemCtx, dofIdx, timeIdx);
453
454 permFactor_ = permfactTable.eval(porosityFactor);
455 }
456 }
457
458 const Evaluation& brineRefDensity() const
459 { return refDensity_; }
460
461 Scalar saltSolubility() const
462 { return saltSolubility_; }
463
464 Scalar saltDensity() const
465 { return saltDensity_; }
466
467 const Evaluation& permFactor() const
468 { return permFactor_; }
469
470protected:
471 Implementation& asImp_()
472 { return *static_cast<Implementation*>(this); }
473
475 Evaluation refDensity_;
476 Evaluation saltSaturation_;
477 Evaluation permFactor_;
480};
481
482template <class TypeTag>
484{
488
489public:
490 void updateSaltConcentration_(const ElementContext&,
491 unsigned,
492 unsigned)
493 {}
494
495 void saltPropertiesUpdate_(const ElementContext&,
496 unsigned,
497 unsigned)
498 {}
499
500 const Evaluation& brineRefDensity() const
501 { throw std::runtime_error("brineRefDensity() called but brine are disabled"); }
502
503 const Scalar saltSolubility() const
504 { throw std::logic_error("saltSolubility() called but salt precipitation is disabled"); }
505
506 const Scalar saltDensity() const
507 { throw std::logic_error("saltDensity() called but salt precipitation is disabled"); }
508
509 const Evaluation& permFactor() const
510 { throw std::logic_error("permFactor() called but salt precipitation is disabled"); }
511};
512
513} // namespace Opm
514
515#endif
Defines a type tags and some fundamental properties all models.
Contains the parameters required to extend the black-oil model by brine.
Declares the properties required by the black oil model.
void updateSaltConcentration_(const ElementContext &, unsigned, unsigned)
Definition: blackoilbrinemodules.hh:490
const Scalar saltDensity() const
Definition: blackoilbrinemodules.hh:506
const Evaluation & permFactor() const
Definition: blackoilbrinemodules.hh:509
void saltPropertiesUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilbrinemodules.hh:495
const Scalar saltSolubility() const
Definition: blackoilbrinemodules.hh:503
const Evaluation & brineRefDensity() const
Definition: blackoilbrinemodules.hh:500
void updateSaltConcentration_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the intensive properties needed to handle brine from the primary variables.
Definition: blackoilbrinemodules.hh:408
Evaluation refDensity_
Definition: blackoilbrinemodules.hh:475
const Evaluation & brineRefDensity() const
Definition: blackoilbrinemodules.hh:458
Evaluation permFactor_
Definition: blackoilbrinemodules.hh:477
Evaluation saltConcentration_
Definition: blackoilbrinemodules.hh:474
Scalar saltSolubility_
Definition: blackoilbrinemodules.hh:478
const Evaluation & permFactor() const
Definition: blackoilbrinemodules.hh:467
Scalar saltDensity_
Definition: blackoilbrinemodules.hh:479
void updateSaltConcentration_(const PrimaryVariables &priVars, const unsigned timeIdx, const LinearizationType lintype)
Definition: blackoilbrinemodules.hh:417
Evaluation saltSaturation_
Definition: blackoilbrinemodules.hh:476
Scalar saltSolubility() const
Definition: blackoilbrinemodules.hh:461
Implementation & asImp_()
Definition: blackoilbrinemodules.hh:471
void saltPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:445
Scalar saltDensity() const
Definition: blackoilbrinemodules.hh:464
Definition: blackoilbrinemodules.hh:368
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:58
static bool hasBDensityTables()
Definition: blackoilbrinemodules.hh:340
static const TabulatedFunction & bdensityTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:292
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilbrinemodules.hh:260
static Scalar saltSol(unsigned regionIdx)
Definition: blackoilbrinemodules.hh:356
static std::string eqName(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:148
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:101
static const TabulatedFunction & pcfactTable(unsigned satnumRegionIdx)
Definition: blackoilbrinemodules.hh:300
static const TabulatedFunction & permfactTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:303
static bool hasPcfactTables()
Definition: blackoilbrinemodules.hh:346
static Scalar referencePressure(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:284
static Scalar saltsolTable(const unsigned pvtnumRegionIdx)
Definition: blackoilbrinemodules.hh:322
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:123
static void setParams(BlackOilBrineParams< Scalar > &&params)
Set parameters.
Definition: blackoilbrinemodules.hh:90
static void registerParameters()
Register all run-time parameters for the black-oil brine module.
Definition: blackoilbrinemodules.hh:98
static const TabulatedFunction & permfactTable(unsigned satnumRegionIdx)
Definition: blackoilbrinemodules.hh:311
static bool hasSaltsolTables()
Definition: blackoilbrinemodules.hh:343
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:155
static Scalar saltdenTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:327
static Scalar saltsolTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:314
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:130
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:200
static bool eqApplies(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:138
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:165
static Scalar saltdenTable(const unsigned pvtnumRegionIdx)
Definition: blackoilbrinemodules.hh:335
static void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the brine specific primary variables to a PrimaryVariables object.
Definition: blackoilbrinemodules.hh:115
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilbrinemodules.hh:270
static void addBrineFluxes_(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:218
static void addBrineFluxes_(RateVector &flux, unsigned phaseIdx, const Evaluation &volFlux, const FluidState &upFs)
Definition: blackoilbrinemodules.hh:232
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilbrinemodules.hh:250
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
Struct holding the parameters for the BlackoilBrineModule class.
Definition: blackoilbrineparams.hpp:42
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilbrineparams.hpp:46
Definition: linearizationtype.hh:34