28#ifndef EWOMS_BLACK_OIL_BRINE_MODULE_HH
29#define EWOMS_BLACK_OIL_BRINE_MODULE_HH
31#include <dune/common/fvector.hh>
54template <class TypeTag, bool enableBrineV = getPropValue<TypeTag, Properties::EnableBrine>()>
70 using Toolbox = MathToolbox<Evaluation>;
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>();
83 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
84 static constexpr unsigned numPhases = FluidSystem::numPhases;
101 if constexpr (enableBrine) {
102 return pvIdx == saltConcentrationIdx;
112 template <
class Flu
idState>
114 const FluidState& fluidState)
116 if constexpr (enableBrine) {
117 priVars[saltConcentrationIdx] = fluidState.saltConcentration();
125 return "saltConcentration";
133 return static_cast<Scalar
>(1.0);
138 if constexpr (enableBrine) {
139 return eqIdx == contiBrineEqIdx;
146 static std::string
eqName([[maybe_unused]]
unsigned eqIdx)
150 return "conti^brine";
153 static Scalar
eqWeight([[maybe_unused]]
unsigned eqIdx)
158 return static_cast<Scalar
>(1.0);
162 template <
class LhsEval>
163 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
164 const IntensiveQuantities& intQuants)
166 if constexpr (enableBrine) {
167 const auto& fs = intQuants.fluidState();
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()),
177 const LhsEval massBrine = surfaceVolumeWater *
178 Toolbox::template decay<LhsEval>(fs.saltConcentration());
180 if constexpr (enableSaltPrecipitation) {
181 const double saltDensity = intQuants.saltDensity();
182 const LhsEval solidSalt =
183 Toolbox::template decay<LhsEval>(intQuants.porosity()) /
184 (1.0 - Toolbox::template decay<LhsEval>(fs.saltSaturation()) + 1.e-8) *
186 Toolbox::template decay<LhsEval>(fs.saltSaturation());
188 storage[contiBrineEqIdx] += massBrine + solidSalt;
191 storage[contiBrineEqIdx] += massBrine;
197 [[maybe_unused]]
const ElementContext& elemCtx,
198 [[maybe_unused]]
unsigned scvfIdx,
199 [[maybe_unused]]
unsigned timeIdx)
201 if constexpr (enableBrine) {
202 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
203 unsigned focusIdx = elemCtx.focusDofIndex();
204 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
205 flux[contiBrineEqIdx] = 0.0;
206 if (upIdx == focusIdx)
207 addBrineFluxes_<Evaluation>(flux, elemCtx, scvfIdx, timeIdx);
209 addBrineFluxes_<Scalar>(flux, elemCtx, scvfIdx, timeIdx);
213 template <
class UpstreamEval>
215 const ElementContext& elemCtx,
219 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
220 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
221 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
222 const auto& upFs = up.fluidState();
223 const auto& volFlux = extQuants.volumeFlux(waterPhaseIdx);
224 addBrineFluxes_<UpstreamEval>(flux, waterPhaseIdx, volFlux, upFs);
227 template <
class UpEval,
class Flu
idState>
230 const Evaluation& volFlux,
231 const FluidState& upFs)
233 if constexpr (enableBrine) {
234 if (phaseIdx == waterPhaseIdx) {
235 flux[contiBrineEqIdx] =
236 decay<UpEval>(upFs.saltConcentration())
237 * decay<UpEval>(upFs.invB(waterPhaseIdx))
252 return static_cast<Scalar
>(0.0);
255 template <
class DofEntity>
256 static void serializeEntity(
const Model& model, std::ostream& outstream,
const DofEntity& dof)
258 if constexpr (enableBrine) {
259 const unsigned dofIdx = model.dofMapper().index(dof);
260 const PrimaryVariables& priVars = model.solution(0)[dofIdx];
261 outstream << priVars[saltConcentrationIdx];
265 template <
class DofEntity>
268 if constexpr (enableBrine) {
269 const unsigned dofIdx = model.dofMapper().index(dof);
270 PrimaryVariables& priVars0 = model.solution(0)[dofIdx];
271 PrimaryVariables& priVars1 = model.solution(1)[dofIdx];
273 instream >> priVars0[saltConcentrationIdx];
276 priVars1[saltConcentrationIdx] = priVars0[saltConcentrationIdx];
284 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
285 return params_.referencePressure_[pvtnumRegionIdx];
288 static const TabulatedFunction&
bdensityTable(
const ElementContext& elemCtx,
292 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
293 return params_.bdensityTable_[pvtnumRegionIdx];
296 static const TabulatedFunction&
pcfactTable(
unsigned satnumRegionIdx)
297 {
return params_.pcfactTable_[satnumRegionIdx]; }
299 static const TabulatedFunction&
permfactTable(
const ElementContext& elemCtx,
303 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
304 return params_.permfactTable_[satnumRegionIdx];
308 {
return params_.permfactTable_[satnumRegionIdx]; }
314 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
315 return params_.saltsolTable_[pvtnumRegionIdx];
320 return params_.saltsolTable_[pvtnumRegionIdx];
327 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
328 return params_.saltdenTable_[pvtnumRegionIdx];
333 return params_.saltdenTable_[pvtnumRegionIdx];
337 {
return !params_.bdensityTable_.empty(); }
340 {
return !params_.saltsolTable_.empty(); }
344 if constexpr (enableSaltPrecipitation) {
345 return !params_.pcfactTable_.empty();
353 {
return params_.saltsolTable_[regionIdx]; }
359template <
class TypeTag,
bool enableBrineV>
360BlackOilBrineParams<typename BlackOilBrineModule<TypeTag, enableBrineV>::Scalar>
361BlackOilBrineModule<TypeTag, enableBrineV>::params_;
363template <
class TypeTag,
bool enableBrineV>
373template <
class TypeTag>
388 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
389 static constexpr int saltConcentrationIdx = Indices::saltConcentrationIdx;
390 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
391 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
392 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
393 static constexpr bool enableBrine =
true;
394 static constexpr bool enableSaltPrecipitation =
395 getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
396 static constexpr int contiBrineEqIdx = Indices::contiBrineEqIdx;
408 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
410 updateSaltConcentration_(priVars, timeIdx, lintype);
414 const unsigned timeIdx,
417 const unsigned pvtnumRegionIdx = priVars.pvtRegionIndex();
418 auto& fs = asImp_().fluidState_;
420 if constexpr (enableSaltPrecipitation) {
421 saltSolubility_ = BrineModule::saltsolTable(pvtnumRegionIdx);
422 saltDensity_ = BrineModule::saltdenTable(pvtnumRegionIdx);
424 if (priVars.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
425 saltSaturation_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
426 fs.setSaltConcentration(saltSolubility_);
429 saltConcentration_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
430 fs.setSaltConcentration(saltConcentration_);
431 saltSaturation_ = 0.0;
433 fs.setSaltSaturation(saltSaturation_);
436 saltConcentration_ = priVars.makeEvaluation(saltConcentrationIdx, timeIdx, lintype);
437 fs.setSaltConcentration(saltConcentration_);
442 [[maybe_unused]]
unsigned dofIdx,
443 [[maybe_unused]]
unsigned timeIdx)
445 if constexpr (enableSaltPrecipitation) {
446 const Evaluation porosityFactor = min(1.0 - asImp_().fluidState_.saltSaturation(), 1.0);
448 const auto& permfactTable = BrineModule::permfactTable(elemCtx, dofIdx, timeIdx);
450 permFactor_ = permfactTable.eval(porosityFactor);
451 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
452 if (!FluidSystem::phaseIsActive(phaseIdx)) {
456 asImp_().mobility_[phaseIdx] *= permFactor_;
462 {
return refDensity_; }
465 {
return saltSolubility_; }
468 {
return saltDensity_; }
471 {
return permFactor_; }
475 {
return *
static_cast<Implementation*
>(
this); }
485template <
class TypeTag>
504 {
throw std::runtime_error(
"brineRefDensity() called but brine are disabled"); }
507 {
throw std::logic_error(
"saltSolubility() called but salt precipitation is disabled"); }
510 {
throw std::logic_error(
"saltDensity() called but salt precipitation is disabled"); }
513 {
throw std::logic_error(
"permFactor() called but salt precipitation is disabled"); }
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:493
const Scalar saltDensity() const
Definition: blackoilbrinemodules.hh:509
const Evaluation & permFactor() const
Definition: blackoilbrinemodules.hh:512
void saltPropertiesUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilbrinemodules.hh:498
const Scalar saltSolubility() const
Definition: blackoilbrinemodules.hh:506
const Evaluation & brineRefDensity() const
Definition: blackoilbrinemodules.hh:503
void updateSaltConcentration_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the intensive properties needed to handle brine from the primary variables.
Definition: blackoilbrinemodules.hh:404
Evaluation refDensity_
Definition: blackoilbrinemodules.hh:478
const Evaluation & brineRefDensity() const
Definition: blackoilbrinemodules.hh:461
Evaluation permFactor_
Definition: blackoilbrinemodules.hh:480
Evaluation saltConcentration_
Definition: blackoilbrinemodules.hh:477
Scalar saltSolubility_
Definition: blackoilbrinemodules.hh:481
const Evaluation & permFactor() const
Definition: blackoilbrinemodules.hh:470
Scalar saltDensity_
Definition: blackoilbrinemodules.hh:482
void updateSaltConcentration_(const PrimaryVariables &priVars, const unsigned timeIdx, const LinearizationType lintype)
Definition: blackoilbrinemodules.hh:413
Evaluation saltSaturation_
Definition: blackoilbrinemodules.hh:479
Scalar saltSolubility() const
Definition: blackoilbrinemodules.hh:464
Implementation & asImp_()
Definition: blackoilbrinemodules.hh:474
void saltPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:441
Scalar saltDensity() const
Definition: blackoilbrinemodules.hh:467
Definition: blackoilbrinemodules.hh:364
Contains the high level supplements required to extend the black oil model by brine.
Definition: blackoilbrinemodules.hh:56
static bool hasBDensityTables()
Definition: blackoilbrinemodules.hh:336
static const TabulatedFunction & bdensityTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:288
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilbrinemodules.hh:256
static Scalar saltSol(unsigned regionIdx)
Definition: blackoilbrinemodules.hh:352
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:296
static const TabulatedFunction & permfactTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:299
static bool hasPcfactTables()
Definition: blackoilbrinemodules.hh:342
static Scalar referencePressure(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:280
static Scalar saltsolTable(const unsigned pvtnumRegionIdx)
Definition: blackoilbrinemodules.hh:318
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilbrinemodules.hh:121
static void setParams(BlackOilBrineParams< Scalar > &¶ms)
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:307
static bool hasSaltsolTables()
Definition: blackoilbrinemodules.hh:339
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilbrinemodules.hh:153
static Scalar saltdenTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:323
static Scalar saltsolTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:310
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:331
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:266
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilbrinemodules.hh:163
static void addBrineFluxes_(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilbrinemodules.hh:214
static void addBrineFluxes_(RateVector &flux, unsigned phaseIdx, const Evaluation &volFlux, const FluidState &upFs)
Definition: blackoilbrinemodules.hh:228
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilbrinemodules.hh:246
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilbioeffectsmodules.hh:43
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