-
Notifications
You must be signed in to change notification settings - Fork 67
/
Copy pathblackoilpolymermodules.hh
1349 lines (1129 loc) · 57.7 KB
/
blackoilpolymermodules.hh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 2 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
/*!
* \file
*
* \brief Contains the classes required to extend the black-oil model by polymer.
*/
#ifndef EWOMS_BLACK_OIL_POLYMER_MODULE_HH
#define EWOMS_BLACK_OIL_POLYMER_MODULE_HH
#include "blackoilproperties.hh"
#include <opm/models/io/vtkblackoilpolymermodule.hh>
#include <opm/models/common/quantitycallbacks.hh>
#include <opm/material/common/Tabulated1DFunction.hpp>
#include <opm/material/common/IntervalTabulated2DFunction.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PlyadsTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PlymaxTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PlyrockTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PlyshlogTable.hpp>
#include <opm/parser/eclipse/EclipseState/Tables/PlyviscTable.hpp>
#include <opm/material/common/Valgrind.hpp>
#include <opm/material/common/Unused.hpp>
#include <opm/material/common/Exceptions.hpp>
#include <dune/common/fvector.hh>
#include <string>
namespace Opm {
/*!
* \ingroup BlackOil
* \brief Contains the high level supplements required to extend the black oil
* model by polymer.
*/
template <class TypeTag, bool enablePolymerV = GET_PROP_VALUE(TypeTag, EnablePolymer)>
class BlackOilPolymerModule
{
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
typedef typename GET_PROP_TYPE(TypeTag, Model) Model;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
typedef Opm::MathToolbox<Evaluation> Toolbox;
typedef typename Opm::Tabulated1DFunction<Scalar> TabulatedFunction;
typedef typename Opm::IntervalTabulated2DFunction<Scalar> TabulatedTwoDFunction;
static constexpr unsigned polymerConcentrationIdx = Indices::polymerConcentrationIdx;
static constexpr unsigned polymerMoleWeightIdx = Indices::polymerMoleWeightIdx;
static constexpr unsigned contiPolymerEqIdx = Indices::contiPolymerEqIdx;
static constexpr unsigned contiPolymerMolarWeightEqIdx = Indices::contiPolymerMWEqIdx;
static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
static constexpr unsigned enablePolymer = enablePolymerV;
static constexpr bool enablePolymerMolarWeight = GET_PROP_VALUE(TypeTag, EnablePolymerMW);
static constexpr unsigned numEq = GET_PROP_VALUE(TypeTag, NumEq);
static constexpr unsigned numPhases = FluidSystem::numPhases;
struct SkprpolyTable {
double refConcentration;
TabulatedTwoDFunction table_func;
};
public:
enum AdsorptionBehaviour { Desorption = 1, NoDesorption = 2 };
// a struct containing the constants to calculate polymer viscosity
// based on Mark-Houwink equation and Huggins equation, the constants are provided
// by the keyword PLYVMH
struct PlyvmhCoefficients {
Scalar k_mh;
Scalar a_mh;
Scalar gamma;
Scalar kappa;
};
/*!
* \brief Initialize all internal data structures needed by the polymer module
*/
static void initFromDeck(const Opm::Deck& deck, const Opm::EclipseState& eclState)
{
// some sanity checks: if polymers are enabled, the POLYMER keyword must be
// present, if polymers are disabled the keyword must not be present.
if (enablePolymer && !deck.hasKeyword("POLYMER")) {
throw std::runtime_error("Non-trivial polymer treatment requested at compile time, but "
"the deck does not contain the POLYMER keyword");
}
else if (!enablePolymer && deck.hasKeyword("POLYMER")) {
throw std::runtime_error("Polymer treatment disabled at compile time, but the deck "
"contains the POLYMER keyword");
}
if (enablePolymerMolarWeight && !deck.hasKeyword("POLYMW")) {
throw std::runtime_error("Polymer molecular weight tracking is enabled at compile time, but "
"the deck does not contain the POLYMW keyword");
}
else if (!enablePolymerMolarWeight && deck.hasKeyword("POLYMW")) {
throw std::runtime_error("Polymer molecular weight tracking is disabled at compile time, but the deck "
"contains the POLYMW keyword");
}
if (enablePolymerMolarWeight && !enablePolymer) {
throw std::runtime_error("Polymer molecular weight tracking is enabled while polymer treatment "
"is disabled at compile time");
}
if (!deck.hasKeyword("POLYMER"))
return; // polymer treatment is supposed to be disabled
const auto& tableManager = eclState.getTableManager();
unsigned numSatRegions = tableManager.getTabdims().getNumSatTables();
setNumSatRegions(numSatRegions);
// initialize the objects which deal with the PLYROCK keyword
const auto& plyrockTables = tableManager.getPlyrockTables();
if (!plyrockTables.empty()) {
assert(numSatRegions == plyrockTables.size());
for (unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
const auto& plyrockTable = plyrockTables.template getTable<Opm::PlyrockTable>(satRegionIdx);
setPlyrock(satRegionIdx,
plyrockTable.getDeadPoreVolumeColumn()[0],
plyrockTable.getResidualResistanceFactorColumn()[0],
plyrockTable.getRockDensityFactorColumn()[0],
static_cast<AdsorptionBehaviour>(plyrockTable.getAdsorbtionIndexColumn()[0]),
plyrockTable.getMaxAdsorbtionColumn()[0]);
}
}
else {
throw std::runtime_error("PLYROCK must be specified in POLYMER runs\n");
}
// initialize the objects which deal with the PLYADS keyword
const auto& plyadsTables = tableManager.getPlyadsTables();
if (!plyadsTables.empty()) {
assert(numSatRegions == plyadsTables.size());
for (unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
const auto& plyadsTable = plyadsTables.template getTable<Opm::PlyadsTable>(satRegionIdx);
// Copy data
const auto& c = plyadsTable.getPolymerConcentrationColumn();
const auto& ads = plyadsTable.getAdsorbedPolymerColumn();
plyadsAdsorbedPolymer_[satRegionIdx].setXYContainers(c, ads);
}
}
else {
throw std::runtime_error("PLYADS must be specified in POLYMER runs\n");
}
unsigned numPvtRegions = tableManager.getTabdims().getNumPVTTables();
setNumPvtRegions(numPvtRegions);
// initialize the objects which deal with the PLYVISC keyword
const auto& plyviscTables = tableManager.getPlyviscTables();
if (!plyviscTables.empty()) {
// different viscosity model is used for POLYMW
if (enablePolymerMolarWeight) {
Opm::OpmLog::warning("PLYVISC should not be used in POLYMW runs, "
"it will have no effect. A viscosity model based on PLYVMH is used instead.\n");
}
else {
assert(numPvtRegions == plyviscTables.size());
for (unsigned pvtRegionIdx = 0; pvtRegionIdx < numPvtRegions; ++ pvtRegionIdx) {
const auto& plyadsTable = plyviscTables.template getTable<Opm::PlyviscTable>(pvtRegionIdx);
// Copy data
const auto& c = plyadsTable.getPolymerConcentrationColumn();
const auto& visc = plyadsTable.getViscosityMultiplierColumn();
plyviscViscosityMultiplierTable_[pvtRegionIdx].setXYContainers(c, visc);
}
}
}
else if (!enablePolymerMolarWeight) {
throw std::runtime_error("PLYVISC must be specified in POLYMER runs\n");
}
// initialize the objects which deal with the PLYMAX keyword
const auto& plymaxTables = tableManager.getPlymaxTables();
const unsigned numMixRegions = plymaxTables.size();
setNumMixRegions(numMixRegions);
if (!plymaxTables.empty()) {
for (unsigned mixRegionIdx = 0; mixRegionIdx < numMixRegions; ++ mixRegionIdx) {
const auto& plymaxTable = plymaxTables.template getTable<Opm::PlymaxTable>(mixRegionIdx);
setPlymax(mixRegionIdx, plymaxTable.getPolymerConcentrationColumn()[mixRegionIdx]);
}
}
else {
throw std::runtime_error("PLYMAX must be specified in POLYMER runs\n");
}
if (deck.hasKeyword("PLMIXPAR")) {
if (enablePolymerMolarWeight) {
Opm::OpmLog::warning("PLMIXPAR should not be used in POLYMW runs, it will have no effect.\n");
}
else {
// initialize the objects which deal with the PLMIXPAR keyword
for (unsigned mixRegionIdx = 0; mixRegionIdx < numMixRegions; ++ mixRegionIdx) {
const auto& plmixparRecord = deck.getKeyword("PLMIXPAR").getRecord(mixRegionIdx);
setPlmixpar(mixRegionIdx, plmixparRecord.getItem("TODD_LONGSTAFF").getSIDouble(0));
}
}
}
else if (!enablePolymerMolarWeight) {
throw std::runtime_error("PLMIXPAR must be specified in POLYMER runs\n");
}
hasPlyshlog_ = deck.hasKeyword("PLYSHLOG");
hasShrate_ = deck.hasKeyword("SHRATE");
if ((hasPlyshlog_ || hasShrate_) && enablePolymerMolarWeight) {
Opm::OpmLog::warning("PLYSHLOG and SHRATE should not be used in POLYMW runs, they will have no effect.\n");
}
if (hasPlyshlog_ && !enablePolymerMolarWeight) {
const auto& plyshlogTables = tableManager.getPlyshlogTables();
assert(numPvtRegions == plyshlogTables.size());
plyshlogShearEffectRefMultiplier_.resize(numPvtRegions);
plyshlogShearEffectRefLogVelocity_.resize(numPvtRegions);
for (unsigned pvtRegionIdx = 0; pvtRegionIdx < numPvtRegions; ++ pvtRegionIdx) {
const auto& plyshlogTable = plyshlogTables.template getTable<Opm::PlyshlogTable>(pvtRegionIdx);
Scalar plyshlogRefPolymerConcentration = plyshlogTable.getRefPolymerConcentration();
auto waterVelocity = plyshlogTable.getWaterVelocityColumn().vectorCopy();
auto shearMultiplier = plyshlogTable.getShearMultiplierColumn().vectorCopy();
// do the unit version here for the waterVelocity
Opm::UnitSystem unitSystem = deck.getActiveUnitSystem();
double siFactor = hasShrate_? unitSystem.parse("1/Time").getSIScaling() : unitSystem.parse("Length/Time").getSIScaling();
for (size_t i = 0; i < waterVelocity.size(); ++i) {
waterVelocity[i] *= siFactor;
// for plyshlog the input must be stored as logarithms
// the interpolation is then done the log-space.
waterVelocity[i] = std::log(waterVelocity[i]);
}
Scalar refViscMult = plyviscViscosityMultiplierTable_[pvtRegionIdx].eval(plyshlogRefPolymerConcentration, /*extrapolate=*/true);
// convert the table using referece conditions
for (size_t i = 0; i < waterVelocity.size(); ++i) {
shearMultiplier[i] *= refViscMult;
shearMultiplier[i] -= 1;
shearMultiplier[i] /= (refViscMult - 1);
shearMultiplier[i] = shearMultiplier[i];
}
plyshlogShearEffectRefMultiplier_[pvtRegionIdx].resize(waterVelocity.size());
plyshlogShearEffectRefLogVelocity_[pvtRegionIdx].resize(waterVelocity.size());
for (size_t i = 0; i < waterVelocity.size(); ++i) {
plyshlogShearEffectRefMultiplier_[pvtRegionIdx][i] = shearMultiplier[i];
plyshlogShearEffectRefLogVelocity_[pvtRegionIdx][i] = waterVelocity[i];
}
}
}
if (hasShrate_ && !enablePolymerMolarWeight) {
if(!hasPlyshlog_) {
throw std::runtime_error("PLYSHLOG must be specified if SHRATE is used in POLYMER runs\n");
}
const auto& shrateKeyword = deck.getKeyword("SHRATE");
const std::vector<double>& shrateFromDeck = shrateKeyword.getSIDoubleData();
shrate_.resize(numPvtRegions);
for (unsigned pvtRegionIdx = 0; pvtRegionIdx < numPvtRegions; ++ pvtRegionIdx) {
if (shrateFromDeck.empty()) {
shrate_[pvtRegionIdx] = 4.8; //default;
}
else if (shrateFromDeck.size() == numPvtRegions) {
shrate_[pvtRegionIdx] = shrateKeyword.getSIDoubleData()[pvtRegionIdx];
}
else {
throw std::runtime_error("SHRATE must either have 0 or number of NUMPVT entries\n");
}
}
}
if (enablePolymerMolarWeight) {
const Opm::DeckKeyword& plyvmhKeyword = deck.getKeyword("PLYVMH");
assert(plyvmhKeyword.size() == numMixRegions);
if (plyvmhKeyword.size() > 0) {
for (size_t regionIdx = 0; regionIdx < plyvmhKeyword.size(); ++regionIdx) {
const Opm::DeckRecord& record = plyvmhKeyword.getRecord(regionIdx);
plyvmhCoefficients_[regionIdx].k_mh = record.getItem("K_MH").getSIDouble(0);
plyvmhCoefficients_[regionIdx].a_mh = record.getItem("A_MH").getSIDouble(0);
plyvmhCoefficients_[regionIdx].gamma = record.getItem("GAMMA").getSIDouble(0);
plyvmhCoefficients_[regionIdx].kappa = record.getItem("KAPPA").getSIDouble(0);
}
}
else {
throw std::runtime_error("PLYVMH keyword must be specified in POLYMW rus \n");
}
// handling PLYMWINJ keyword
const auto& plymwinjTables = tableManager.getPlymwinjTables();
for (const auto& table : plymwinjTables) {
const int tableNumber = table.first;
const auto& plymwinjtable = table.second;
const std::vector<double>& throughput = plymwinjtable.getThroughputs();
const std::vector<double>& watervelocity = plymwinjtable.getVelocities();
const std::vector<std::vector<double>>& molecularweight = plymwinjtable.getMoleWeights();
TabulatedTwoDFunction tablefunc(throughput, watervelocity, molecularweight, true, false);
plymwinjTables_[tableNumber] = std::move(tablefunc);
}
// handling SKPRWAT keyword
const auto& skprwatTables = tableManager.getSkprwatTables();
for (const auto& table : skprwatTables) {
const int tableNumber = table.first;
const auto& skprwattable = table.second;
const std::vector<double>& throughput = skprwattable.getThroughputs();
const std::vector<double>& watervelocity = skprwattable.getVelocities();
const std::vector<std::vector<double>>& skinpressure = skprwattable.getSkinPressures();
TabulatedTwoDFunction tablefunc(throughput, watervelocity, skinpressure, true, false);
skprwatTables_[tableNumber] = std::move(tablefunc);
}
// handling SKPRPOLY keyword
const auto& skprpolyTables = tableManager.getSkprpolyTables();
for (const auto& table : skprpolyTables) {
const int tableNumber = table.first;
const auto& skprpolytable = table.second;
const std::vector<double>& throughput = skprpolytable.getThroughputs();
const std::vector<double>& watervelocity = skprpolytable.getVelocities();
const std::vector<std::vector<double>>& skinpressure = skprpolytable.getSkinPressures();
const double refPolymerConcentration = skprpolytable.referenceConcentration();
SkprpolyTable tablefunc = {refPolymerConcentration,
TabulatedTwoDFunction(throughput, watervelocity, skinpressure, true, false)};
skprpolyTables_[tableNumber] = std::move(tablefunc);
}
}
}
/*!
* \brief Specify the number of satuation regions.
*
* This must be called before setting the PLYROCK and PLYADS of any region.
*/
static void setNumSatRegions(unsigned numRegions)
{
plyrockDeadPoreVolume_.resize(numRegions);
plyrockResidualResistanceFactor_.resize(numRegions);
plyrockRockDensityFactor_.resize(numRegions);
plyrockAdsorbtionIndex_.resize(numRegions);
plyrockMaxAdsorbtion_.resize(numRegions);
plyadsAdsorbedPolymer_.resize(numRegions);
}
/*!
* \brief Specify the polymer rock properties a single region.
*
* The index of specified here must be in range [0, numSatRegions)
*/
static void setPlyrock(unsigned satRegionIdx,
const Scalar& plyrockDeadPoreVolume,
const Scalar& plyrockResidualResistanceFactor,
const Scalar& plyrockRockDensityFactor,
const Scalar& plyrockAdsorbtionIndex,
const Scalar& plyrockMaxAdsorbtion)
{
plyrockDeadPoreVolume_[satRegionIdx] = plyrockDeadPoreVolume;
plyrockResidualResistanceFactor_[satRegionIdx] = plyrockResidualResistanceFactor;
plyrockRockDensityFactor_[satRegionIdx] = plyrockRockDensityFactor;
plyrockAdsorbtionIndex_[satRegionIdx] = plyrockAdsorbtionIndex;
plyrockMaxAdsorbtion_[satRegionIdx] = plyrockMaxAdsorbtion;
}
/*!
* \brief Specify the number of pvt regions.
*
* This must be called before setting the PLYVISC of any region.
*/
static void setNumPvtRegions(unsigned numRegions)
{
plyviscViscosityMultiplierTable_.resize(numRegions);
}
/*!
* \brief Specify the polymer viscosity a single region.
*
* The index of specified here must be in range [0, numSatRegions)
*/
static void setPlyvisc(unsigned satRegionIdx,
const TabulatedFunction& plyviscViscosityMultiplierTable)
{
plyviscViscosityMultiplierTable_[satRegionIdx] = plyviscViscosityMultiplierTable;
}
/*!
* \brief Specify the number of mix regions.
*
* This must be called before setting the PLYMAC and PLMIXPAR of any region.
*/
static void setNumMixRegions(unsigned numRegions)
{
plymaxMaxConcentration_.resize(numRegions);
plymixparToddLongstaff_.resize(numRegions);
if (enablePolymerMolarWeight) {
plyvmhCoefficients_.resize(numRegions);
}
}
/*!
* \brief Specify the maximum polymer concentration a single region.
*
* The index of specified here must be in range [0, numMixRegionIdx)
*/
static void setPlymax(unsigned mixRegionIdx,
const Scalar& plymaxMaxConcentration)
{
plymaxMaxConcentration_[mixRegionIdx] = plymaxMaxConcentration;
}
/*!
* \brief Specify the maximum polymer concentration a single region.
*
* The index of specified here must be in range [0, numMixRegionIdx)
*/
static void setPlmixpar(unsigned mixRegionIdx,
const Scalar& plymixparToddLongstaff)
{
plymixparToddLongstaff_[mixRegionIdx] = plymixparToddLongstaff;
}
/*!
* \brief get the PLYMWINJ table
*/
static TabulatedTwoDFunction& getPlymwinjTable(const int tableNumber)
{
const auto iterTable = plymwinjTables_.find(tableNumber);
if (iterTable != plymwinjTables_.end()) {
return iterTable->second;
}
else {
throw std::runtime_error(" the PLYMWINJ table " + std::to_string(tableNumber) + " does not exist\n");
}
}
/*!
* \brief get the SKPRWAT table
*/
static TabulatedTwoDFunction& getSkprwatTable(const int tableNumber)
{
const auto iterTable = skprwatTables_.find(tableNumber);
if (iterTable != skprwatTables_.end()) {
return iterTable->second;
}
else {
throw std::runtime_error(" the SKPRWAT table " + std::to_string(tableNumber) + " does not exist\n");
}
}
/*!
* \brief get the SKPRPOLY table
*/
static SkprpolyTable& getSkprpolyTable(const int tableNumber)
{
const auto iterTable = skprpolyTables_.find(tableNumber);
if (iterTable != skprpolyTables_.end()) {
return iterTable->second;
}
else {
throw std::runtime_error(" the SKPRPOLY table " + std::to_string(tableNumber) + " does not exist\n");
}
}
/*!
* \brief Register all run-time parameters for the black-oil polymer module.
*/
static void registerParameters()
{
if (!enablePolymer)
// polymers have been disabled at compile time
return;
Opm::VtkBlackOilPolymerModule<TypeTag>::registerParameters();
}
/*!
* \brief Register all polymer specific VTK and ECL output modules.
*/
static void registerOutputModules(Model& model,
Simulator& simulator)
{
if (!enablePolymer)
// polymers have been disabled at compile time
return;
model.addOutputModule(new Opm::VtkBlackOilPolymerModule<TypeTag>(simulator));
}
static bool primaryVarApplies(unsigned pvIdx)
{
if (!enablePolymer)
// polymers have been disabled at compile time
return false;
if (!enablePolymerMolarWeight)
return pvIdx == polymerConcentrationIdx;
// both enablePolymer and enablePolymerMolarWeight are true here
return pvIdx == polymerConcentrationIdx || pvIdx == polymerMoleWeightIdx;
}
static std::string primaryVarName(unsigned pvIdx)
{
assert(primaryVarApplies(pvIdx));
if (pvIdx == polymerConcentrationIdx) {
return "polymer_waterconcentration";
}
else {
return "polymer_molecularweight";
}
}
static Scalar primaryVarWeight(unsigned pvIdx OPM_OPTIM_UNUSED)
{
assert(primaryVarApplies(pvIdx));
// TODO: it may be beneficial to chose this differently.
return static_cast<Scalar>(1.0);
}
static bool eqApplies(unsigned eqIdx)
{
if (!enablePolymer)
return false;
if (!enablePolymerMolarWeight)
return eqIdx == contiPolymerEqIdx;
// both enablePolymer and enablePolymerMolarWeight are true here
return eqIdx == contiPolymerEqIdx || eqIdx == contiPolymerMolarWeightEqIdx;
}
static std::string eqName(unsigned eqIdx)
{
assert(eqApplies(eqIdx));
if (eqIdx == contiPolymerEqIdx)
return "conti^polymer";
else
return "conti^polymer_molecularweight";
}
static Scalar eqWeight(unsigned eqIdx OPM_OPTIM_UNUSED)
{
assert(eqApplies(eqIdx));
// TODO: it may be beneficial to chose this differently.
return static_cast<Scalar>(1.0);
}
// must be called after water storage is computed
template <class LhsEval>
static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
const IntensiveQuantities& intQuants)
{
if (!enablePolymer)
return;
const auto& fs = intQuants.fluidState();
LhsEval surfaceVolumeWater =
Toolbox::template decay<LhsEval>(fs.saturation(waterPhaseIdx))
* Toolbox::template decay<LhsEval>(fs.invB(waterPhaseIdx))
* Toolbox::template decay<LhsEval>(intQuants.porosity());
// avoid singular matrix if no water is present.
surfaceVolumeWater = Opm::max(surfaceVolumeWater, 1e-10);
// polymer in water phase
const LhsEval massPolymer = surfaceVolumeWater
* Toolbox::template decay<LhsEval>(intQuants.polymerConcentration())
* (1.0 - Toolbox::template decay<LhsEval>(intQuants.polymerDeadPoreVolume()));
// polymer in solid phase
const LhsEval adsorptionPolymer =
Toolbox::template decay<LhsEval>(1.0 - intQuants.porosity())
* Toolbox::template decay<LhsEval>(intQuants.polymerRockDensity())
* Toolbox::template decay<LhsEval>(intQuants.polymerAdsorption());
LhsEval accumulationPolymer = massPolymer + adsorptionPolymer;
storage[contiPolymerEqIdx] += accumulationPolymer;
// tracking the polymer molecular weight
if (enablePolymerMolarWeight) {
accumulationPolymer = Opm::max(accumulationPolymer, 1e-10);
storage[contiPolymerMolarWeightEqIdx] += accumulationPolymer
* Toolbox::template decay<LhsEval> (intQuants.polymerMoleWeight());
}
}
static void computeFlux(RateVector& flux,
const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx)
{
if (!enablePolymer)
return;
const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
const unsigned upIdx = extQuants.upstreamIndex(FluidSystem::waterPhaseIdx);
const unsigned inIdx = extQuants.interiorIndex();
const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
const unsigned contiWaterEqIdx = Indices::conti0EqIdx + Indices::canonicalToActiveComponentIndex(FluidSystem::waterCompIdx);
if (upIdx == inIdx) {
flux[contiPolymerEqIdx] =
extQuants.volumeFlux(waterPhaseIdx)
*up.fluidState().invB(waterPhaseIdx)
*up.polymerViscosityCorrection()
/extQuants.polymerShearFactor()
*up.polymerConcentration();
// modify water
flux[contiWaterEqIdx] /=
extQuants.waterShearFactor();
}
else {
flux[contiPolymerEqIdx] =
extQuants.volumeFlux(waterPhaseIdx)
*Opm::decay<Scalar>(up.fluidState().invB(waterPhaseIdx))
*Opm::decay<Scalar>(up.polymerViscosityCorrection())
/Opm::decay<Scalar>(extQuants.polymerShearFactor())
*Opm::decay<Scalar>(up.polymerConcentration());
// modify water
flux[contiWaterEqIdx] /=
Opm::decay<Scalar>(extQuants.waterShearFactor());
}
// flux related to transport of polymer molecular weight
if (enablePolymerMolarWeight) {
if (upIdx == inIdx)
flux[contiPolymerMolarWeightEqIdx] =
flux[contiPolymerEqIdx]*up.polymerMoleWeight();
else
flux[contiPolymerMolarWeightEqIdx] =
flux[contiPolymerEqIdx]*Opm::decay<Scalar>(up.polymerMoleWeight());
}
}
/*!
* \brief Return how much a Newton-Raphson update is considered an error
*/
static Scalar computeUpdateError(const PrimaryVariables& oldPv OPM_UNUSED,
const EqVector& delta OPM_UNUSED)
{
// do not consider consider the cange of polymer primary variables for
// convergence
// TODO: maybe this should be changed
return static_cast<Scalar>(0.0);
}
template <class DofEntity>
static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
{
if (!enablePolymer)
return;
unsigned dofIdx = model.dofMapper().index(dof);
const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
outstream << priVars[polymerConcentrationIdx];
outstream << priVars[polymerMoleWeightIdx];
}
template <class DofEntity>
static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
{
if (!enablePolymer)
return;
unsigned dofIdx = model.dofMapper().index(dof);
PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
instream >> priVars0[polymerConcentrationIdx];
instream >> priVars0[polymerMoleWeightIdx];
// set the primary variables for the beginning of the current time step.
priVars1[polymerConcentrationIdx] = priVars0[polymerConcentrationIdx];
priVars1[polymerMoleWeightIdx] = priVars0[polymerMoleWeightIdx];
}
static const Scalar plyrockDeadPoreVolume(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyrockDeadPoreVolume_[satnumRegionIdx];
}
static const Scalar plyrockResidualResistanceFactor(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyrockResidualResistanceFactor_[satnumRegionIdx];
}
static const Scalar plyrockRockDensityFactor(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyrockRockDensityFactor_[satnumRegionIdx];
}
static const Scalar plyrockAdsorbtionIndex(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyrockAdsorbtionIndex_[satnumRegionIdx];
}
static const Scalar plyrockMaxAdsorbtion(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyrockMaxAdsorbtion_[satnumRegionIdx];
}
static const TabulatedFunction& plyadsAdsorbedPolymer(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyadsAdsorbedPolymer_[satnumRegionIdx];
}
static const TabulatedFunction& plyviscViscosityMultiplierTable(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned pvtnumRegionIdx = elemCtx.problem().pvtRegionIndex(elemCtx, scvIdx, timeIdx);
return plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
}
static const TabulatedFunction& plyviscViscosityMultiplierTable(unsigned pvtnumRegionIdx)
{
return plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
}
static const Scalar plymaxMaxConcentration(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned polymerMixRegionIdx = elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plymaxMaxConcentration_[polymerMixRegionIdx];
}
static const Scalar plymixparToddLongstaff(const ElementContext& elemCtx,
unsigned scvIdx,
unsigned timeIdx)
{
unsigned polymerMixRegionIdx = elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plymixparToddLongstaff_[polymerMixRegionIdx];
}
static const PlyvmhCoefficients& plyvmhCoefficients(const ElementContext& elemCtx,
const unsigned scvIdx,
const unsigned timeIdx)
{
const unsigned polymerMixRegionIdx = elemCtx.problem().plmixnumRegionIndex(elemCtx, scvIdx, timeIdx);
return plyvmhCoefficients_[polymerMixRegionIdx];
}
static bool hasPlyshlog()
{
return hasPlyshlog_;
}
static bool hasShrate()
{
return hasShrate_;
}
static const Scalar shrate(unsigned pvtnumRegionIdx)
{
return shrate_[pvtnumRegionIdx];
}
/*!
* \brief Computes the shear factor
*
* Input is polymer concentration and either the water velocity or the shrate if hasShrate_ is true.
* The pvtnumRegionIdx is needed to make sure the right table is used.
*/
template <class Evaluation>
static Evaluation computeShearFactor(const Evaluation& polymerConcentration,
unsigned pvtnumRegionIdx,
const Evaluation& v0)
{
typedef Opm::MathToolbox<Evaluation> ToolboxLocal;
const auto& viscosityMultiplierTable = plyviscViscosityMultiplierTable_[pvtnumRegionIdx];
Scalar viscosityMultiplier = viscosityMultiplierTable.eval(Opm::scalarValue(polymerConcentration), /*extrapolate=*/true);
const Scalar eps = 1e-14;
// return 1.0 if the polymer has no effect on the water.
if (std::abs((viscosityMultiplier - 1.0)) < eps)
return ToolboxLocal::createConstant(v0, 1.0);
const std::vector<Scalar>& shearEffectRefLogVelocity = plyshlogShearEffectRefLogVelocity_[pvtnumRegionIdx];
auto v0AbsLog = Opm::log(Opm::abs(v0));
// return 1.0 if the velocity /sharte is smaller than the first velocity entry.
if (v0AbsLog < shearEffectRefLogVelocity[0])
return ToolboxLocal::createConstant(v0, 1.0);
// compute shear factor from input
// Z = (1 + (P - 1) * M(v)) / P
// where M(v) is computed from user input
// and P = viscosityMultiplier
const std::vector<Scalar>& shearEffectRefMultiplier = plyshlogShearEffectRefMultiplier_[pvtnumRegionIdx];
size_t numTableEntries = shearEffectRefLogVelocity.size();
assert(shearEffectRefMultiplier.size() == numTableEntries);
std::vector<Scalar> shearEffectMultiplier(numTableEntries, 1.0);
for (size_t i = 0; i < numTableEntries; ++i) {
shearEffectMultiplier[i] = (1.0 + (viscosityMultiplier - 1.0)*shearEffectRefMultiplier[i]) / viscosityMultiplier;
shearEffectMultiplier[i] = Opm::log(shearEffectMultiplier[i]);
}
// store the logarithmic velocity and logarithmic multipliers in a table for easy look up and
// linear interpolation in the logarithmic space.
TabulatedFunction logShearEffectMultiplier = TabulatedFunction(numTableEntries, shearEffectRefLogVelocity, shearEffectMultiplier, /*bool sortInputs =*/ false);
// Find sheared velocity (v) that satisfies
// F = log(v) + log (Z) - log(v0) = 0;
// Set up the function
// u = log(v)
auto F = [&logShearEffectMultiplier, &v0AbsLog](const Evaluation& u) {
return u + logShearEffectMultiplier.eval(u, true) - v0AbsLog;
};
// and its derivative
auto dF = [&logShearEffectMultiplier](const Evaluation& u) {
return 1 + logShearEffectMultiplier.evalDerivative(u, true);
};
// Solve F = 0 using Newton
// Use log(v0) as initial value for u
auto u = v0AbsLog;
bool converged = false;
for (int i = 0; i < 20; ++i) {
auto f = F(u);
auto df = dF(u);
u -= f/df;
if (std::abs(Opm::scalarValue(f)) < 1e-12) {
converged = true;
break;
}
}
if (!converged) {
throw std::runtime_error("Not able to compute shear velocity. \n");
}
// return the shear factor
return Opm::exp(logShearEffectMultiplier.eval(u, /*extrapolate=*/true));
}
const Scalar molarMass() const
{
return 0.25; // kg/mol
}
private:
static std::vector<Scalar> plyrockDeadPoreVolume_;
static std::vector<Scalar> plyrockResidualResistanceFactor_;
static std::vector<Scalar> plyrockRockDensityFactor_;
static std::vector<Scalar> plyrockAdsorbtionIndex_;
static std::vector<Scalar> plyrockMaxAdsorbtion_;
static std::vector<TabulatedFunction> plyadsAdsorbedPolymer_;
static std::vector<TabulatedFunction> plyviscViscosityMultiplierTable_;
static std::vector<Scalar> plymaxMaxConcentration_;
static std::vector<Scalar> plymixparToddLongstaff_;
static std::vector<std::vector<Scalar>> plyshlogShearEffectRefMultiplier_;
static std::vector<std::vector<Scalar>> plyshlogShearEffectRefLogVelocity_;
static std::vector<Scalar> shrate_;
static bool hasShrate_;
static bool hasPlyshlog_;
static std::vector<PlyvmhCoefficients> plyvmhCoefficients_;
static std::map<int, TabulatedTwoDFunction> plymwinjTables_;
static std::map<int, TabulatedTwoDFunction> skprwatTables_;
static std::map<int, SkprpolyTable> skprpolyTables_;
};
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyrockDeadPoreVolume_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyrockResidualResistanceFactor_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyrockRockDensityFactor_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyrockAdsorbtionIndex_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyrockMaxAdsorbtion_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::TabulatedFunction>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyadsAdsorbedPolymer_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::TabulatedFunction>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyviscViscosityMultiplierTable_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plymaxMaxConcentration_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plymixparToddLongstaff_;
template <class TypeTag, bool enablePolymerV>
std::vector<std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyshlogShearEffectRefMultiplier_;
template <class TypeTag, bool enablePolymerV>
std::vector<std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyshlogShearEffectRefLogVelocity_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::Scalar>
BlackOilPolymerModule<TypeTag, enablePolymerV>::shrate_;
template <class TypeTag, bool enablePolymerV>
bool
BlackOilPolymerModule<TypeTag, enablePolymerV>::hasShrate_;
template <class TypeTag, bool enablePolymerV>
bool
BlackOilPolymerModule<TypeTag, enablePolymerV>::hasPlyshlog_;
template <class TypeTag, bool enablePolymerV>
std::vector<typename BlackOilPolymerModule<TypeTag, enablePolymerV>::PlyvmhCoefficients>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plyvmhCoefficients_;
template <class TypeTag, bool enablePolymerV>
std::map<int, typename BlackOilPolymerModule<TypeTag, enablePolymerV>::TabulatedTwoDFunction>
BlackOilPolymerModule<TypeTag, enablePolymerV>::plymwinjTables_;