28#ifndef EWOMS_BLACK_OIL_POLYMER_MODULE_HH
29#define EWOMS_BLACK_OIL_POLYMER_MODULE_HH
31#include <dune/common/fvector.hh>
33#include <opm/common/OpmLog/OpmLog.hpp>
35#include <opm/material/common/MathToolbox.hpp>
62template <class TypeTag, bool enablePolymerV = getPropValue<TypeTag, Properties::EnablePolymer>()>
77 using Toolbox = MathToolbox<Evaluation>;
82 static constexpr unsigned polymerConcentrationIdx = Indices::polymerConcentrationIdx;
83 static constexpr unsigned polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
84 static constexpr unsigned contiPolymerEqIdx = Indices::contiPolymerEqIdx;
85 static constexpr unsigned contiPolymerMolarWeightEqIdx = Indices::contiPolymerMWEqIdx;
86 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
88 static constexpr unsigned enablePolymer = enablePolymerV;
89 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
91 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
105 const auto iterTable = params_.plymwinjTables_.find(tableNumber);
106 if (iterTable != params_.plymwinjTables_.end()) {
107 return iterTable->second;
110 throw std::runtime_error(
" the PLYMWINJ table " +
std::to_string(tableNumber) +
" does not exist\n");
119 const auto iterTable = params_.skprwatTables_.find(tableNumber);
120 if (iterTable != params_.skprwatTables_.end()) {
121 return iterTable->second;
124 throw std::runtime_error(
" the SKPRWAT table " +
std::to_string(tableNumber) +
" does not exist\n");
134 const auto iterTable = params_.skprpolyTables_.find(tableNumber);
135 if (iterTable != params_.skprpolyTables_.end()) {
136 return iterTable->second;
139 throw std::runtime_error(
" the SKPRPOLY table " +
std::to_string(tableNumber) +
" does not exist\n");
148 if constexpr (enablePolymer) {
157 Simulator& simulator)
159 if constexpr (enablePolymer) {
166 if constexpr (enablePolymer) {
167 if constexpr (enablePolymerMolarWeight) {
168 return pvIdx == polymerConcentrationIdx || pvIdx == polymerMoleWeightIdx;
171 return pvIdx == polymerConcentrationIdx;
183 if (pvIdx == polymerConcentrationIdx) {
184 return "polymer_waterconcentration";
187 return "polymer_molecularweight";
196 return static_cast<Scalar
>(1.0);
201 if constexpr (enablePolymer) {
202 if constexpr (enablePolymerMolarWeight) {
203 return eqIdx == contiPolymerEqIdx || eqIdx == contiPolymerMolarWeightEqIdx;
206 return eqIdx == contiPolymerEqIdx;
214 static std::string
eqName(
unsigned eqIdx)
218 if (eqIdx == contiPolymerEqIdx) {
219 return "conti^polymer";
222 return "conti^polymer_molecularweight";
226 static Scalar
eqWeight([[maybe_unused]]
unsigned eqIdx)
231 return static_cast<Scalar
>(1.0);
235 template <
class LhsEval>
236 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
237 const IntensiveQuantities& intQuants)
239 if constexpr (enablePolymer) {
240 const auto& fs = intQuants.fluidState();
243 const LhsEval surfaceVolumeWater =
244 max(Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx)) *
245 Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx)) *
246 Toolbox::template decay<LhsEval>(intQuants.porosity()),
250 const LhsEval massPolymer =
252 Toolbox::template decay<LhsEval>(intQuants.polymerConcentration()) *
253 (1.0 - Toolbox::template decay<LhsEval>(intQuants.polymerDeadPoreVolume()));
256 const LhsEval adsorptionPolymer =
257 Toolbox::template decay<LhsEval>(1.0 - intQuants.porosity()) *
258 Toolbox::template decay<LhsEval>(intQuants.polymerRockDensity()) *
259 Toolbox::template decay<LhsEval>(intQuants.polymerAdsorption());
261 LhsEval accumulationPolymer = massPolymer + adsorptionPolymer;
263 storage[contiPolymerEqIdx] += accumulationPolymer;
266 if constexpr (enablePolymerMolarWeight) {
267 accumulationPolymer = max(accumulationPolymer, 1e-10);
269 storage[contiPolymerMolarWeightEqIdx] +=
270 accumulationPolymer * Toolbox::template decay<LhsEval>(intQuants.polymerMoleWeight());
276 [[maybe_unused]]
const ElementContext& elemCtx,
277 [[maybe_unused]]
unsigned scvfIdx,
278 [[maybe_unused]]
unsigned timeIdx)
280 if constexpr (enablePolymer) {
281 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
283 const unsigned upIdx = extQuants.upstreamIndex(FluidSystem::waterPhaseIdx);
284 const unsigned inIdx = extQuants.interiorIndex();
285 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
286 const unsigned contiWaterEqIdx =
287 Indices::conti0EqIdx + FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
289 if (upIdx == inIdx) {
290 flux[contiPolymerEqIdx] =
291 extQuants.volumeFlux(waterPhaseIdx) *
292 up.fluidState().invB(waterPhaseIdx) *
293 up.polymerViscosityCorrection() /
294 extQuants.polymerShearFactor() *
295 up.polymerConcentration();
298 flux[contiWaterEqIdx] /= extQuants.waterShearFactor();
301 flux[contiPolymerEqIdx] =
302 extQuants.volumeFlux(waterPhaseIdx) *
303 decay<Scalar>(up.fluidState().invB(waterPhaseIdx)) *
304 decay<Scalar>(up.polymerViscosityCorrection()) /
305 decay<Scalar>(extQuants.polymerShearFactor()) *
306 decay<Scalar>(up.polymerConcentration());
309 flux[contiWaterEqIdx] /= decay<Scalar>(extQuants.waterShearFactor());
313 if constexpr (enablePolymerMolarWeight) {
314 if (upIdx == inIdx) {
315 flux[contiPolymerMolarWeightEqIdx] =
316 flux[contiPolymerEqIdx] * up.polymerMoleWeight();
319 flux[contiPolymerMolarWeightEqIdx] =
320 flux[contiPolymerEqIdx] * decay<Scalar>(up.polymerMoleWeight());
335 return static_cast<Scalar
>(0.0);
338 template <
class DofEntity>
339 static void serializeEntity(
const Model& model, std::ostream& outstream,
const DofEntity& dof)
341 if constexpr (enablePolymer) {
342 const unsigned dofIdx = model.dofMapper().index(dof);
343 const PrimaryVariables& priVars = model.solution(0)[dofIdx];
344 outstream << priVars[polymerConcentrationIdx];
345 outstream << priVars[polymerMoleWeightIdx];
349 template <
class DofEntity>
352 if constexpr (enablePolymer) {
353 const unsigned dofIdx = model.dofMapper().index(dof);
354 PrimaryVariables& priVars0 = model.solution(0)[dofIdx];
355 PrimaryVariables& priVars1 = model.solution(1)[dofIdx];
357 instream >> priVars0[polymerConcentrationIdx];
358 instream >> priVars0[polymerMoleWeightIdx];
361 priVars1[polymerConcentrationIdx] = priVars0[polymerConcentrationIdx];
362 priVars1[polymerMoleWeightIdx] = priVars0[polymerMoleWeightIdx];
370 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
371 return params_.plyrockDeadPoreVolume_[satnumRegionIdx];
378 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
379 return params_.plyrockResidualResistanceFactor_[satnumRegionIdx];
386 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
387 return params_.plyrockRockDensityFactor_[satnumRegionIdx];
394 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
395 return params_.plyrockAdsorbtionIndex_[satnumRegionIdx];
402 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
403 return params_.plyrockMaxAdsorbtion_[satnumRegionIdx];
406 static const TabulatedFunction&
411 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
412 return params_.plyadsAdsorbedPolymer_[satnumRegionIdx];
415 static const TabulatedFunction&
420 unsigned pvtnumRegionIdx =
421 elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
422 return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
425 static const TabulatedFunction&
427 {
return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx]; }
433 const unsigned polymerMixRegionIdx =
434 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
435 return params_.plymaxMaxConcentration_[polymerMixRegionIdx];
442 const unsigned polymerMixRegionIdx =
443 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
444 return params_.plymixparToddLongstaff_[polymerMixRegionIdx];
449 const unsigned scvIdx,
450 const unsigned timeIdx)
452 const unsigned polymerMixRegionIdx =
453 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
454 return params_.plyvmhCoefficients_[polymerMixRegionIdx];
458 {
return params_.hasPlyshlog_; }
461 {
return params_.hasShrate_; }
463 static Scalar
shrate(
unsigned pvtnumRegionIdx)
464 {
return params_.shrate_[pvtnumRegionIdx]; }
472 template <
class Evaluation>
474 unsigned pvtnumRegionIdx,
475 const Evaluation& v0)
477 using ToolboxLocal = MathToolbox<Evaluation>;
479 const auto& viscosityMultiplierTable = params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
480 const Scalar viscosityMultiplier =
481 viscosityMultiplierTable.eval(scalarValue(polymerConcentration),
true);
483 const Scalar eps = 1e-14;
485 if (std::abs((viscosityMultiplier - 1.0)) < eps) {
486 return ToolboxLocal::createConstant(v0, 1.0);
489 const std::vector<Scalar>& shearEffectRefLogVelocity =
490 params_.plyshlogShearEffectRefLogVelocity_[pvtnumRegionIdx];
491 const auto v0AbsLog = log(abs(v0));
493 if (v0AbsLog < shearEffectRefLogVelocity[0]) {
494 return ToolboxLocal::createConstant(v0, 1.0);
501 const std::vector<Scalar>& shearEffectRefMultiplier =
502 params_.plyshlogShearEffectRefMultiplier_[pvtnumRegionIdx];
503 const std::size_t numTableEntries = shearEffectRefLogVelocity.size();
504 assert(shearEffectRefMultiplier.size() == numTableEntries);
506 std::vector<Scalar> shearEffectMultiplier(numTableEntries, 1.0);
507 for (std::size_t i = 0; i < numTableEntries; ++i) {
508 shearEffectMultiplier[i] = (1.0 + (viscosityMultiplier - 1.0) *
509 shearEffectRefMultiplier[i]) / viscosityMultiplier;
510 shearEffectMultiplier[i] = log(shearEffectMultiplier[i]);
514 const TabulatedFunction logShearEffectMultiplier =
515 TabulatedFunction(numTableEntries, shearEffectRefLogVelocity,
516 shearEffectMultiplier,
false);
523 auto F = [&logShearEffectMultiplier, &v0AbsLog](
const Evaluation& u) {
524 return u + logShearEffectMultiplier.eval(u,
true) - v0AbsLog;
527 auto dF = [&logShearEffectMultiplier](
const Evaluation& u) {
528 return 1 + logShearEffectMultiplier.evalDerivative(u,
true);
534 bool converged =
false;
536 for (
int i = 0; i < 20; ++i) {
538 const auto df = dF(u);
540 if (std::abs(scalarValue(f)) < 1e-12) {
546 throw std::runtime_error(
"Not able to compute shear velocity. \n");
550 return exp(logShearEffectMultiplier.eval(u,
true));
560template <
class TypeTag,
bool enablePolymerV>
561BlackOilPolymerParams<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
562BlackOilPolymerModule<TypeTag, enablePolymerV>::params_;
564template <
class TypeTag,
bool enablePolymerV>
574template <
class TypeTag>
589 static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
590 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
591 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
592 static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
604 const auto linearizationType = elemCtx.linearizationType();
605 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
606 polymerConcentration_ = priVars.makeEvaluation(polymerConcentrationIdx, timeIdx, linearizationType);
607 if constexpr (enablePolymerMolarWeight) {
608 polymerMoleWeight_ = priVars.makeEvaluation(polymerMoleWeightIdx, timeIdx, linearizationType);
612 const Scalar& maxAdsorbtion = PolymerModule::plyrockMaxAdsorbtion(elemCtx, dofIdx, timeIdx);
613 const auto& plyadsAdsorbedPolymer = PolymerModule::plyadsAdsorbedPolymer(elemCtx, dofIdx, timeIdx);
614 polymerAdsorption_ = plyadsAdsorbedPolymer.eval(polymerConcentration_,
true);
615 if (
static_cast<int>(PolymerModule::plyrockAdsorbtionIndex(elemCtx, dofIdx, timeIdx)) ==
618 const auto maxPolymerAdsorption =
619 elemCtx.problem().maxPolymerAdsorption(elemCtx, dofIdx, timeIdx);
620 polymerAdsorption_ = std::max(Evaluation(maxPolymerAdsorption), polymerAdsorption_);
624 const Scalar& residualResistanceFactor =
625 PolymerModule::plyrockResidualResistanceFactor(elemCtx, dofIdx, timeIdx);
626 const Evaluation resistanceFactor = 1.0 + (residualResistanceFactor - 1.0) *
627 polymerAdsorption_ / maxAdsorbtion;
630 if constexpr (!enablePolymerMolarWeight) {
631 const Scalar cmax = PolymerModule::plymaxMaxConcentration(elemCtx, dofIdx, timeIdx);
632 const auto& fs = asImp_().fluidState_;
633 const Evaluation& muWater = fs.viscosity(waterPhaseIdx);
634 const auto& viscosityMultiplier =
635 PolymerModule::plyviscViscosityMultiplierTable(elemCtx, dofIdx, timeIdx);
636 const Evaluation viscosityMixture =
637 viscosityMultiplier.eval(polymerConcentration_,
true) * muWater;
640 const Scalar plymixparToddLongstaff = PolymerModule::plymixparToddLongstaff(elemCtx, dofIdx, timeIdx);
641 const Evaluation viscosityPolymer = viscosityMultiplier.eval(cmax,
true) * muWater;
642 const Evaluation viscosityPolymerEffective =
643 pow(viscosityMixture, plymixparToddLongstaff) * pow(viscosityPolymer, 1.0 - plymixparToddLongstaff);
644 const Evaluation viscosityWaterEffective =
645 pow(viscosityMixture, plymixparToddLongstaff) * pow(muWater, 1.0 - plymixparToddLongstaff);
647 const Evaluation cbar = polymerConcentration_ / cmax;
649 waterViscosityCorrection_ = muWater * ((1.0 - cbar) / viscosityWaterEffective + cbar / viscosityPolymerEffective);
651 polymerViscosityCorrection_ = (muWater / waterViscosityCorrection_) / viscosityPolymerEffective;
654 const auto& plyvmhCoefficients = PolymerModule::plyvmhCoefficients(elemCtx, dofIdx, timeIdx);
655 const Scalar k_mh = plyvmhCoefficients.k_mh;
656 const Scalar a_mh = plyvmhCoefficients.a_mh;
657 const Scalar gamma = plyvmhCoefficients.gamma;
658 const Scalar kappa = plyvmhCoefficients.kappa;
662 const Evaluation intrinsicViscosity = k_mh * pow(polymerMoleWeight_ * 1000., a_mh);
663 const Evaluation x = polymerConcentration_ * intrinsicViscosity;
664 waterViscosityCorrection_ = 1.0 / (1.0 + gamma * (x + kappa * x * x));
665 polymerViscosityCorrection_ = 1.0;
669 asImp_().mobility_[waterPhaseIdx] *= waterViscosityCorrection_ / resistanceFactor;
672 polymerDeadPoreVolume_ = PolymerModule::plyrockDeadPoreVolume(elemCtx, dofIdx, timeIdx);
673 polymerRockDensity_ = PolymerModule::plyrockRockDensityFactor(elemCtx, dofIdx, timeIdx);
677 {
return polymerConcentration_; }
681 if constexpr (enablePolymerMolarWeight) {
682 return polymerMoleWeight_;
685 throw std::logic_error(
"polymerMoleWeight() is called but polymer milecular weight is disabled");
690 {
return polymerDeadPoreVolume_; }
693 {
return polymerAdsorption_; }
696 {
return polymerRockDensity_; }
700 {
return polymerViscosityCorrection_; }
704 {
return waterViscosityCorrection_; }
708 {
return *
static_cast<Implementation*
>(
this); }
720template <
class TypeTag>
733 {
throw std::logic_error(
"polymerMoleWeight() called but polymer molecular weight is disabled"); }
736 {
throw std::runtime_error(
"polymerConcentration() called but polymers are disabled"); }
739 {
throw std::runtime_error(
"polymerDeadPoreVolume() called but polymers are disabled"); }
742 {
throw std::runtime_error(
"polymerAdsorption() called but polymers are disabled"); }
745 {
throw std::runtime_error(
"polymerRockDensity() called but polymers are disabled"); }
748 {
throw std::runtime_error(
"polymerViscosityCorrection() called but polymers are disabled"); }
751 {
throw std::runtime_error(
"waterViscosityCorrection() called but polymers are disabled"); }
754template <
class TypeTag,
bool enablePolymerV>
764template <
class TypeTag>
775 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
786 template <
class Dummy =
bool>
792 throw std::runtime_error(
"The extension of the blackoil model for polymers is not yet "
793 "implemented for problems specified using permeabilities.");
802 template <
class Dummy =
bool>
808 waterShearFactor_ = 1.0;
809 polymerShearFactor_ = 1.0;
811 if (!PolymerModule::hasPlyshlog()) {
815 const ExtensiveQuantities& extQuants = asImp_();
816 const unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
817 const unsigned interiorDofIdx = extQuants.interiorIndex();
818 const unsigned exteriorDofIdx = extQuants.exteriorIndex();
819 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
820 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
821 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
824 const Evaluation poroAvg = intQuantsIn.porosity() * 0.5 + intQuantsEx.porosity() * 0.5;
825 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvfIdx, timeIdx);
826 const Evaluation& Sw = up.fluidState().saturation(waterPhaseIdx);
827 const unsigned cellIdx = elemCtx.globalSpaceIndex(scvfIdx, timeIdx);
828 const auto& materialLawManager = elemCtx.problem().materialLawManager();
829 const auto& scaledDrainageInfo =
830 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
831 const Scalar& Swcr = scaledDrainageInfo.Swcr;
834 const Evaluation denom = max(poroAvg * (Sw - Swcr), 1e-12);
835 Evaluation waterVolumeVelocity = extQuants.volumeFlux(waterPhaseIdx) / denom;
838 if (PolymerModule::hasShrate()) {
839 const Evaluation& relWater = up.relativePermeability(waterPhaseIdx);
840 const Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
842 const Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
843 const auto dist = elemCtx.pos(interiorDofIdx, timeIdx) - elemCtx.pos(exteriorDofIdx, timeIdx);
845 const Scalar absPerm = trans / faceArea * dist.two_norm();
846 waterVolumeVelocity *=
847 PolymerModule::shrate(pvtnumRegionIdx) * sqrt(poroAvg * Sw / (relWater * absPerm));
848 assert(isfinite(waterVolumeVelocity));
854 PolymerModule::computeShearFactor(up.polymerConcentration(),
856 waterVolumeVelocity);
857 polymerShearFactor_ =
858 PolymerModule::computeShearFactor(up.polymerConcentration(),
860 waterVolumeVelocity * up.polymerViscosityCorrection());
864 {
return polymerShearFactor_; }
867 {
return waterShearFactor_; }
870 Implementation& asImp_()
871 {
return *
static_cast<Implementation*
>(
this); }
873 Evaluation polymerShearFactor_;
874 Evaluation waterShearFactor_;
877template <
class TypeTag>
895 {
throw std::runtime_error(
"polymerShearFactor() called but polymers are disabled"); }
898 {
throw std::runtime_error(
"waterShearFactor() called but polymers are disabled"); }
Contains the parameters required to extend the black-oil model by polymer.
Declares the properties required by the black oil model.
const Evaluation & polymerShearFactor() const
Definition: blackoilpolymermodules.hh:894
void updateShearMultipliersPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:889
void updateShearMultipliers(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:884
const Evaluation & waterShearFactor() const
Definition: blackoilpolymermodules.hh:897
const Evaluation & polymerShearFactor() const
Definition: blackoilpolymermodules.hh:863
void updateShearMultipliersPerm(const ElementContext &, unsigned, unsigned)
Method which calculates the shear factor based on flow velocity.
Definition: blackoilpolymermodules.hh:788
const Evaluation & waterShearFactor() const
Definition: blackoilpolymermodules.hh:866
void updateShearMultipliers(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the shear factor based on flow velocity.
Definition: blackoilpolymermodules.hh:804
Provides the polymer specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilpolymermodules.hh:755
const Evaluation & polymerConcentration() const
Definition: blackoilpolymermodules.hh:735
const Evaluation & waterViscosityCorrection() const
Definition: blackoilpolymermodules.hh:750
const Evaluation & polymerRockDensity() const
Definition: blackoilpolymermodules.hh:744
const Evaluation & polymerDeadPoreVolume() const
Definition: blackoilpolymermodules.hh:738
const Evaluation & polymerMoleWeight() const
Definition: blackoilpolymermodules.hh:732
void polymerPropertiesUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:727
const Evaluation & polymerAdsorption() const
Definition: blackoilpolymermodules.hh:741
const Evaluation & polymerViscosityCorrection() const
Definition: blackoilpolymermodules.hh:747
Evaluation polymerViscosityCorrection_
Definition: blackoilpolymermodules.hh:716
Scalar polymerDeadPoreVolume_
Definition: blackoilpolymermodules.hh:713
void polymerPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the intensive properties needed to handle polymers from the primary variables.
Definition: blackoilpolymermodules.hh:600
Implementation & asImp_()
Definition: blackoilpolymermodules.hh:707
Scalar polymerRockDensity_
Definition: blackoilpolymermodules.hh:714
const Evaluation & polymerMoleWeight() const
Definition: blackoilpolymermodules.hh:679
Evaluation polymerConcentration_
Definition: blackoilpolymermodules.hh:710
Scalar polymerDeadPoreVolume() const
Definition: blackoilpolymermodules.hh:689
Evaluation polymerMoleWeight_
Definition: blackoilpolymermodules.hh:712
const Evaluation & waterViscosityCorrection() const
Definition: blackoilpolymermodules.hh:703
const Evaluation & polymerAdsorption() const
Definition: blackoilpolymermodules.hh:692
const Evaluation & polymerConcentration() const
Definition: blackoilpolymermodules.hh:676
const Evaluation & polymerViscosityCorrection() const
Definition: blackoilpolymermodules.hh:699
Evaluation waterViscosityCorrection_
Definition: blackoilpolymermodules.hh:717
Scalar polymerRockDensity() const
Definition: blackoilpolymermodules.hh:695
Evaluation polymerAdsorption_
Definition: blackoilpolymermodules.hh:715
Provides the volumetric quantities required for the equations needed by the polymers extension of the...
Definition: blackoilpolymermodules.hh:565
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:64
static const TabulatedFunction & plyviscViscosityMultiplierTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:416
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:179
static const TabulatedFunction & plyadsAdsorbedPolymer(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:407
static TabulatedTwoDFunction & getPlymwinjTable(const int tableNumber)
get the PLYMWINJ table
Definition: blackoilpolymermodules.hh:103
static Scalar plyrockAdsorbtionIndex(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:390
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:339
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:226
static void registerParameters()
Register all run-time parameters for the black-oil polymer module.
Definition: blackoilpolymermodules.hh:146
static bool eqApplies(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:199
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:191
static BlackOilPolymerParams< Scalar >::SkprpolyTable & getSkprpolyTable(const int tableNumber)
get the SKPRPOLY table
Definition: blackoilpolymermodules.hh:132
static bool hasShrate()
Definition: blackoilpolymermodules.hh:460
static Scalar plyrockRockDensityFactor(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:382
static Scalar plyrockDeadPoreVolume(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:366
static const TabulatedFunction & plyviscViscosityMultiplierTable(unsigned pvtnumRegionIdx)
Definition: blackoilpolymermodules.hh:426
static void registerOutputModules(Model &model, Simulator &simulator)
Register all polymer specific VTK and ECL output modules.
Definition: blackoilpolymermodules.hh:156
static Scalar plymaxMaxConcentration(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:429
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilpolymermodules.hh:329
Scalar molarMass() const
Definition: blackoilpolymermodules.hh:553
static const BlackOilPolymerParams< Scalar >::PlyvmhCoefficients & plyvmhCoefficients(const ElementContext &elemCtx, const unsigned scvIdx, const unsigned timeIdx)
Definition: blackoilpolymermodules.hh:448
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:164
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilpolymermodules.hh:236
static bool hasPlyshlog()
Definition: blackoilpolymermodules.hh:457
static Scalar shrate(unsigned pvtnumRegionIdx)
Definition: blackoilpolymermodules.hh:463
static std::string eqName(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:214
static void setParams(BlackOilPolymerParams< Scalar > &¶ms)
Set parameters.
Definition: blackoilpolymermodules.hh:95
static TabulatedTwoDFunction & getSkprwatTable(const int tableNumber)
get the SKPRWAT table
Definition: blackoilpolymermodules.hh:117
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:350
static Scalar plymixparToddLongstaff(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:438
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:275
static Scalar plyrockMaxAdsorbtion(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:398
static Scalar plyrockResidualResistanceFactor(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:374
static Evaluation computeShearFactor(const Evaluation &polymerConcentration, unsigned pvtnumRegionIdx, const Evaluation &v0)
Computes the shear factor.
Definition: blackoilpolymermodules.hh:473
VTK output module for the black oil model's polymer related quantities.
Definition: vtkblackoilpolymermodule.hpp:61
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilpolymermodule.hpp:91
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
The Opm property system, traits with inheritance.
Definition: blackoilpolymerparams.hpp:85
Definition: blackoilpolymerparams.hpp:92
Struct holding the parameters for the BlackOilPolymerModule class.
Definition: blackoilpolymerparams.hpp:45
IntervalTabulated2DFunction< Scalar > TabulatedTwoDFunction
Definition: blackoilpolymerparams.hpp:47
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilpolymerparams.hpp:46