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