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