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
35#include <opm/material/common/MathToolbox.hpp>
36
39
41
43
44#include <algorithm>
45#include <cassert>
46#include <cmath>
47#include <cstddef>
48#include <istream>
49#include <memory>
50#include <ostream>
51#include <stdexcept>
52#include <string>
53#include <vector>
54
55namespace Opm {
56
62template <class TypeTag, bool enablePolymerV = getPropValue<TypeTag, Properties::EnablePolymer>()>
64{
76
77 using Toolbox = MathToolbox<Evaluation>;
78
79 using TabulatedFunction = typename BlackOilPolymerParams<Scalar>::TabulatedFunction;
80 using TabulatedTwoDFunction = typename BlackOilPolymerParams<Scalar>::TabulatedTwoDFunction;
81
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;
87
88 static constexpr unsigned enablePolymer = enablePolymerV;
89 static constexpr bool enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>();
90
91 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
92
93public:
96 {
97 params_ = params;
98 }
99
103 static TabulatedTwoDFunction& getPlymwinjTable(const int tableNumber)
104 {
105 const auto iterTable = params_.plymwinjTables_.find(tableNumber);
106 if (iterTable != params_.plymwinjTables_.end()) {
107 return iterTable->second;
108 }
109 else {
110 throw std::runtime_error(" the PLYMWINJ table " + std::to_string(tableNumber) + " does not exist\n");
111 }
112 }
113
117 static TabulatedTwoDFunction& getSkprwatTable(const int tableNumber)
118 {
119 const auto iterTable = params_.skprwatTables_.find(tableNumber);
120 if (iterTable != params_.skprwatTables_.end()) {
121 return iterTable->second;
122 }
123 else {
124 throw std::runtime_error(" the SKPRWAT table " + std::to_string(tableNumber) + " does not exist\n");
125 }
126 }
127
132 getSkprpolyTable(const int tableNumber)
133 {
134 const auto iterTable = params_.skprpolyTables_.find(tableNumber);
135 if (iterTable != params_.skprpolyTables_.end()) {
136 return iterTable->second;
137 }
138 else {
139 throw std::runtime_error(" the SKPRPOLY table " + std::to_string(tableNumber) + " does not exist\n");
140 }
141 }
142
146 static void registerParameters()
147 {
148 if constexpr (enablePolymer) {
150 }
151 }
152
156 static void registerOutputModules(Model& model,
157 Simulator& simulator)
158 {
159 if constexpr (enablePolymer) {
160 model.addOutputModule(std::make_unique<VtkBlackOilPolymerModule<TypeTag>>(simulator));
161 }
162 }
163
164 static bool primaryVarApplies(unsigned pvIdx)
165 {
166 if constexpr (enablePolymer) {
167 if constexpr (enablePolymerMolarWeight) {
168 return pvIdx == polymerConcentrationIdx || pvIdx == polymerMoleWeightIdx;
169 }
170 else {
171 return pvIdx == polymerConcentrationIdx;
172 }
173 }
174 else {
175 return false;
176 }
177 }
178
179 static std::string primaryVarName(unsigned pvIdx)
180 {
181 assert(primaryVarApplies(pvIdx));
182
183 if (pvIdx == polymerConcentrationIdx) {
184 return "polymer_waterconcentration";
185 }
186 else {
187 return "polymer_molecularweight";
188 }
189 }
190
191 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
192 {
193 assert(primaryVarApplies(pvIdx));
194
195 // TODO: it may be beneficial to chose this differently.
196 return static_cast<Scalar>(1.0);
197 }
198
199 static bool eqApplies(unsigned eqIdx)
200 {
201 if constexpr (enablePolymer) {
202 if constexpr (enablePolymerMolarWeight) {
203 return eqIdx == contiPolymerEqIdx || eqIdx == contiPolymerMolarWeightEqIdx;
204 }
205 else {
206 return eqIdx == contiPolymerEqIdx;
207 }
208 }
209 else {
210 return false;
211 }
212 }
213
214 static std::string eqName(unsigned eqIdx)
215 {
216 assert(eqApplies(eqIdx));
217
218 if (eqIdx == contiPolymerEqIdx) {
219 return "conti^polymer";
220 }
221 else {
222 return "conti^polymer_molecularweight";
223 }
224 }
225
226 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
227 {
228 assert(eqApplies(eqIdx));
229
230 // TODO: it may be beneficial to chose this differently.
231 return static_cast<Scalar>(1.0);
232 }
233
234 // must be called after water storage is computed
235 template <class LhsEval>
236 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
237 const IntensiveQuantities& intQuants)
238 {
239 if constexpr (enablePolymer) {
240 const auto& fs = intQuants.fluidState();
241
242 // avoid singular matrix if no water is present.
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()),
247 1e-10);
248
249 // polymer in water phase
250 const LhsEval massPolymer =
251 surfaceVolumeWater *
252 Toolbox::template decay<LhsEval>(intQuants.polymerConcentration()) *
253 (1.0 - Toolbox::template decay<LhsEval>(intQuants.polymerDeadPoreVolume()));
254
255 // polymer in solid phase
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());
260
261 LhsEval accumulationPolymer = massPolymer + adsorptionPolymer;
262
263 storage[contiPolymerEqIdx] += accumulationPolymer;
264
265 // tracking the polymer molecular weight
266 if constexpr (enablePolymerMolarWeight) {
267 accumulationPolymer = max(accumulationPolymer, 1e-10);
268
269 storage[contiPolymerMolarWeightEqIdx] +=
270 accumulationPolymer * Toolbox::template decay<LhsEval>(intQuants.polymerMoleWeight());
271 }
272 }
273 }
274
275 static void computeFlux([[maybe_unused]] RateVector& flux,
276 [[maybe_unused]] const ElementContext& elemCtx,
277 [[maybe_unused]] unsigned scvfIdx,
278 [[maybe_unused]] unsigned timeIdx)
279 {
280 if constexpr (enablePolymer) {
281 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
282
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);
288
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();
296
297 // modify water
298 flux[contiWaterEqIdx] /= extQuants.waterShearFactor();
299 }
300 else {
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());
307
308 // modify water
309 flux[contiWaterEqIdx] /= decay<Scalar>(extQuants.waterShearFactor());
310 }
311
312 // flux related to transport of polymer molecular weight
313 if constexpr (enablePolymerMolarWeight) {
314 if (upIdx == inIdx) {
315 flux[contiPolymerMolarWeightEqIdx] =
316 flux[contiPolymerEqIdx] * up.polymerMoleWeight();
317 }
318 else {
319 flux[contiPolymerMolarWeightEqIdx] =
320 flux[contiPolymerEqIdx] * decay<Scalar>(up.polymerMoleWeight());
321 }
322 }
323 }
324 }
325
329 static Scalar computeUpdateError(const PrimaryVariables&,
330 const EqVector&)
331 {
332 // do not consider consider the cange of polymer primary variables for
333 // convergence
334 // TODO: maybe this should be changed
335 return static_cast<Scalar>(0.0);
336 }
337
338 template <class DofEntity>
339 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
340 {
341 if constexpr (enablePolymer) {
342 const unsigned dofIdx = model.dofMapper().index(dof);
343 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
344 outstream << priVars[polymerConcentrationIdx];
345 outstream << priVars[polymerMoleWeightIdx];
346 }
347 }
348
349 template <class DofEntity>
350 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
351 {
352 if constexpr (enablePolymer) {
353 const unsigned dofIdx = model.dofMapper().index(dof);
354 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
355 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
356
357 instream >> priVars0[polymerConcentrationIdx];
358 instream >> priVars0[polymerMoleWeightIdx];
359
360 // set the primary variables for the beginning of the current time step.
361 priVars1[polymerConcentrationIdx] = priVars0[polymerConcentrationIdx];
362 priVars1[polymerMoleWeightIdx] = priVars0[polymerMoleWeightIdx];
363 }
364 }
365
366 static Scalar plyrockDeadPoreVolume(const ElementContext& elemCtx,
367 unsigned scvIdx,
368 unsigned timeIdx)
369 {
370 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
371 return params_.plyrockDeadPoreVolume_[satnumRegionIdx];
372 }
373
374 static Scalar plyrockResidualResistanceFactor(const ElementContext& elemCtx,
375 unsigned scvIdx,
376 unsigned timeIdx)
377 {
378 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
379 return params_.plyrockResidualResistanceFactor_[satnumRegionIdx];
380 }
381
382 static Scalar plyrockRockDensityFactor(const ElementContext& elemCtx,
383 unsigned scvIdx,
384 unsigned timeIdx)
385 {
386 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
387 return params_.plyrockRockDensityFactor_[satnumRegionIdx];
388 }
389
390 static Scalar plyrockAdsorbtionIndex(const ElementContext& elemCtx,
391 unsigned scvIdx,
392 unsigned timeIdx)
393 {
394 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
395 return params_.plyrockAdsorbtionIndex_[satnumRegionIdx];
396 }
397
398 static Scalar plyrockMaxAdsorbtion(const ElementContext& elemCtx,
399 unsigned scvIdx,
400 unsigned timeIdx)
401 {
402 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
403 return params_.plyrockMaxAdsorbtion_[satnumRegionIdx];
404 }
405
406 static const TabulatedFunction&
407 plyadsAdsorbedPolymer(const ElementContext& elemCtx,
408 unsigned scvIdx,
409 unsigned timeIdx)
410 {
411 const unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
412 return params_.plyadsAdsorbedPolymer_[satnumRegionIdx];
413 }
414
415 static const TabulatedFunction&
416 plyviscViscosityMultiplierTable(const ElementContext& elemCtx,
417 unsigned scvIdx,
418 unsigned timeIdx)
419 {
420 unsigned pvtnumRegionIdx =
421 elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
422 return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
423 }
424
425 static const TabulatedFunction&
426 plyviscViscosityMultiplierTable(unsigned pvtnumRegionIdx)
427 { return params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx]; }
428
429 static Scalar plymaxMaxConcentration(const ElementContext& elemCtx,
430 unsigned scvIdx,
431 unsigned timeIdx)
432 {
433 const unsigned polymerMixRegionIdx =
434 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
435 return params_.plymaxMaxConcentration_[polymerMixRegionIdx];
436 }
437
438 static Scalar plymixparToddLongstaff(const ElementContext& elemCtx,
439 unsigned scvIdx,
440 unsigned timeIdx)
441 {
442 const unsigned polymerMixRegionIdx =
443 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
444 return params_.plymixparToddLongstaff_[polymerMixRegionIdx];
445 }
446
448 plyvmhCoefficients(const ElementContext& elemCtx,
449 const unsigned scvIdx,
450 const unsigned timeIdx)
451 {
452 const unsigned polymerMixRegionIdx =
453 elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
454 return params_.plyvmhCoefficients_[polymerMixRegionIdx];
455 }
456
457 static bool hasPlyshlog()
458 { return params_.hasPlyshlog_; }
459
460 static bool hasShrate()
461 { return params_.hasShrate_; }
462
463 static Scalar shrate(unsigned pvtnumRegionIdx)
464 { return params_.shrate_[pvtnumRegionIdx]; }
465
472 template <class Evaluation>
473 static Evaluation computeShearFactor(const Evaluation& polymerConcentration,
474 unsigned pvtnumRegionIdx,
475 const Evaluation& v0)
476 {
477 using ToolboxLocal = MathToolbox<Evaluation>;
478
479 const auto& viscosityMultiplierTable = params_.plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
480 const Scalar viscosityMultiplier =
481 viscosityMultiplierTable.eval(scalarValue(polymerConcentration), /*extrapolate=*/true);
482
483 const Scalar eps = 1e-14;
484 // return 1.0 if the polymer has no effect on the water.
485 if (std::abs((viscosityMultiplier - 1.0)) < eps) {
486 return ToolboxLocal::createConstant(v0, 1.0);
487 }
488
489 const std::vector<Scalar>& shearEffectRefLogVelocity =
490 params_.plyshlogShearEffectRefLogVelocity_[pvtnumRegionIdx];
491 const auto v0AbsLog = log(abs(v0));
492 // return 1.0 if the velocity /sharte is smaller than the first velocity entry.
493 if (v0AbsLog < shearEffectRefLogVelocity[0]) {
494 return ToolboxLocal::createConstant(v0, 1.0);
495 }
496
497 // compute shear factor from input
498 // Z = (1 + (P - 1) * M(v)) / P
499 // where M(v) is computed from user input
500 // and P = viscosityMultiplier
501 const std::vector<Scalar>& shearEffectRefMultiplier =
502 params_.plyshlogShearEffectRefMultiplier_[pvtnumRegionIdx];
503 const std::size_t numTableEntries = shearEffectRefLogVelocity.size();
504 assert(shearEffectRefMultiplier.size() == numTableEntries);
505
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]);
511 }
512 // store the logarithmic velocity and logarithmic multipliers in a table for easy look up and
513 // linear interpolation in the logarithmic space.
514 const TabulatedFunction logShearEffectMultiplier =
515 TabulatedFunction(numTableEntries, shearEffectRefLogVelocity,
516 shearEffectMultiplier, /*bool sortInputs =*/ false);
517
518 // Find sheared velocity (v) that satisfies
519 // F = log(v) + log (Z) - log(v0) = 0;
520
521 // Set up the function
522 // u = log(v)
523 auto F = [&logShearEffectMultiplier, &v0AbsLog](const Evaluation& u) {
524 return u + logShearEffectMultiplier.eval(u, true) - v0AbsLog;
525 };
526 // and its derivative
527 auto dF = [&logShearEffectMultiplier](const Evaluation& u) {
528 return 1 + logShearEffectMultiplier.evalDerivative(u, true);
529 };
530
531 // Solve F = 0 using Newton
532 // Use log(v0) as initial value for u
533 auto u = v0AbsLog;
534 bool converged = false;
535 // TODO make this into parameters
536 for (int i = 0; i < 20; ++i) {
537 const auto f = F(u);
538 const auto df = dF(u);
539 u -= f / df;
540 if (std::abs(scalarValue(f)) < 1e-12) {
541 converged = true;
542 break;
543 }
544 }
545 if (!converged) {
546 throw std::runtime_error("Not able to compute shear velocity. \n");
547 }
548
549 // return the shear factor
550 return exp(logShearEffectMultiplier.eval(u, /*extrapolate=*/true));
551 }
552
553 Scalar molarMass() const
554 { return 0.25; } // kg/mol
555
556private:
557 static BlackOilPolymerParams<Scalar> params_;
558};
559
560template <class TypeTag, bool enablePolymerV>
561BlackOilPolymerParams<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
562BlackOilPolymerModule<TypeTag, enablePolymerV>::params_;
563
564template <class TypeTag, bool enablePolymerV>
566
574template <class TypeTag>
575class BlackOilPolymerIntensiveQuantities<TypeTag, /*enablePolymerV=*/true>
576{
578
586
588
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;
593
594public:
600 void polymerPropertiesUpdate_(const ElementContext& elemCtx,
601 unsigned dofIdx,
602 unsigned timeIdx)
603 {
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);
609 }
610
611 // permeability reduction due to polymer
612 const Scalar& maxAdsorbtion = PolymerModule::plyrockMaxAdsorbtion(elemCtx, dofIdx, timeIdx);
613 const auto& plyadsAdsorbedPolymer = PolymerModule::plyadsAdsorbedPolymer(elemCtx, dofIdx, timeIdx);
614 polymerAdsorption_ = plyadsAdsorbedPolymer.eval(polymerConcentration_, /*extrapolate=*/true);
615 if (static_cast<int>(PolymerModule::plyrockAdsorbtionIndex(elemCtx, dofIdx, timeIdx)) ==
617 {
618 const auto maxPolymerAdsorption =
619 elemCtx.problem().maxPolymerAdsorption(elemCtx, dofIdx, timeIdx);
620 polymerAdsorption_ = std::max(Evaluation(maxPolymerAdsorption), polymerAdsorption_);
621 }
622
623 // compute resitanceFactor
624 const Scalar& residualResistanceFactor =
625 PolymerModule::plyrockResidualResistanceFactor(elemCtx, dofIdx, timeIdx);
626 const Evaluation resistanceFactor = 1.0 + (residualResistanceFactor - 1.0) *
627 polymerAdsorption_ / maxAdsorbtion;
628
629 // compute effective viscosities
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_, /*extrapolate=*/true) * muWater;
638
639 // Do the Todd-Longstaff mixing
640 const Scalar plymixparToddLongstaff = PolymerModule::plymixparToddLongstaff(elemCtx, dofIdx, timeIdx);
641 const Evaluation viscosityPolymer = viscosityMultiplier.eval(cmax, /*extrapolate=*/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);
646
647 const Evaluation cbar = polymerConcentration_ / cmax;
648 // waterViscosity / effectiveWaterViscosity
649 waterViscosityCorrection_ = muWater * ((1.0 - cbar) / viscosityWaterEffective + cbar / viscosityPolymerEffective);
650 // effectiveWaterViscosity / effectivePolymerViscosity
651 polymerViscosityCorrection_ = (muWater / waterViscosityCorrection_) / viscosityPolymerEffective;
652 }
653 else { // based on PLYVMH
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;
659
660 // viscosity model based on Mark-Houwink equation and Huggins equation
661 // 1000 is a emperical constant, most likely related to unit conversion
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;
666 }
667
668 // adjust water mobility
669 asImp_().mobility_[waterPhaseIdx] *= waterViscosityCorrection_ / resistanceFactor;
670
671 // update rock properties
672 polymerDeadPoreVolume_ = PolymerModule::plyrockDeadPoreVolume(elemCtx, dofIdx, timeIdx);
673 polymerRockDensity_ = PolymerModule::plyrockRockDensityFactor(elemCtx, dofIdx, timeIdx);
674 }
675
676 const Evaluation& polymerConcentration() const
677 { return polymerConcentration_; }
678
679 const Evaluation& polymerMoleWeight() const
680 {
681 if constexpr (enablePolymerMolarWeight) {
682 return polymerMoleWeight_;
683 }
684 else {
685 throw std::logic_error("polymerMoleWeight() is called but polymer milecular weight is disabled");
686 }
687 }
688
690 { return polymerDeadPoreVolume_; }
691
692 const Evaluation& polymerAdsorption() const
693 { return polymerAdsorption_; }
694
695 Scalar polymerRockDensity() const
696 { return polymerRockDensity_; }
697
698 // effectiveWaterViscosity / effectivePolymerViscosity
699 const Evaluation& polymerViscosityCorrection() const
700 { return polymerViscosityCorrection_; }
701
702 // waterViscosity / effectiveWaterViscosity
703 const Evaluation& waterViscosityCorrection() const
704 { return waterViscosityCorrection_; }
705
706protected:
707 Implementation& asImp_()
708 { return *static_cast<Implementation*>(this); }
709
711 // polymer molecular weight
718};
719
720template <class TypeTag>
722{
725
726public:
727 void polymerPropertiesUpdate_(const ElementContext&,
728 unsigned,
729 unsigned)
730 {}
731
732 const Evaluation& polymerMoleWeight() const
733 { throw std::logic_error("polymerMoleWeight() called but polymer molecular weight is disabled"); }
734
735 const Evaluation& polymerConcentration() const
736 { throw std::runtime_error("polymerConcentration() called but polymers are disabled"); }
737
738 const Evaluation& polymerDeadPoreVolume() const
739 { throw std::runtime_error("polymerDeadPoreVolume() called but polymers are disabled"); }
740
741 const Evaluation& polymerAdsorption() const
742 { throw std::runtime_error("polymerAdsorption() called but polymers are disabled"); }
743
744 const Evaluation& polymerRockDensity() const
745 { throw std::runtime_error("polymerRockDensity() called but polymers are disabled"); }
746
747 const Evaluation& polymerViscosityCorrection() const
748 { throw std::runtime_error("polymerViscosityCorrection() called but polymers are disabled"); }
749
750 const Evaluation& waterViscosityCorrection() const
751 { throw std::runtime_error("waterViscosityCorrection() called but polymers are disabled"); }
752};
753
754template <class TypeTag, bool enablePolymerV>
756
764template <class TypeTag>
765class BlackOilPolymerExtensiveQuantities<TypeTag, /*enablePolymerV=*/true>
766{
768
774
775 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
776
778
779public:
786 template <class Dummy = bool> // we need to make this method a template to avoid
787 // compiler errors if it is not instantiated!
788 void updateShearMultipliersPerm(const ElementContext&,
789 unsigned,
790 unsigned)
791 {
792 throw std::runtime_error("The extension of the blackoil model for polymers is not yet "
793 "implemented for problems specified using permeabilities.");
794 }
795
802 template <class Dummy = bool> // we need to make this method a template to avoid
803 // compiler errors if it is not instantiated!
804 void updateShearMultipliers(const ElementContext& elemCtx,
805 unsigned scvfIdx,
806 unsigned timeIdx)
807 {
808 waterShearFactor_ = 1.0;
809 polymerShearFactor_ = 1.0;
810
811 if (!PolymerModule::hasPlyshlog()) {
812 return;
813 }
814
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);
822
823 // compute water velocity from flux
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;
832
833 // guard against zero porosity and no mobile water
834 const Evaluation denom = max(poroAvg * (Sw - Swcr), 1e-12);
835 Evaluation waterVolumeVelocity = extQuants.volumeFlux(waterPhaseIdx) / denom;
836
837 // if shrate is specified. Compute shrate based on the water velocity
838 if (PolymerModule::hasShrate()) {
839 const Evaluation& relWater = up.relativePermeability(waterPhaseIdx);
840 const Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
841 if (trans > 0.0) {
842 const Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
843 const auto dist = elemCtx.pos(interiorDofIdx, timeIdx) - elemCtx.pos(exteriorDofIdx, timeIdx);
844 // compute permeability from transmissibility.
845 const Scalar absPerm = trans / faceArea * dist.two_norm();
846 waterVolumeVelocity *=
847 PolymerModule::shrate(pvtnumRegionIdx) * sqrt(poroAvg * Sw / (relWater * absPerm));
848 assert(isfinite(waterVolumeVelocity));
849 }
850 }
851
852 // compute share factors for water and polymer
853 waterShearFactor_ =
854 PolymerModule::computeShearFactor(up.polymerConcentration(),
855 pvtnumRegionIdx,
856 waterVolumeVelocity);
857 polymerShearFactor_ =
858 PolymerModule::computeShearFactor(up.polymerConcentration(),
859 pvtnumRegionIdx,
860 waterVolumeVelocity * up.polymerViscosityCorrection());
861 }
862
863 const Evaluation& polymerShearFactor() const
864 { return polymerShearFactor_; }
865
866 const Evaluation& waterShearFactor() const
867 { return waterShearFactor_; }
868
869private:
870 Implementation& asImp_()
871 { return *static_cast<Implementation*>(this); }
872
873 Evaluation polymerShearFactor_;
874 Evaluation waterShearFactor_;
875};
876
877template <class TypeTag>
879{
882
883public:
884 void updateShearMultipliers(const ElementContext&,
885 unsigned,
886 unsigned)
887 {}
888
889 void updateShearMultipliersPerm(const ElementContext&,
890 unsigned,
891 unsigned)
892 {}
893
894 const Evaluation& polymerShearFactor() const
895 { throw std::runtime_error("polymerShearFactor() called but polymers are disabled"); }
896
897 const Evaluation& waterShearFactor() const
898 { throw std::runtime_error("waterShearFactor() called but polymers are disabled"); }
899};
900
901} // namespace Opm
902
903#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: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 > &&params)
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