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