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