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
38
41
43
44#include <cassert>
45#include <istream>
46#include <ostream>
47#include <string>
48
49namespace Opm {
50
56template <class TypeTag>
57class BlackOilBrineModule<TypeTag, true>
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 = true;
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 OPM_HOST_DEVICE bool primaryVarApplies(unsigned pvIdx)
102 {
103 return pvIdx == saltConcentrationIdx;
104 }
105
109 template <class FluidState>
110 static void assignPrimaryVars(PrimaryVariables& priVars,
111 const FluidState& fluidState)
112 {
113 priVars[saltConcentrationIdx] = fluidState.saltConcentration();
114 }
115
116 static std::string primaryVarName([[maybe_unused]] unsigned pvIdx)
117 {
118 assert(primaryVarApplies(pvIdx));
119
120 return "saltConcentration";
121 }
122
123 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
124 {
125 assert(primaryVarApplies(pvIdx));
126
127 // TODO: it may be beneficial to choose this differently.
128 return static_cast<Scalar>(1.0);
129 }
130
131 static OPM_HOST_DEVICE bool eqApplies(unsigned eqIdx)
132 {
133 return eqIdx == contiBrineEqIdx;
134 }
135
136 static std::string eqName([[maybe_unused]] unsigned eqIdx)
137 {
138 assert(eqApplies(eqIdx));
139
140 return "conti^brine";
141 }
142
143 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
144 {
145 assert(eqApplies(eqIdx));
146
147 // TODO: it may be beneficial to choose this differently.
148 return static_cast<Scalar>(1.0);
149 }
150
151 // must be called after water storage is computed
152 template <class StorageType>
153 OPM_HOST_DEVICE static void addStorage(StorageType& storage,
154 const IntensiveQuantities& intQuants)
155 {
156 using LhsEval = typename StorageType::value_type;
157
158 const auto& fs = intQuants.fluidState();
159
160 // avoid singular matrix if no water is present.
161 const LhsEval surfaceVolumeWater =
162 max(Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx)) *
163 Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx)) *
164 Toolbox::template decay<LhsEval>(intQuants.porosity()),
165 1e-10);
166
167 // Brine in water phase
168 const LhsEval massBrine = surfaceVolumeWater *
169 Toolbox::template decay<LhsEval>(fs.saltConcentration());
170
171 if constexpr (enableSaltPrecipitation) {
172 const double saltDensity = intQuants.saltDensity(); // Solid salt density kg/m3
173 const LhsEval solidSalt =
174 Toolbox::template decay<LhsEval>(intQuants.porosity()) /
175 (1.0 - Toolbox::template decay<LhsEval>(fs.saltSaturation()) + 1.e-8) *
176 saltDensity *
177 Toolbox::template decay<LhsEval>(fs.saltSaturation());
178
179 storage[contiBrineEqIdx] += massBrine + solidSalt;
180 }
181 else {
182 storage[contiBrineEqIdx] += massBrine;
183 }
184 }
185
186 static void computeFlux(RateVector& flux,
187 const ElementContext& elemCtx,
188 unsigned scvfIdx,
189 unsigned timeIdx)
190 {
191 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
192 unsigned focusIdx = elemCtx.focusDofIndex();
193 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
194 flux[contiBrineEqIdx] = 0.0;
195 if (upIdx == focusIdx) {
196 addBrineFluxes_<Evaluation>(flux, elemCtx, scvfIdx, timeIdx);
197 }
198 else {
199 addBrineFluxes_<Scalar>(flux, elemCtx, scvfIdx, timeIdx);
200 }
201 }
202
203 template <class UpstreamEval>
204 static void addBrineFluxes_(RateVector& flux,
205 const ElementContext& elemCtx,
206 unsigned scvfIdx,
207 unsigned timeIdx)
208 {
209 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
210 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
211 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
212 const auto& upFs = up.fluidState();
213 const auto& volFlux = extQuants.volumeFlux(waterPhaseIdx);
214 addBrineFluxes_<UpstreamEval>(flux, waterPhaseIdx, volFlux, upFs);
215 }
216
217 template <class UpEval, class FluidState>
218 static void addBrineFluxes_(RateVector& flux,
219 unsigned phaseIdx,
220 const Evaluation& volFlux,
221 const FluidState& upFs)
222 {
223 if (phaseIdx == waterPhaseIdx) {
224 flux[contiBrineEqIdx] =
225 decay<UpEval>(upFs.saltConcentration())
226 * decay<UpEval>(upFs.invB(waterPhaseIdx))
227 * volFlux;
228 }
229 }
230
234 static Scalar computeUpdateError(const PrimaryVariables&,
235 const EqVector&)
236 {
237 // do not consider the change of Brine primary variables for
238 // convergence
239 // TODO: maybe this should be changed
240 return static_cast<Scalar>(0.0);
241 }
242
243 template <class DofEntity>
244 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
245 {
246 const unsigned dofIdx = model.dofMapper().index(dof);
247 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
248 outstream << priVars[saltConcentrationIdx];
249 }
250
251 template <class DofEntity>
252 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
253 {
254 const unsigned dofIdx = model.dofMapper().index(dof);
255 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
256 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
257
258 instream >> priVars0[saltConcentrationIdx];
259
260 // set the primary variables for the beginning of the current time step.
261 priVars1[saltConcentrationIdx] = priVars0[saltConcentrationIdx];
262 }
263
264 static Scalar referencePressure(const ElementContext& elemCtx,
265 unsigned scvIdx,
266 unsigned timeIdx)
267 {
268 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
269 return params_.referencePressure_[pvtnumRegionIdx];
270 }
271
272 static const TabulatedFunction& bdensityTable(const ElementContext& elemCtx,
273 unsigned scvIdx,
274 unsigned timeIdx)
275 {
276 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
277 return params_.bdensityTable_[pvtnumRegionIdx];
278 }
279
280 static const TabulatedFunction& pcfactTable(unsigned satnumRegionIdx)
281 { return params_.pcfactTable_[satnumRegionIdx]; }
282
283 static const TabulatedFunction& permfactTable(const ElementContext& elemCtx,
284 unsigned scvIdx,
285 unsigned timeIdx)
286 {
287 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
288 return params_.permfactTable_[satnumRegionIdx];
289 }
290
291 static const TabulatedFunction& permfactTable(unsigned satnumRegionIdx)
292 { return params_.permfactTable_[satnumRegionIdx]; }
293
294 static Scalar saltsolTable(const ElementContext& elemCtx,
295 unsigned scvIdx,
296 unsigned timeIdx)
297 {
298 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
299 return params_.saltsolTable_[pvtnumRegionIdx];
300 }
301
302 static Scalar saltsolTable(const unsigned pvtnumRegionIdx)
303 {
304 return params_.saltsolTable_[pvtnumRegionIdx];
305 }
306
307 static Scalar saltdenTable(const ElementContext& elemCtx,
308 unsigned scvIdx,
309 unsigned timeIdx)
310 {
311 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
312 return params_.saltdenTable_[pvtnumRegionIdx];
313 }
314
315 static Scalar saltdenTable(const unsigned pvtnumRegionIdx)
316 {
317 return params_.saltdenTable_[pvtnumRegionIdx];
318 }
319
320 static bool hasBDensityTables()
321 { return !params_.bdensityTable_.empty(); }
322
323 static bool hasSaltsolTables()
324 { return !params_.saltsolTable_.empty(); }
325
326 static bool hasPcfactTables()
327 {
328 if constexpr (enableSaltPrecipitation) {
329 return !params_.pcfactTable_.empty();
330 }
331 else {
332 return false;
333 }
334 }
335
336 static Scalar saltSol(unsigned regionIdx)
337 { return params_.saltsolTable_[regionIdx]; }
338
339private:
340 static BlackOilBrineParams<Scalar> params_;
341};
342
343template <class TypeTag>
344BlackOilBrineParams<typename BlackOilBrineModule<TypeTag, true>::Scalar>
345BlackOilBrineModule<TypeTag, true>::params_;
346
354template <class TypeTag>
355class BlackOilBrineIntensiveQuantities<TypeTag, /*enableBrineV=*/true>
356{
358
366
368
369 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
370 static constexpr unsigned 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 OPM_HOST_DEVICE 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 OPM_HOST_DEVICE 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 OPM_HOST_DEVICE 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 - asImp_().fluidState_.saltSaturation(), 1.0); //phi/phi_0
428
429 const auto& permfactTable = BrineModule::permfactTable(elemCtx, dofIdx, timeIdx);
430
431 permFactor_ = permfactTable.eval(porosityFactor);
432 }
433 }
434
435 OPM_HOST_DEVICE const Evaluation& brineRefDensity() const
436 { return refDensity_; }
437
438 OPM_HOST_DEVICE Scalar saltSolubility() const
439 { return saltSolubility_; }
440
441 OPM_HOST_DEVICE Scalar saltDensity() const
442 { return saltDensity_; }
443
444 OPM_HOST_DEVICE const Evaluation& permFactor() const
445 { return permFactor_; }
446
447protected:
448 OPM_HOST_DEVICE Implementation& asImp_()
449 { return *static_cast<Implementation*>(this); }
450
452 Evaluation refDensity_;
453 Evaluation saltSaturation_;
454 Evaluation permFactor_;
457};
458
459} // namespace Opm
460
461#endif
Defines a type tags and some fundamental properties all models.
Contains the parameters required to extend the black-oil model by brine.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
Provides the volumetric quantities required for the equations needed by the brine extension of the bl...
OPM_HOST_DEVICE Implementation & asImp_()
Definition: blackoilbrinemodules.hh:448
OPM_HOST_DEVICE 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:452
Evaluation permFactor_
Definition: blackoilbrinemodules.hh:454
Evaluation saltConcentration_
Definition: blackoilbrinemodules.hh:451
Scalar saltSolubility_
Definition: blackoilbrinemodules.hh:455
OPM_HOST_DEVICE void updateSaltConcentration_(const PrimaryVariables &priVars, const unsigned timeIdx, const LinearizationType lintype)
Definition: blackoilbrinemodules.hh:394
OPM_HOST_DEVICE Scalar saltDensity() const
Definition: blackoilbrinemodules.hh:441
OPM_HOST_DEVICE const Evaluation & permFactor() const
Definition: blackoilbrinemodules.hh:444
Scalar saltDensity_
Definition: blackoilbrinemodules.hh:456
Evaluation saltSaturation_
Definition: blackoilbrinemodules.hh:453
OPM_HOST_DEVICE const Evaluation & brineRefDensity() const
Definition: blackoilbrinemodules.hh:435
OPM_HOST_DEVICE void saltPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:422
OPM_HOST_DEVICE Scalar saltSolubility() const
Definition: blackoilbrinemodules.hh:438
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:58
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilbrinemodules.hh:252
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:123
static void registerParameters()
Register all run-time parameters for the black-oil brine module.
Definition: blackoilbrinemodules.hh:98
static void setParams(BlackOilBrineParams< Scalar > &&params)
Set parameters.
Definition: blackoilbrinemodules.hh:90
static const TabulatedFunction & pcfactTable(unsigned satnumRegionIdx)
Definition: blackoilbrinemodules.hh:280
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:143
static std::string eqName(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:136
static const TabulatedFunction & permfactTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:283
static OPM_HOST_DEVICE bool eqApplies(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:131
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:116
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilbrinemodules.hh:244
static void addBrineFluxes_(RateVector &flux, unsigned phaseIdx, const Evaluation &volFlux, const FluidState &upFs)
Definition: blackoilbrinemodules.hh:218
static const TabulatedFunction & bdensityTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:272
static void assignPrimaryVars(PrimaryVariables &priVars, const FluidState &fluidState)
Assign the brine specific primary variables to a PrimaryVariables object.
Definition: blackoilbrinemodules.hh:110
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilbrinemodules.hh:234
static const TabulatedFunction & permfactTable(unsigned satnumRegionIdx)
Definition: blackoilbrinemodules.hh:291
static Scalar saltsolTable(const unsigned pvtnumRegionIdx)
Definition: blackoilbrinemodules.hh:302
static Scalar saltdenTable(const unsigned pvtnumRegionIdx)
Definition: blackoilbrinemodules.hh:315
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:186
static OPM_HOST_DEVICE bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:101
static bool hasPcfactTables()
Definition: blackoilbrinemodules.hh:326
static bool hasSaltsolTables()
Definition: blackoilbrinemodules.hh:323
static Scalar saltdenTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:307
static Scalar referencePressure(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:264
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:153
static Scalar saltsolTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:294
static bool hasBDensityTables()
Definition: blackoilbrinemodules.hh:320
static Scalar saltSol(unsigned regionIdx)
Definition: blackoilbrinemodules.hh:336
static void addBrineFluxes_(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:204
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