blackoilpolymermodules.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_POLYMER_MODULE_HH
29#define EWOMS_BLACK_OIL_POLYMER_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/common/OpmLog/OpmLog.hpp>
34#include <opm/common/utility/gpuDecorators.hpp>
35
36#include <opm/material/common/MathToolbox.hpp>
37
40
42
44
45#include <algorithm>
46#include <cassert>
47#include <cmath>
48#include <cstddef>
49#include <istream>
50#include <memory>
51#include <ostream>
52#include <stdexcept>
53#include <string>
54#include <vector>
55
56namespace Opm {
57
63template <class TypeTag, bool enablePolymerV = getPropValue<TypeTag, Properties::EnablePolymer>()>
65{
77
78 using Toolbox = MathToolbox<Evaluation>;
79
80 using TabulatedFunction = typename BlackOilPolymerParams<Scalar>::TabulatedFunction;
81 using TabulatedTwoDFunction = typename BlackOilPolymerParams<Scalar>::TabulatedTwoDFunction;
82
83 static constexpr unsigned polymerConcentrationIdx = Indices::polymerConcentrationIdx;
84 static constexpr unsigned polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
85 static constexpr unsigned contiPolymerEqIdx = Indices::contiPolymerEqIdx;
86 static constexpr unsigned contiPolymerMolarWeightEqIdx = Indices::contiPolymerMWEqIdx;
87 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
88
89 static constexpr unsigned enablePolymer = enablePolymerV;
90 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
91
92 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
93
94public:
97 {
98 params_ = params;
99 }
100
104 static TabulatedTwoDFunction& getPlymwinjTable(const int tableNumber)
105 {
106 const auto iterTable = params_.plymwinjTables_.find(tableNumber);
107 if (iterTable != params_.plymwinjTables_.end()) {
108 return iterTable->second;
109 }
110 else {
111 throw std::runtime_error(" the PLYMWINJ table " + std::to_string(tableNumber) + " does not exist\n");
112 }
113 }
114
118 static TabulatedTwoDFunction& getSkprwatTable(const int tableNumber)
119 {
120 const auto iterTable = params_.skprwatTables_.find(tableNumber);
121 if (iterTable != params_.skprwatTables_.end()) {
122 return iterTable->second;
123 }
124 else {
125 throw std::runtime_error(" the SKPRWAT table " + std::to_string(tableNumber) + " does not exist\n");
126 }
127 }
128
133 getSkprpolyTable(const int tableNumber)
134 {
135 const auto iterTable = params_.skprpolyTables_.find(tableNumber);
136 if (iterTable != params_.skprpolyTables_.end()) {
137 return iterTable->second;
138 }
139 else {
140 throw std::runtime_error(" the SKPRPOLY table " + std::to_string(tableNumber) + " does not exist\n");
141 }
142 }
143
147 static void registerParameters()
148 {
149 if constexpr (enablePolymer) {
151 }
152 }
153
157 static void registerOutputModules(Model& model,
158 Simulator& simulator)
159 {
160 if constexpr (enablePolymer) {
161 model.addOutputModule(std::make_unique<VtkBlackOilPolymerModule<TypeTag>>(simulator));
162 }
163 }
164
165 static bool primaryVarApplies(unsigned pvIdx)
166 {
167 if constexpr (enablePolymer) {
168 if constexpr (enablePolymerMolarWeight) {
169 return pvIdx == polymerConcentrationIdx || pvIdx == polymerMoleWeightIdx;
170 }
171 else {
172 return pvIdx == polymerConcentrationIdx;
173 }
174 }
175 else {
176 return false;
177 }
178 }
179
180 static std::string primaryVarName(unsigned pvIdx)
181 {
182 assert(primaryVarApplies(pvIdx));
183
184 if (pvIdx == polymerConcentrationIdx) {
185 return "polymer_waterconcentration";
186 }
187 else {
188 return "polymer_molecularweight";
189 }
190 }
191
192 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
193 {
194 assert(primaryVarApplies(pvIdx));
195
196 // TODO: it may be beneficial to chose this differently.
197 return static_cast<Scalar>(1.0);
198 }
199
200 static bool eqApplies(unsigned eqIdx)
201 {
202 if constexpr (enablePolymer) {
203 if constexpr (enablePolymerMolarWeight) {
204 return eqIdx == contiPolymerEqIdx || eqIdx == contiPolymerMolarWeightEqIdx;
205 }
206 else {
207 return eqIdx == contiPolymerEqIdx;
208 }
209 }
210 else {
211 return false;
212 }
213 }
214
215 static std::string eqName(unsigned eqIdx)
216 {
217 assert(eqApplies(eqIdx));
218
219 if (eqIdx == contiPolymerEqIdx) {
220 return "conti^polymer";
221 }
222 else {
223 return "conti^polymer_molecularweight";
224 }
225 }
226
227 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
228 {
229 assert(eqApplies(eqIdx));
230
231 // TODO: it may be beneficial to chose this differently.
232 return static_cast<Scalar>(1.0);
233 }
234
235 // must be called after water storage is computed
236 template <class StorageType>
237 OPM_HOST_DEVICE static void addStorage(StorageType& storage,
238 const IntensiveQuantities& intQuants)
239 {
240 using LhsEval = typename StorageType::value_type;
241
242 if constexpr (enablePolymer) {
243 const auto& fs = intQuants.fluidState();
244
245 // avoid singular matrix if no water is present.
246 const LhsEval surfaceVolumeWater =
247 max(Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx)) *
248 Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx)) *
249 Toolbox::template decay<LhsEval>(intQuants.porosity()),
250 1e-10);
251
252 // polymer in water phase
253 const LhsEval massPolymer =
254 surfaceVolumeWater *
255 Toolbox::template decay<LhsEval>(intQuants.polymerConcentration()) *
256 (1.0 - Toolbox::template decay<LhsEval>(intQuants.polymerDeadPoreVolume()));
257
258 // polymer in solid phase
259 const LhsEval adsorptionPolymer =
260 Toolbox::template decay<LhsEval>(1.0 - intQuants.porosity()) *
261 Toolbox::template decay<LhsEval>(intQuants.polymerRockDensity()) *
262 Toolbox::template decay<LhsEval>(intQuants.polymerAdsorption());
263
264 LhsEval accumulationPolymer = massPolymer + adsorptionPolymer;
265
266 storage[contiPolymerEqIdx] += accumulationPolymer;
267
268 // tracking the polymer molecular weight
269 if constexpr (enablePolymerMolarWeight) {
270 accumulationPolymer = max(accumulationPolymer, 1e-10);
271
272 storage[contiPolymerMolarWeightEqIdx] +=
273 accumulationPolymer * Toolbox::template decay<LhsEval>(intQuants.polymerMoleWeight());
274 }
275 }
276 }
277
278 static void computeFlux([[maybe_unused]] RateVector& flux,
279 [[maybe_unused]] const ElementContext& elemCtx,
280 [[maybe_unused]] unsigned scvfIdx,
281 [[maybe_unused]] unsigned timeIdx)
282 {
283 if constexpr (enablePolymer) {
284 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
285
286 const unsigned upIdx = extQuants.upstreamIndex(FluidSystem::waterPhaseIdx);
287 const unsigned inIdx = extQuants.interiorIndex();
288 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
289 const unsigned contiWaterEqIdx =
290 Indices::conti0EqIdx + FluidSystem::canonicalToActiveCompIdx(FluidSystem::waterCompIdx);
291
292 if (upIdx == inIdx) {
293 flux[contiPolymerEqIdx] =
294 extQuants.volumeFlux(waterPhaseIdx) *
295 up.fluidState().invB(waterPhaseIdx) *
296 up.polymerViscosityCorrection() /
297 extQuants.polymerShearFactor() *
298 up.polymerConcentration();
299
300 // modify water
301 flux[contiWaterEqIdx] /= extQuants.waterShearFactor();
302 }
303 else {
304 flux[contiPolymerEqIdx] =
305 extQuants.volumeFlux(waterPhaseIdx) *
306 decay<Scalar>(up.fluidState().invB(waterPhaseIdx)) *
307 decay<Scalar>(up.polymerViscosityCorrection()) /
308 decay<Scalar>(extQuants.polymerShearFactor()) *
309 decay<Scalar>(up.polymerConcentration());
310
311 // modify water
312 flux[contiWaterEqIdx] /= decay<Scalar>(extQuants.waterShearFactor());
313 }
314
315 // flux related to transport of polymer molecular weight
316 if constexpr (enablePolymerMolarWeight) {
317 if (upIdx == inIdx) {
318 flux[contiPolymerMolarWeightEqIdx] =
319 flux[contiPolymerEqIdx] * up.polymerMoleWeight();
320 }
321 else {
322 flux[contiPolymerMolarWeightEqIdx] =
323 flux[contiPolymerEqIdx] * decay<Scalar>(up.polymerMoleWeight());
324 }
325 }
326 }
327 }
328
332 static Scalar computeUpdateError(const PrimaryVariables&,
333 const EqVector&)
334 {
335 // do not consider consider the cange of polymer primary variables for
336 // convergence
337 // TODO: maybe this should be changed
338 return static_cast<Scalar>(0.0);
339 }
340
341 template <class DofEntity>
342 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
343 {
344 if constexpr (enablePolymer) {
345 const unsigned dofIdx = model.dofMapper().index(dof);
346 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
347 outstream << priVars[polymerConcentrationIdx];
348 outstream << priVars[polymerMoleWeightIdx];
349 }
350 }
351
352 template <class DofEntity>
353 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
354 {
355 if constexpr (enablePolymer) {
356 const unsigned dofIdx = model.dofMapper().index(dof);
357 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
358 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
359
360 instream >> priVars0[polymerConcentrationIdx];
361 instream >> priVars0[polymerMoleWeightIdx];
362
363 // set the primary variables for the beginning of the current time step.
364 priVars1[polymerConcentrationIdx] = priVars0[polymerConcentrationIdx];
365 priVars1[polymerMoleWeightIdx] = priVars0[polymerMoleWeightIdx];
366 }
367 }
368
369 static Scalar plyrockDeadPoreVolume(const ElementContext& elemCtx,
370 unsigned scvIdx,
371 unsigned timeIdx)
372 {
373 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
374 return params_.plyrockDeadPoreVolume_[satnumRegionIdx];
375 }
376
377 static Scalar plyrockResidualResistanceFactor(const ElementContext& elemCtx,
378 unsigned scvIdx,
379 unsigned timeIdx)
380 {
381 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
382 return params_.plyrockResidualResistanceFactor_[satnumRegionIdx];
383 }
384
385 static Scalar plyrockRockDensityFactor(const ElementContext& elemCtx,
386 unsigned scvIdx,
387 unsigned timeIdx)
388 {
389 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
390 return params_.plyrockRockDensityFactor_[satnumRegionIdx];
391 }
392
393 static Scalar plyrockAdsorbtionIndex(const ElementContext& elemCtx,
394 unsigned scvIdx,
395 unsigned timeIdx)
396 {
397 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
398 return params_.plyrockAdsorbtionIndex_[satnumRegionIdx];
399 }
400
401 static Scalar plyrockMaxAdsorbtion(const ElementContext& elemCtx,
402 unsigned scvIdx,
403 unsigned timeIdx)
404 {
405 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
406 return params_.plyrockMaxAdsorbtion_[satnumRegionIdx];
407 }
408
409 static const TabulatedFunction&
410 plyadsAdsorbedPolymer(const ElementContext& elemCtx,
411 unsigned scvIdx,
412 unsigned timeIdx)
413 {
414 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
415 return params_.plyadsAdsorbedPolymer_[satnumRegionIdx];
416 }
417
418 static const TabulatedFunction&
419 plyviscViscosityMultiplierTable(const ElementContext& elemCtx,
420 unsigned scvIdx,
421 unsigned timeIdx)
422 {
423 unsigned pvtnumRegionIdx =
424 elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
425 return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
426 }
427
428 static const TabulatedFunction&
429 plyviscViscosityMultiplierTable(unsigned pvtnumRegionIdx)
430 { return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx]; }
431
432 static Scalar plymaxMaxConcentration(const ElementContext& elemCtx,
433 unsigned scvIdx,
434 unsigned timeIdx)
435 {
436 const unsigned polymerMixRegionIdx =
437 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
438 return params_.plymaxMaxConcentration_[polymerMixRegionIdx];
439 }
440
441 static Scalar plymixparToddLongstaff(const ElementContext& elemCtx,
442 unsigned scvIdx,
443 unsigned timeIdx)
444 {
445 const unsigned polymerMixRegionIdx =
446 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
447 return params_.plymixparToddLongstaff_[polymerMixRegionIdx];
448 }
449
451 plyvmhCoefficients(const ElementContext& elemCtx,
452 const unsigned scvIdx,
453 const unsigned timeIdx)
454 {
455 const unsigned polymerMixRegionIdx =
456 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
457 return params_.plyvmhCoefficients_[polymerMixRegionIdx];
458 }
459
460 static bool hasPlyshlog()
461 { return params_.hasPlyshlog_; }
462
463 static bool hasShrate()
464 { return params_.hasShrate_; }
465
466 static Scalar shrate(unsigned pvtnumRegionIdx)
467 { return params_.shrate_[pvtnumRegionIdx]; }
468
475 template <class Evaluation>
476 static Evaluation computeShearFactor(const Evaluation& polymerConcentration,
477 unsigned pvtnumRegionIdx,
478 const Evaluation& v0)
479 {
480 using ToolboxLocal = MathToolbox<Evaluation>;
481
482 const auto& viscosityMultiplierTable = params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
483 const Scalar viscosityMultiplier =
484 viscosityMultiplierTable.eval(scalarValue(polymerConcentration), /*extrapolate=*/true);
485
486 const Scalar eps = 1e-14;
487 // return 1.0 if the polymer has no effect on the water.
488 if (std::abs((viscosityMultiplier - 1.0)) < eps) {
489 return ToolboxLocal::createConstant(v0, 1.0);
490 }
491
492 const std::vector<Scalar>& shearEffectRefLogVelocity =
493 params_.plyshlogShearEffectRefLogVelocity_[pvtnumRegionIdx];
494 const auto v0AbsLog = log(abs(v0));
495 // return 1.0 if the velocity /sharte is smaller than the first velocity entry.
496 if (v0AbsLog < shearEffectRefLogVelocity[0]) {
497 return ToolboxLocal::createConstant(v0, 1.0);
498 }
499
500 // compute shear factor from input
501 // Z = (1 + (P - 1) * M(v)) / P
502 // where M(v) is computed from user input
503 // and P = viscosityMultiplier
504 const std::vector<Scalar>& shearEffectRefMultiplier =
505 params_.plyshlogShearEffectRefMultiplier_[pvtnumRegionIdx];
506 const std::size_t numTableEntries = shearEffectRefLogVelocity.size();
507 assert(shearEffectRefMultiplier.size() == numTableEntries);
508
509 std::vector<Scalar> shearEffectMultiplier(numTableEntries, 1.0);
510 for (std::size_t i = 0; i < numTableEntries; ++i) {
511 shearEffectMultiplier[i] = (1.0 + (viscosityMultiplier - 1.0) *
512 shearEffectRefMultiplier[i]) / viscosityMultiplier;
513 shearEffectMultiplier[i] = log(shearEffectMultiplier[i]);
514 }
515 // store the logarithmic velocity and logarithmic multipliers in a table for easy look up and
516 // linear interpolation in the logarithmic space.
517 const TabulatedFunction logShearEffectMultiplier =
518 TabulatedFunction(numTableEntries, shearEffectRefLogVelocity,
519 shearEffectMultiplier, /*bool sortInputs =*/ false);
520
521 // Find sheared velocity (v) that satisfies
522 // F = log(v) + log (Z) - log(v0) = 0;
523
524 // Set up the function
525 // u = log(v)
526 auto F = [&logShearEffectMultiplier, &v0AbsLog](const Evaluation& u) {
527 return u + logShearEffectMultiplier.eval(u, true) - v0AbsLog;
528 };
529 // and its derivative
530 auto dF = [&logShearEffectMultiplier](const Evaluation& u) {
531 return 1 + logShearEffectMultiplier.evalDerivative(u, true);
532 };
533
534 // Solve F = 0 using Newton
535 // Use log(v0) as initial value for u
536 auto u = v0AbsLog;
537 bool converged = false;
538 // TODO make this into parameters
539 for (int i = 0; i < 20; ++i) {
540 const auto f = F(u);
541 const auto df = dF(u);
542 u -= f / df;
543 if (std::abs(scalarValue(f)) < 1e-12) {
544 converged = true;
545 break;
546 }
547 }
548 if (!converged) {
549 throw std::runtime_error("Not able to compute shear velocity. \n");
550 }
551
552 // return the shear factor
553 return exp(logShearEffectMultiplier.eval(u, /*extrapolate=*/true));
554 }
555
556 Scalar molarMass() const
557 { return 0.25; } // kg/mol
558
559private:
560 static BlackOilPolymerParams<Scalar> params_;
561};
562
563template <class TypeTag, bool enablePolymerV>
564BlackOilPolymerParams<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
565BlackOilPolymerModule<TypeTag, enablePolymerV>::params_;
566
567template <class TypeTag, bool enablePolymerV>
569
577template <class TypeTag>
578class BlackOilPolymerIntensiveQuantities<TypeTag, /*enablePolymerV=*/true>
579{
581
589
591
592 static constexpr int polymerConcentrationIdx = Indices::polymerConcentrationIdx;
593 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
594 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
595 static constexpr int polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
596
597public:
603 void polymerPropertiesUpdate_(const ElementContext& elemCtx,
604 unsigned dofIdx,
605 unsigned timeIdx)
606 {
607 const auto linearizationType = elemCtx.linearizationType();
608 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
609 polymerConcentration_ = priVars.makeEvaluation(polymerConcentrationIdx, timeIdx, linearizationType);
610 if constexpr (enablePolymerMolarWeight) {
611 polymerMoleWeight_ = priVars.makeEvaluation(polymerMoleWeightIdx, timeIdx, linearizationType);
612 }
613
614 // permeability reduction due to polymer
615 const Scalar& maxAdsorbtion = PolymerModule::plyrockMaxAdsorbtion(elemCtx, dofIdx, timeIdx);
616 const auto& plyadsAdsorbedPolymer = PolymerModule::plyadsAdsorbedPolymer(elemCtx, dofIdx, timeIdx);
617 polymerAdsorption_ = plyadsAdsorbedPolymer.eval(polymerConcentration_, /*extrapolate=*/true);
618 if (static_cast<int>(PolymerModule::plyrockAdsorbtionIndex(elemCtx, dofIdx, timeIdx)) ==
620 {
621 const auto maxPolymerAdsorption =
622 elemCtx.problem().maxPolymerAdsorption(elemCtx, dofIdx, timeIdx);
623 polymerAdsorption_ = std::max(Evaluation(maxPolymerAdsorption), polymerAdsorption_);
624 }
625
626 // compute resitanceFactor
627 const Scalar& residualResistanceFactor =
628 PolymerModule::plyrockResidualResistanceFactor(elemCtx, dofIdx, timeIdx);
629 const Evaluation resistanceFactor = 1.0 + (residualResistanceFactor - 1.0) *
630 polymerAdsorption_ / maxAdsorbtion;
631
632 // compute effective viscosities
633 if constexpr (!enablePolymerMolarWeight) {
634 const Scalar cmax = PolymerModule::plymaxMaxConcentration(elemCtx, dofIdx, timeIdx);
635 const auto& fs = asImp_().fluidState_;
636 const Evaluation& muWater = fs.viscosity(waterPhaseIdx);
637 const auto& viscosityMultiplier =
638 PolymerModule::plyviscViscosityMultiplierTable(elemCtx, dofIdx, timeIdx);
639 const Evaluation viscosityMixture =
640 viscosityMultiplier.eval(polymerConcentration_, /*extrapolate=*/true) * muWater;
641
642 // Do the Todd-Longstaff mixing
643 const Scalar plymixparToddLongstaff = PolymerModule::plymixparToddLongstaff(elemCtx, dofIdx, timeIdx);
644 const Evaluation viscosityPolymer = viscosityMultiplier.eval(cmax, /*extrapolate=*/true) * muWater;
645 const Evaluation viscosityPolymerEffective =
646 pow(viscosityMixture, plymixparToddLongstaff) * pow(viscosityPolymer, 1.0 - plymixparToddLongstaff);
647 const Evaluation viscosityWaterEffective =
648 pow(viscosityMixture, plymixparToddLongstaff) * pow(muWater, 1.0 - plymixparToddLongstaff);
649
650 const Evaluation cbar = polymerConcentration_ / cmax;
651 // waterViscosity / effectiveWaterViscosity
652 waterViscosityCorrection_ = muWater * ((1.0 - cbar) / viscosityWaterEffective + cbar / viscosityPolymerEffective);
653 // effectiveWaterViscosity / effectivePolymerViscosity
654 polymerViscosityCorrection_ = (muWater / waterViscosityCorrection_) / viscosityPolymerEffective;
655 }
656 else { // based on PLYVMH
657 const auto& plyvmhCoefficients = PolymerModule::plyvmhCoefficients(elemCtx, dofIdx, timeIdx);
658 const Scalar k_mh = plyvmhCoefficients.k_mh;
659 const Scalar a_mh = plyvmhCoefficients.a_mh;
660 const Scalar gamma = plyvmhCoefficients.gamma;
661 const Scalar kappa = plyvmhCoefficients.kappa;
662
663 // viscosity model based on Mark-Houwink equation and Huggins equation
664 // 1000 is a emperical constant, most likely related to unit conversion
665 const Evaluation intrinsicViscosity = k_mh * pow(polymerMoleWeight_ * 1000., a_mh);
666 const Evaluation x = polymerConcentration_ * intrinsicViscosity;
667 waterViscosityCorrection_ = 1.0 / (1.0 + gamma * (x + kappa * x * x));
668 polymerViscosityCorrection_ = 1.0;
669 }
670
671 // adjust water mobility
672 asImp_().mobility_[waterPhaseIdx] *= waterViscosityCorrection_ / resistanceFactor;
673
674 // update rock properties
675 polymerDeadPoreVolume_ = PolymerModule::plyrockDeadPoreVolume(elemCtx, dofIdx, timeIdx);
676 polymerRockDensity_ = PolymerModule::plyrockRockDensityFactor(elemCtx, dofIdx, timeIdx);
677 }
678
679 const Evaluation& polymerConcentration() const
680 { return polymerConcentration_; }
681
682 const Evaluation& polymerMoleWeight() const
683 {
684 if constexpr (enablePolymerMolarWeight) {
685 return polymerMoleWeight_;
686 }
687 else {
688 throw std::logic_error("polymerMoleWeight() is called but polymer milecular weight is disabled");
689 }
690 }
691
693 { return polymerDeadPoreVolume_; }
694
695 const Evaluation& polymerAdsorption() const
696 { return polymerAdsorption_; }
697
698 Scalar polymerRockDensity() const
699 { return polymerRockDensity_; }
700
701 // effectiveWaterViscosity / effectivePolymerViscosity
702 const Evaluation& polymerViscosityCorrection() const
703 { return polymerViscosityCorrection_; }
704
705 // waterViscosity / effectiveWaterViscosity
706 const Evaluation& waterViscosityCorrection() const
707 { return waterViscosityCorrection_; }
708
709protected:
710 Implementation& asImp_()
711 { return *static_cast<Implementation*>(this); }
712
714 // polymer molecular weight
721};
722
723template <class TypeTag>
725{
728
729public:
730 void polymerPropertiesUpdate_(const ElementContext&,
731 unsigned,
732 unsigned)
733 {}
734
735 const Evaluation& polymerMoleWeight() const
736 { throw std::logic_error("polymerMoleWeight() called but polymer molecular weight is disabled"); }
737
738 const Evaluation& polymerConcentration() const
739 { throw std::runtime_error("polymerConcentration() called but polymers are disabled"); }
740
741 const Evaluation& polymerDeadPoreVolume() const
742 { throw std::runtime_error("polymerDeadPoreVolume() called but polymers are disabled"); }
743
744 const Evaluation& polymerAdsorption() const
745 { throw std::runtime_error("polymerAdsorption() called but polymers are disabled"); }
746
747 const Evaluation& polymerRockDensity() const
748 { throw std::runtime_error("polymerRockDensity() called but polymers are disabled"); }
749
750 const Evaluation& polymerViscosityCorrection() const
751 { throw std::runtime_error("polymerViscosityCorrection() called but polymers are disabled"); }
752
753 const Evaluation& waterViscosityCorrection() const
754 { throw std::runtime_error("waterViscosityCorrection() called but polymers are disabled"); }
755};
756
757template <class TypeTag, bool enablePolymerV>
759
767template <class TypeTag>
768class BlackOilPolymerExtensiveQuantities<TypeTag, /*enablePolymerV=*/true>
769{
771
777
778 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
779
781
782public:
789 template <class Dummy = bool> // we need to make this method a template to avoid
790 // compiler errors if it is not instantiated!
791 void updateShearMultipliersPerm(const ElementContext&,
792 unsigned,
793 unsigned)
794 {
795 throw std::runtime_error("The extension of the blackoil model for polymers is not yet "
796 "implemented for problems specified using permeabilities.");
797 }
798
805 template <class Dummy = bool> // we need to make this method a template to avoid
806 // compiler errors if it is not instantiated!
807 void updateShearMultipliers(const ElementContext& elemCtx,
808 unsigned scvfIdx,
809 unsigned timeIdx)
810 {
811 waterShearFactor_ = 1.0;
812 polymerShearFactor_ = 1.0;
813
814 if (!PolymerModule::hasPlyshlog()) {
815 return;
816 }
817
818 const ExtensiveQuantities& extQuants = asImp_();
819 const unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
820 const unsigned interiorDofIdx = extQuants.interiorIndex();
821 const unsigned exteriorDofIdx = extQuants.exteriorIndex();
822 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
823 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
824 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
825
826 // compute water velocity from flux
827 const Evaluation poroAvg = intQuantsIn.porosity() * 0.5 + intQuantsEx.porosity() * 0.5;
828 const unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvfIdx, timeIdx);
829 const Evaluation& Sw = up.fluidState().saturation(waterPhaseIdx);
830 const unsigned cellIdx = elemCtx.globalSpaceIndex(scvfIdx, timeIdx);
831 const auto& materialLawManager = elemCtx.problem().materialLawManager();
832 const auto& scaledDrainageInfo =
833 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
834 const Scalar& Swcr = scaledDrainageInfo.Swcr;
835
836 // guard against zero porosity and no mobile water
837 const Evaluation denom = max(poroAvg * (Sw - Swcr), 1e-12);
838 Evaluation waterVolumeVelocity = extQuants.volumeFlux(waterPhaseIdx) / denom;
839
840 // if shrate is specified. Compute shrate based on the water velocity
841 if (PolymerModule::hasShrate()) {
842 const Evaluation& relWater = up.relativePermeability(waterPhaseIdx);
843 const Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
844 if (trans > 0.0) {
845 const Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
846 const auto dist = elemCtx.pos(interiorDofIdx, timeIdx) - elemCtx.pos(exteriorDofIdx, timeIdx);
847 // compute permeability from transmissibility.
848 const Scalar absPerm = trans / faceArea * dist.two_norm();
849 waterVolumeVelocity *=
850 PolymerModule::shrate(pvtnumRegionIdx) * sqrt(poroAvg * Sw / (relWater * absPerm));
851 assert(isfinite(waterVolumeVelocity));
852 }
853 }
854
855 // compute share factors for water and polymer
856 waterShearFactor_ =
857 PolymerModule::computeShearFactor(up.polymerConcentration(),
858 pvtnumRegionIdx,
859 waterVolumeVelocity);
860 polymerShearFactor_ =
861 PolymerModule::computeShearFactor(up.polymerConcentration(),
862 pvtnumRegionIdx,
863 waterVolumeVelocity * up.polymerViscosityCorrection());
864 }
865
866 const Evaluation& polymerShearFactor() const
867 { return polymerShearFactor_; }
868
869 const Evaluation& waterShearFactor() const
870 { return waterShearFactor_; }
871
872private:
873 Implementation& asImp_()
874 { return *static_cast<Implementation*>(this); }
875
876 Evaluation polymerShearFactor_;
877 Evaluation waterShearFactor_;
878};
879
880template <class TypeTag>
882{
885
886public:
887 void updateShearMultipliers(const ElementContext&,
888 unsigned,
889 unsigned)
890 {}
891
892 void updateShearMultipliersPerm(const ElementContext&,
893 unsigned,
894 unsigned)
895 {}
896
897 const Evaluation& polymerShearFactor() const
898 { throw std::runtime_error("polymerShearFactor() called but polymers are disabled"); }
899
900 const Evaluation& waterShearFactor() const
901 { throw std::runtime_error("waterShearFactor() called but polymers are disabled"); }
902};
903
904} // namespace Opm
905
906#endif
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:897
void updateShearMultipliersPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:892
void updateShearMultipliers(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:887
const Evaluation & waterShearFactor() const
Definition: blackoilpolymermodules.hh:900
const Evaluation & polymerShearFactor() const
Definition: blackoilpolymermodules.hh:866
void updateShearMultipliersPerm(const ElementContext &, unsigned, unsigned)
Method which calculates the shear factor based on flow velocity.
Definition: blackoilpolymermodules.hh:791
const Evaluation & waterShearFactor() const
Definition: blackoilpolymermodules.hh:869
void updateShearMultipliers(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the shear factor based on flow velocity.
Definition: blackoilpolymermodules.hh:807
Provides the polymer specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilpolymermodules.hh:758
const Evaluation & polymerConcentration() const
Definition: blackoilpolymermodules.hh:738
const Evaluation & waterViscosityCorrection() const
Definition: blackoilpolymermodules.hh:753
const Evaluation & polymerRockDensity() const
Definition: blackoilpolymermodules.hh:747
const Evaluation & polymerDeadPoreVolume() const
Definition: blackoilpolymermodules.hh:741
const Evaluation & polymerMoleWeight() const
Definition: blackoilpolymermodules.hh:735
void polymerPropertiesUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilpolymermodules.hh:730
const Evaluation & polymerAdsorption() const
Definition: blackoilpolymermodules.hh:744
const Evaluation & polymerViscosityCorrection() const
Definition: blackoilpolymermodules.hh:750
Evaluation polymerViscosityCorrection_
Definition: blackoilpolymermodules.hh:719
Scalar polymerDeadPoreVolume_
Definition: blackoilpolymermodules.hh:716
void polymerPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the intensive properties needed to handle polymers from the primary variables.
Definition: blackoilpolymermodules.hh:603
Implementation & asImp_()
Definition: blackoilpolymermodules.hh:710
Scalar polymerRockDensity_
Definition: blackoilpolymermodules.hh:717
const Evaluation & polymerMoleWeight() const
Definition: blackoilpolymermodules.hh:682
Evaluation polymerConcentration_
Definition: blackoilpolymermodules.hh:713
Scalar polymerDeadPoreVolume() const
Definition: blackoilpolymermodules.hh:692
Evaluation polymerMoleWeight_
Definition: blackoilpolymermodules.hh:715
const Evaluation & waterViscosityCorrection() const
Definition: blackoilpolymermodules.hh:706
const Evaluation & polymerAdsorption() const
Definition: blackoilpolymermodules.hh:695
const Evaluation & polymerConcentration() const
Definition: blackoilpolymermodules.hh:679
const Evaluation & polymerViscosityCorrection() const
Definition: blackoilpolymermodules.hh:702
Evaluation waterViscosityCorrection_
Definition: blackoilpolymermodules.hh:720
Scalar polymerRockDensity() const
Definition: blackoilpolymermodules.hh:698
Evaluation polymerAdsorption_
Definition: blackoilpolymermodules.hh:718
Provides the volumetric quantities required for the equations needed by the polymers extension of the...
Definition: blackoilpolymermodules.hh:568
Contains the high level supplements required to extend the black oil model by polymer.
Definition: blackoilpolymermodules.hh:65
static const TabulatedFunction & plyviscViscosityMultiplierTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:419
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:180
static const TabulatedFunction & plyadsAdsorbedPolymer(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:410
static TabulatedTwoDFunction & getPlymwinjTable(const int tableNumber)
get the PLYMWINJ table
Definition: blackoilpolymermodules.hh:104
static Scalar plyrockAdsorbtionIndex(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:393
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:342
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:227
static void registerParameters()
Register all run-time parameters for the black-oil polymer module.
Definition: blackoilpolymermodules.hh:147
static bool eqApplies(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:200
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:192
static BlackOilPolymerParams< Scalar >::SkprpolyTable & getSkprpolyTable(const int tableNumber)
get the SKPRPOLY table
Definition: blackoilpolymermodules.hh:133
static bool hasShrate()
Definition: blackoilpolymermodules.hh:463
static Scalar plyrockRockDensityFactor(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:385
static Scalar plyrockDeadPoreVolume(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:369
static const TabulatedFunction & plyviscViscosityMultiplierTable(unsigned pvtnumRegionIdx)
Definition: blackoilpolymermodules.hh:429
static void registerOutputModules(Model &model, Simulator &simulator)
Register all polymer specific VTK and ECL output modules.
Definition: blackoilpolymermodules.hh:157
static Scalar plymaxMaxConcentration(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:432
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilpolymermodules.hh:332
Scalar molarMass() const
Definition: blackoilpolymermodules.hh:556
static const BlackOilPolymerParams< Scalar >::PlyvmhCoefficients & plyvmhCoefficients(const ElementContext &elemCtx, const unsigned scvIdx, const unsigned timeIdx)
Definition: blackoilpolymermodules.hh:451
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:165
static bool hasPlyshlog()
Definition: blackoilpolymermodules.hh:460
static Scalar shrate(unsigned pvtnumRegionIdx)
Definition: blackoilpolymermodules.hh:466
static std::string eqName(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:215
static OPM_HOST_DEVICE void addStorage(StorageType &storage, const IntensiveQuantities &intQuants)
Definition: blackoilpolymermodules.hh:237
static void setParams(BlackOilPolymerParams< Scalar > &&params)
Set parameters.
Definition: blackoilpolymermodules.hh:96
static TabulatedTwoDFunction & getSkprwatTable(const int tableNumber)
get the SKPRWAT table
Definition: blackoilpolymermodules.hh:118
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:353
static Scalar plymixparToddLongstaff(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:441
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:278
static Scalar plyrockMaxAdsorbtion(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:401
static Scalar plyrockResidualResistanceFactor(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilpolymermodules.hh:377
static Evaluation computeShearFactor(const Evaluation &polymerConcentration, unsigned pvtnumRegionIdx, const Evaluation &v0)
Computes the shear factor.
Definition: blackoilpolymermodules.hh:476
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: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
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