blackoilsolventmodules.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_SOLVENT_MODULE_HH
29#define EWOMS_BLACK_OIL_SOLVENT_MODULE_HH
30
31#include "blackoilproperties.hh"
32
33#include <opm/common/Exceptions.hpp>
34
39
40#include <opm/material/fluidsystems/blackoilpvt/SolventPvt.hpp>
41
42#if HAVE_ECL_INPUT
43#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
44#include <opm/input/eclipse/EclipseState/Tables/SsfnTable.hpp>
45#include <opm/input/eclipse/EclipseState/Tables/Sof2Table.hpp>
46#include <opm/input/eclipse/EclipseState/Tables/MsfnTable.hpp>
47#include <opm/input/eclipse/EclipseState/Tables/PmiscTable.hpp>
48#include <opm/input/eclipse/EclipseState/Tables/MiscTable.hpp>
49#include <opm/input/eclipse/EclipseState/Tables/SorwmisTable.hpp>
50#include <opm/input/eclipse/EclipseState/Tables/SgcwmisTable.hpp>
51#include <opm/input/eclipse/EclipseState/Tables/TlpmixpaTable.hpp>
52#endif
53
54#include <opm/material/common/Valgrind.hpp>
55
56#include <dune/common/fvector.hh>
57
58#include <string>
59
60namespace Opm {
61
67template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
69{
82 using Toolbox = MathToolbox<Evaluation>;
83 using SolventPvt = typename BlackOilSolventParams<Scalar>::SolventPvt;
84 using Co2GasPvt = typename BlackOilSolventParams<Scalar>::Co2GasPvt;
85 using H2GasPvt = typename BlackOilSolventParams<Scalar>::H2GasPvt;
86 using BrineCo2Pvt = typename BlackOilSolventParams<Scalar>::BrineCo2Pvt;
87 using BrineH2Pvt = typename BlackOilSolventParams<Scalar>::BrineH2Pvt;
88
89 using TabulatedFunction = typename BlackOilSolventParams<Scalar>::TabulatedFunction;
90
91 static constexpr unsigned solventSaturationIdx = Indices::solventSaturationIdx;
92 static constexpr unsigned contiSolventEqIdx = Indices::contiSolventEqIdx;
93 static constexpr unsigned enableSolvent = enableSolventV;
94 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
95 static constexpr unsigned numPhases = FluidSystem::numPhases;
96 static constexpr bool blackoilConserveSurfaceVolume = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
97 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
98
99
100public:
101#if HAVE_ECL_INPUT
105 static void initFromState(const EclipseState& eclState, const Schedule& schedule)
106 {
107 // some sanity checks: if solvents are enabled, the SOLVENT keyword must be
108 // present, if solvents are disabled the keyword must not be present.
109 if (enableSolvent && !eclState.runspec().phases().active(Phase::SOLVENT))
110 throw std::runtime_error("Non-trivial solvent treatment requested at compile "
111 "time, but the deck does not contain the SOLVENT keyword");
112 else if (!enableSolvent && eclState.runspec().phases().active(Phase::SOLVENT))
113 throw std::runtime_error("Solvent treatment disabled at compile time, but the deck "
114 "contains the SOLVENT keyword");
115
116 if (!eclState.runspec().phases().active(Phase::SOLVENT))
117 return; // solvent treatment is supposed to be disabled
118
119 params_.co2sol_ = eclState.runspec().co2Sol();
120 params_.h2sol_ = eclState.runspec().h2Sol();
121
122 if (isCO2Sol() && isH2Sol()) {
123 throw std::runtime_error("CO2SOL and H2SOL can not be used together");
124 }
125
126 if (isCO2Sol() || isH2Sol() ) {
127 if (isCO2Sol()) {
128 params_.co2GasPvt_.initFromState(eclState, schedule);
129 params_.brineCo2Pvt_.initFromState(eclState, schedule);
130 } else {
131 params_.h2GasPvt_.initFromState(eclState, schedule);
132 params_.brineH2Pvt_.initFromState(eclState, schedule);
133 }
134 if (eclState.getSimulationConfig().hasDISGASW()) {
135 params_.rsSolw_active_ = true;
136 }
137 } else
138 params_.solventPvt_.initFromState(eclState, schedule);
139
140 const auto& tableManager = eclState.getTableManager();
141 // initialize the objects which deal with the SSFN keyword
142 const auto& ssfnTables = tableManager.getSsfnTables();
143 unsigned numSatRegions = tableManager.getTabdims().getNumSatTables();
144 params_.setNumSatRegions(numSatRegions);
145 for (unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
146 const auto& ssfnTable = ssfnTables.template getTable<SsfnTable>(satRegionIdx);
147 params_.ssfnKrg_[satRegionIdx].setXYContainers(ssfnTable.getSolventFractionColumn(),
148 ssfnTable.getGasRelPermMultiplierColumn(),
149 /*sortInput=*/true);
150 params_.ssfnKrs_[satRegionIdx].setXYContainers(ssfnTable.getSolventFractionColumn(),
151 ssfnTable.getSolventRelPermMultiplierColumn(),
152 /*sortInput=*/true);
153 }
154
155 // initialize the objects needed for miscible solvent and oil simulations
156 params_.isMiscible_ = false;
157 if (!eclState.getTableManager().getMiscTables().empty()) {
158 params_.isMiscible_ = true;
159
160 unsigned numMiscRegions = 1;
161
162 // misicible hydrocabon relative permeability wrt water
163 const auto& sof2Tables = tableManager.getSof2Tables();
164 if (!sof2Tables.empty()) {
165 // resize the attributes of the object
166 params_.sof2Krn_.resize(numSatRegions);
167 for (unsigned satRegionIdx = 0; satRegionIdx < numSatRegions; ++ satRegionIdx) {
168 const auto& sof2Table = sof2Tables.template getTable<Sof2Table>(satRegionIdx);
169 params_.sof2Krn_[satRegionIdx].setXYContainers(sof2Table.getSoColumn(),
170 sof2Table.getKroColumn(),
171 /*sortInput=*/true);
172 }
173 }
174 else if(eclState.runspec().phases().active(Phase::OIL))
175 throw std::runtime_error("SOF2 must be specified in MISCIBLE (SOLVENT and OIL) runs\n");
176
177 const auto& miscTables = tableManager.getMiscTables();
178 if (!miscTables.empty()) {
179 assert(numMiscRegions == miscTables.size());
180
181 // resize the attributes of the object
182 params_.misc_.resize(numMiscRegions);
183 for (unsigned miscRegionIdx = 0; miscRegionIdx < numMiscRegions; ++miscRegionIdx) {
184 const auto& miscTable = miscTables.template getTable<MiscTable>(miscRegionIdx);
185
186 // solventFraction = Ss / (Ss + Sg);
187 const auto& solventFraction = miscTable.getSolventFractionColumn();
188 const auto& misc = miscTable.getMiscibilityColumn();
189 params_.misc_[miscRegionIdx].setXYContainers(solventFraction, misc);
190 }
191 }
192 else
193 throw std::runtime_error("MISC must be specified in MISCIBLE (SOLVENT) runs\n");
194
195 // resize the attributes of the object
196 params_.pmisc_.resize(numMiscRegions);
197 const auto& pmiscTables = tableManager.getPmiscTables();
198 if (!pmiscTables.empty()) {
199 assert(numMiscRegions == pmiscTables.size());
200
201 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
202 const auto& pmiscTable = pmiscTables.template getTable<PmiscTable>(regionIdx);
203
204 // Copy data
205 const auto& po = pmiscTable.getOilPhasePressureColumn();
206 const auto& pmisc = pmiscTable.getMiscibilityColumn();
207
208 params_.pmisc_[regionIdx].setXYContainers(po, pmisc);
209 }
210 }
211 else {
212 std::vector<double> x = {0.0,1.0e20};
213 std::vector<double> y = {1.0,1.0};
214 TabulatedFunction constant = TabulatedFunction(2, x, y);
215 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
216 params_.pmisc_[regionIdx] = constant;
217 }
218 }
219
220 // miscible relative permeability multipleiers
221 params_.msfnKrsg_.resize(numSatRegions);
222 params_.msfnKro_.resize(numSatRegions);
223 const auto& msfnTables = tableManager.getMsfnTables();
224 if (!msfnTables.empty()) {
225 assert(numSatRegions == msfnTables.size());
226
227 for (unsigned regionIdx = 0; regionIdx < numSatRegions; ++regionIdx) {
228 const MsfnTable& msfnTable = msfnTables.template getTable<MsfnTable>(regionIdx);
229
230 // Copy data
231 // Ssg = Ss + Sg;
232 const auto& Ssg = msfnTable.getGasPhaseFractionColumn();
233 const auto& krsg = msfnTable.getGasSolventRelpermMultiplierColumn();
234 const auto& kro = msfnTable.getOilRelpermMultiplierColumn();
235
236 params_.msfnKrsg_[regionIdx].setXYContainers(Ssg, krsg);
237 params_.msfnKro_[regionIdx].setXYContainers(Ssg, kro);
238 }
239 }
240 else {
241 std::vector<double> x = {0.0,1.0};
242 std::vector<double> y = {1.0,0.0};
243 TabulatedFunction unit = TabulatedFunction(2, x, x);
244 TabulatedFunction invUnit = TabulatedFunction(2, x, y);
245
246 for (unsigned regionIdx = 0; regionIdx < numSatRegions; ++regionIdx) {
247 params_.setMsfn(regionIdx, unit, invUnit);
248 }
249 }
250 // resize the attributes of the object
251 params_.sorwmis_.resize(numMiscRegions);
252 const auto& sorwmisTables = tableManager.getSorwmisTables();
253 if (!sorwmisTables.empty()) {
254 assert(numMiscRegions == sorwmisTables.size());
255
256 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
257 const auto& sorwmisTable = sorwmisTables.template getTable<SorwmisTable>(regionIdx);
258
259 // Copy data
260 const auto& sw = sorwmisTable.getWaterSaturationColumn();
261 const auto& sorwmis = sorwmisTable.getMiscibleResidualOilColumn();
262
263 params_.sorwmis_[regionIdx].setXYContainers(sw, sorwmis);
264 }
265 }
266 else {
267 // default
268 std::vector<double> x = {0.0,1.0};
269 std::vector<double> y = {0.0,0.0};
270 TabulatedFunction zero = TabulatedFunction(2, x, y);
271 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
272 params_.sorwmis_[regionIdx] = zero;
273 }
274 }
275
276 // resize the attributes of the object
277 params_.sgcwmis_.resize(numMiscRegions);
278 const auto& sgcwmisTables = tableManager.getSgcwmisTables();
279 if (!sgcwmisTables.empty()) {
280 assert(numMiscRegions == sgcwmisTables.size());
281
282 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
283 const auto& sgcwmisTable = sgcwmisTables.template getTable<SgcwmisTable>(regionIdx);
284
285 // Copy data
286 const auto& sw = sgcwmisTable.getWaterSaturationColumn();
287 const auto& sgcwmis = sgcwmisTable.getMiscibleResidualGasColumn();
288
289 params_.sgcwmis_[regionIdx].setXYContainers(sw, sgcwmis);
290 }
291 }
292 else {
293 // default
294 std::vector<double> x = {0.0,1.0};
295 std::vector<double> y = {0.0,0.0};
296 TabulatedFunction zero = TabulatedFunction(2, x, y);
297 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx)
298 params_.sgcwmis_[regionIdx] = zero;
299 }
300
301 const auto& tlmixpar = eclState.getTableManager().getTLMixpar();
302 if (!tlmixpar.empty()) {
303 // resize the attributes of the object
304 params_.tlMixParamViscosity_.resize(numMiscRegions);
305 params_.tlMixParamDensity_.resize(numMiscRegions);
306
307 assert(numMiscRegions == tlmixpar.size());
308 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
309 const auto& tlp = tlmixpar[regionIdx];
310 params_.tlMixParamViscosity_[regionIdx] = tlp.viscosity_parameter;
311 params_.tlMixParamDensity_[regionIdx] = tlp.density_parameter;
312 }
313 }
314 else
315 throw std::runtime_error("TLMIXPAR must be specified in MISCIBLE (SOLVENT) runs\n");
316
317 // resize the attributes of the object
318 params_.tlPMixTable_.resize(numMiscRegions);
319 if (!eclState.getTableManager().getTlpmixpaTables().empty()) {
320 const auto& tlpmixparTables = tableManager.getTlpmixpaTables();
321 if (!tlpmixparTables.empty()) {
322 assert(numMiscRegions == tlpmixparTables.size());
323 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx) {
324 const auto& tlpmixparTable = tlpmixparTables.template getTable<TlpmixpaTable>(regionIdx);
325
326 // Copy data
327 const auto& po = tlpmixparTable.getOilPhasePressureColumn();
328 const auto& tlpmixpa = tlpmixparTable.getMiscibilityColumn();
329
330 params_.tlPMixTable_[regionIdx].setXYContainers(po, tlpmixpa);
331 }
332 }
333 else {
334 // if empty keyword. Try to use the pmisc table as default.
335 if (params_.pmisc_.size() > 0)
336 params_.tlPMixTable_ = params_.pmisc_;
337 else
338 throw std::invalid_argument("If the pressure dependent TL values in "
339 "TLPMIXPA is defaulted (no entries), then "
340 "the PMISC tables must be specified.");
341 }
342 }
343 else {
344 // default
345 std::vector<double> x = {0.0,1.0e20};
346 std::vector<double> y = {1.0,1.0};
347 TabulatedFunction ones = TabulatedFunction(2, x, y);
348 for (unsigned regionIdx = 0; regionIdx < numMiscRegions; ++regionIdx)
349 params_.tlPMixTable_[regionIdx] = ones;
350 }
351 }
352 }
353#endif
354
358 static void setSolventPvt(const SolventPvt& value)
359 { params_.solventPvt_ = value; }
360
361
362 static void setIsMiscible(const bool isMiscible)
363 { params_.isMiscible_ = isMiscible; }
364
368 static void registerParameters()
369 {
370 if constexpr (enableSolvent)
372 }
373
377 static void registerOutputModules(Model& model,
378 Simulator& simulator)
379 {
380 if constexpr (enableSolvent)
381 model.addOutputModule(new VtkBlackOilSolventModule<TypeTag>(simulator));
382 }
383
384 static bool primaryVarApplies(unsigned pvIdx)
385 {
386 if constexpr (enableSolvent)
387 return pvIdx == solventSaturationIdx;
388 else
389 return false;
390 }
391
392 static std::string primaryVarName([[maybe_unused]] unsigned pvIdx)
393 {
394 assert(primaryVarApplies(pvIdx));
395
396 return "saturation_solvent";
397 }
398
399 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
400 {
401 assert(primaryVarApplies(pvIdx));
402
403 // TODO: it may be beneficial to chose this differently.
404 return static_cast<Scalar>(1.0);
405 }
406
407 static bool eqApplies(unsigned eqIdx)
408 {
409 if constexpr (enableSolvent)
410 return eqIdx == contiSolventEqIdx;
411 else
412 return false;
413 }
414
415 static std::string eqName([[maybe_unused]] unsigned eqIdx)
416 {
417 assert(eqApplies(eqIdx));
418
419 return "conti^solvent";
420 }
421
422 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
423 {
424 assert(eqApplies(eqIdx));
425
426 // TODO: it may be beneficial to chose this differently.
427 return static_cast<Scalar>(1.0);
428 }
429
430 template <class LhsEval>
431 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
432 const IntensiveQuantities& intQuants)
433 {
434 if constexpr (enableSolvent) {
435 if constexpr (blackoilConserveSurfaceVolume) {
436 storage[contiSolventEqIdx] +=
437 Toolbox::template decay<LhsEval>(intQuants.porosity())
438 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
439 * Toolbox::template decay<LhsEval>(intQuants.solventInverseFormationVolumeFactor());
440 if (isSolubleInWater()) {
441 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
442 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
443 * Toolbox::template decay<LhsEval>(intQuants.fluidState().invB(waterPhaseIdx))
444 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
445 }
446 }
447 else {
448 storage[contiSolventEqIdx] +=
449 Toolbox::template decay<LhsEval>(intQuants.porosity())
450 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
451 * Toolbox::template decay<LhsEval>(intQuants.solventDensity());
452 if (isSolubleInWater()) {
453 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
454 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
455 * Toolbox::template decay<LhsEval>(intQuants.fluidState().density(waterPhaseIdx))
456 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
457 }
458
459 }
460 }
461 }
462
463 static void computeFlux([[maybe_unused]] RateVector& flux,
464 [[maybe_unused]] const ElementContext& elemCtx,
465 [[maybe_unused]] unsigned scvfIdx,
466 [[maybe_unused]] unsigned timeIdx)
467
468 {
469 if constexpr (enableSolvent) {
470 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
471
472 unsigned upIdx = extQuants.solventUpstreamIndex();
473 unsigned inIdx = extQuants.interiorIndex();
474 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
475
476 if constexpr (blackoilConserveSurfaceVolume) {
477 if (upIdx == inIdx)
478 flux[contiSolventEqIdx] =
479 extQuants.solventVolumeFlux()
480 *up.solventInverseFormationVolumeFactor();
481 else
482 flux[contiSolventEqIdx] =
483 extQuants.solventVolumeFlux()
484 *decay<Scalar>(up.solventInverseFormationVolumeFactor());
485
486
487 if (isSolubleInWater()) {
488 if (upIdx == inIdx)
489 flux[contiSolventEqIdx] +=
490 extQuants.volumeFlux(waterPhaseIdx)
491 * up.fluidState().invB(waterPhaseIdx)
492 * up.rsSolw();
493 else
494 flux[contiSolventEqIdx] +=
495 extQuants.volumeFlux(waterPhaseIdx)
496 *decay<Scalar>(up.fluidState().invB(waterPhaseIdx))
497 *decay<Scalar>(up.rsSolw());
498 }
499 }
500 else {
501 if (upIdx == inIdx)
502 flux[contiSolventEqIdx] =
503 extQuants.solventVolumeFlux()
504 *up.solventDensity();
505 else
506 flux[contiSolventEqIdx] =
507 extQuants.solventVolumeFlux()
508 *decay<Scalar>(up.solventDensity());
509
510
511 if (isSolubleInWater()) {
512 if (upIdx == inIdx)
513 flux[contiSolventEqIdx] +=
514 extQuants.volumeFlux(waterPhaseIdx)
515 * up.fluidState().density(waterPhaseIdx)
516 * up.rsSolw();
517 else
518 flux[contiSolventEqIdx] +=
519 extQuants.volumeFlux(waterPhaseIdx)
520 *decay<Scalar>(up.fluidState().density(waterPhaseIdx))
521 *decay<Scalar>(up.rsSolw());
522 }
523 }
524 }
525 }
526
530 static void assignPrimaryVars(PrimaryVariables& priVars,
531 Scalar solventSaturation,
532 Scalar solventRsw)
533 {
534 if constexpr (!enableSolvent) {
535 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Disabled);
536 return;
537 }
538 // Determine the meaning of the solvent primary variables
539 if (solventSaturation > 0 || !isSolubleInWater()) {
540 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Ss);
541 priVars[solventSaturationIdx] = solventSaturation;
542 } else {
543 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Rsolw);
544 priVars[solventSaturationIdx] = solventRsw;
545 }
546 }
547
551 static void updatePrimaryVars(PrimaryVariables& newPv,
552 const PrimaryVariables& oldPv,
553 const EqVector& delta)
554 {
555 if constexpr (enableSolvent)
556 // do a plain unchopped Newton update
557 newPv[solventSaturationIdx] = oldPv[solventSaturationIdx] - delta[solventSaturationIdx];
558 }
559
563 static Scalar computeUpdateError(const PrimaryVariables&,
564 const EqVector&)
565 {
566 // do not consider consider the cange of solvent primary variables for
567 // convergence
568 // TODO: maybe this should be changed
569 return static_cast<Scalar>(0.0);
570 }
571
575 static Scalar computeResidualError(const EqVector& resid)
576 {
577 // do not weight the residual of solvents when it comes to convergence
578 return std::abs(Toolbox::scalarValue(resid[contiSolventEqIdx]));
579 }
580
581 template <class DofEntity>
582 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
583 {
584 if constexpr (enableSolvent) {
585 unsigned dofIdx = model.dofMapper().index(dof);
586
587 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
588 outstream << priVars[solventSaturationIdx];
589 }
590 }
591
592 template <class DofEntity>
593 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
594 {
595 if constexpr (enableSolvent) {
596 unsigned dofIdx = model.dofMapper().index(dof);
597
598 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
599 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
600
601 instream >> priVars0[solventSaturationIdx];
602
603 // set the primary variables for the beginning of the current time step.
604 priVars1 = priVars0[solventSaturationIdx];
605 }
606 }
607
608 static const SolventPvt& solventPvt()
609 {
610 return params_.solventPvt_;
611 }
612
613
614 static const Co2GasPvt& co2GasPvt()
615 {
616 return params_.co2GasPvt_;
617 }
618
619 static const H2GasPvt& h2GasPvt()
620 {
621 return params_.h2GasPvt_;
622 }
623
624 static const BrineCo2Pvt& brineCo2Pvt()
625 {
626 return params_.brineCo2Pvt_;
627 }
628
629 static const BrineH2Pvt& brineH2Pvt()
630 {
631 return params_.brineH2Pvt_;
632 }
633
634 static const TabulatedFunction& ssfnKrg(const ElementContext& elemCtx,
635 unsigned scvIdx,
636 unsigned timeIdx)
637 {
638 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
639 return params_.ssfnKrg_[satnumRegionIdx];
640 }
641
642 static const TabulatedFunction& ssfnKrs(const ElementContext& elemCtx,
643 unsigned scvIdx,
644 unsigned timeIdx)
645 {
646 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
647 return params_.ssfnKrs_[satnumRegionIdx];
648 }
649
650 static const TabulatedFunction& sof2Krn(const ElementContext& elemCtx,
651 unsigned scvIdx,
652 unsigned timeIdx)
653 {
654 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
655 return params_.sof2Krn_[satnumRegionIdx];
656 }
657
658 static const TabulatedFunction& misc(const ElementContext& elemCtx,
659 unsigned scvIdx,
660 unsigned timeIdx)
661 {
662 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
663 return params_.misc_[miscnumRegionIdx];
664 }
665
666 static const TabulatedFunction& pmisc(const ElementContext& elemCtx,
667 unsigned scvIdx,
668 unsigned timeIdx)
669 {
670 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
671 return params_.pmisc_[miscnumRegionIdx];
672 }
673
674 static const TabulatedFunction& msfnKrsg(const ElementContext& elemCtx,
675 unsigned scvIdx,
676 unsigned timeIdx)
677 {
678 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
679 return params_.msfnKrsg_[satnumRegionIdx];
680 }
681
682 static const TabulatedFunction& msfnKro(const ElementContext& elemCtx,
683 unsigned scvIdx,
684 unsigned timeIdx)
685 {
686 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
687 return params_.msfnKro_[satnumRegionIdx];
688 }
689
690 static const TabulatedFunction& sorwmis(const ElementContext& elemCtx,
691 unsigned scvIdx,
692 unsigned timeIdx)
693 {
694 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
695 return params_.sorwmis_[miscnumRegionIdx];
696 }
697
698 static const TabulatedFunction& sgcwmis(const ElementContext& elemCtx,
699 unsigned scvIdx,
700 unsigned timeIdx)
701 {
702 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
703 return params_.sgcwmis_[miscnumRegionIdx];
704 }
705
706 static const TabulatedFunction& tlPMixTable(const ElementContext& elemCtx,
707 unsigned scvIdx,
708 unsigned timeIdx)
709 {
710 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
711 return params_.tlPMixTable_[miscnumRegionIdx];
712 }
713
714 static const Scalar& tlMixParamViscosity(const ElementContext& elemCtx,
715 unsigned scvIdx,
716 unsigned timeIdx)
717 {
718 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
719 return params_.tlMixParamViscosity_[miscnumRegionIdx];
720 }
721
722 static const Scalar& tlMixParamDensity(const ElementContext& elemCtx,
723 unsigned scvIdx,
724 unsigned timeIdx)
725 {
726 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
727 return params_.tlMixParamDensity_[miscnumRegionIdx];
728 }
729
730 static bool isMiscible()
731 {
732 return params_.isMiscible_;
733 }
734
735 template <class Value>
736 static const Value solubilityLimit(unsigned pvtIdx, const Value& temperature, const Value& pressure, const Value& saltConcentration)
737 {
738 if (!isSolubleInWater())
739 return 0.0;
740
741 assert(isCO2Sol() || isH2Sol());
742 if (isCO2Sol())
743 return brineCo2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
744 else
745 return brineH2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
746 }
747
748 static bool isSolubleInWater()
749 {
750 return params_.rsSolw_active_;
751 }
752
753 static bool isCO2Sol()
754 {
755 return params_.co2sol_;
756 }
757
758 static bool isH2Sol()
759 {
760 return params_.h2sol_;
761 }
762
763private:
764 static BlackOilSolventParams<Scalar> params_; // the krg(Fs) column of the SSFN table
765};
766
767template <class TypeTag, bool enableSolventV>
768BlackOilSolventParams<typename BlackOilSolventModule<TypeTag, enableSolventV>::Scalar>
769BlackOilSolventModule<TypeTag, enableSolventV>::params_;
770
778template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
780{
782
790
792
793 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
794 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
795 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
796 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
797 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
798 static constexpr double cutOff = 1e-12;
799
800
801public:
808 void solventPreSatFuncUpdate_(const ElementContext& elemCtx,
809 unsigned dofIdx,
810 unsigned timeIdx)
811 {
812 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
813
814 auto& fs = asImp_().fluidState_;
815 solventSaturation_ = 0.0;
816 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
817 solventSaturation_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
818 }
819
820 hydrocarbonSaturation_ = fs.saturation(gasPhaseIdx);
821
822 // apply a cut-off. Don't waste calculations if no solvent
823 if (solventSaturation().value() < cutOff)
824 return;
825
826 // make the saturation of the gas phase which is used by the saturation functions
827 // the sum of the solvent "saturation" and the saturation the hydrocarbon gas.
828 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_ + solventSaturation_);
829 }
830
838 void solventPostSatFuncUpdate_(const ElementContext& elemCtx,
839 unsigned dofIdx,
840 unsigned timeIdx)
841 {
842 // revert the gas "saturation" of the fluid state back to the saturation of the
843 // hydrocarbon gas.
844 auto& fs = asImp_().fluidState_;
845 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_);
846
847
848 // update rsSolw. This needs to be done after the pressure is defined in the fluid state.
849 rsSolw_ = 0.0;
850 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
851 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
852 rsSolw_ = SolventModule::solubilityLimit(asImp_().pvtRegionIndex(), fs.temperature(waterPhaseIdx), fs.pressure(waterPhaseIdx), fs.saltConcentration());
853 } else if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Rsolw) {
854 rsSolw_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
855 }
856
857 solventMobility_ = 0.0;
858
859 // apply a cut-off. Don't waste calculations if no solvent
860 if (solventSaturation().value() < cutOff)
861 return;
862
863 // Pressure effects on capillary pressure miscibility
865 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
866 const Evaluation pmisc = SolventModule::pmisc(elemCtx, dofIdx, timeIdx).eval(p, /*extrapolate=*/true);
867 const Evaluation& pgImisc = fs.pressure(gasPhaseIdx);
868
869 // compute capillary pressure for miscible fluid
870 const auto& problem = elemCtx.problem();
871 Evaluation pgMisc = 0.0;
872 std::array<Evaluation, numPhases> pC;
873 const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
874 MaterialLaw::capillaryPressures(pC, materialParams, fs);
875
876 //oil is the reference phase for pressure
877 const auto linearizationType = elemCtx.linearizationType();
878 if (priVars.primaryVarsMeaningPressure() == PrimaryVariables::PressureMeaning::Pg)
879 pgMisc = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
880 else {
881 const Evaluation& po = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
882 pgMisc = po + (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
883 }
884
885 fs.setPressure(gasPhaseIdx, pmisc * pgMisc + (1.0 - pmisc) * pgImisc);
886 }
887
888
889 Evaluation gasSolventSat = hydrocarbonSaturation_ + solventSaturation_;
890
891 if (gasSolventSat.value() < cutOff) // avoid division by zero
892 return;
893
894 Evaluation Fhydgas = hydrocarbonSaturation_/gasSolventSat;
895 Evaluation Fsolgas = solventSaturation_/gasSolventSat;
896
897 // account for miscibility of oil and solvent
898 if (SolventModule::isMiscible() && FluidSystem::phaseIsActive(oilPhaseIdx)) {
899 const auto& misc = SolventModule::misc(elemCtx, dofIdx, timeIdx);
900 const auto& pmisc = SolventModule::pmisc(elemCtx, dofIdx, timeIdx);
901 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
902 const Evaluation miscibility = misc.eval(Fsolgas, /*extrapolate=*/true) * pmisc.eval(p, /*extrapolate=*/true);
903
904 // TODO adjust endpoints of sn and ssg
905 unsigned cellIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
906 const auto& materialLawManager = elemCtx.problem().materialLawManager();
907 const auto& scaledDrainageInfo =
908 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
909
910 const Scalar& sogcr = scaledDrainageInfo.Sogcr;
911 Evaluation sor = sogcr;
912 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
913 const Evaluation& sw = fs.saturation(waterPhaseIdx);
914 const auto& sorwmis = SolventModule::sorwmis(elemCtx, dofIdx, timeIdx);
915 sor = miscibility * sorwmis.eval(sw, /*extrapolate=*/true) + (1.0 - miscibility) * sogcr;
916 }
917 const Scalar& sgcr = scaledDrainageInfo.Sgcr;
918 Evaluation sgc = sgcr;
919 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
920 const Evaluation& sw = fs.saturation(waterPhaseIdx);
921 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, dofIdx, timeIdx);
922 sgc = miscibility * sgcwmis.eval(sw, /*extrapolate=*/true) + (1.0 - miscibility) * sgcr;
923 }
924
925 Evaluation oilGasSolventSat = gasSolventSat;
926 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
927 oilGasSolventSat += fs.saturation(oilPhaseIdx);
928 }
929 const Evaluation zero = 0.0;
930 const Evaluation oilGasSolventEffSat = std::max(oilGasSolventSat - sor - sgc, zero);
931
932 Evaluation F_totalGas = 0.0;
933 if (oilGasSolventEffSat.value() > cutOff) {
934 const Evaluation gasSolventEffSat = std::max(gasSolventSat - sgc, zero);
935 F_totalGas = gasSolventEffSat / oilGasSolventEffSat;
936 }
937 const auto& msfnKro = SolventModule::msfnKro(elemCtx, dofIdx, timeIdx);
938 const auto& msfnKrsg = SolventModule::msfnKrsg(elemCtx, dofIdx, timeIdx);
939 const auto& sof2Krn = SolventModule::sof2Krn(elemCtx, dofIdx, timeIdx);
940
941 const Evaluation mkrgt = msfnKrsg.eval(F_totalGas, /*extrapolate=*/true) * sof2Krn.eval(oilGasSolventSat, /*extrapolate=*/true);
942 const Evaluation mkro = msfnKro.eval(F_totalGas, /*extrapolate=*/true) * sof2Krn.eval(oilGasSolventSat, /*extrapolate=*/true);
943
944 Evaluation& kro = asImp_().mobility_[oilPhaseIdx];
945 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
946
947 // combine immiscible and miscible part of the relperm
948 krg *= (1.0 - miscibility);
949 krg += miscibility * mkrgt;
950 kro *= (1.0 - miscibility);
951 kro += miscibility * mkro;
952 }
953
954 // compute the mobility of the solvent "phase" and modify the gas phase
955 const auto& ssfnKrg = SolventModule::ssfnKrg(elemCtx, dofIdx, timeIdx);
956 const auto& ssfnKrs = SolventModule::ssfnKrs(elemCtx, dofIdx, timeIdx);
957
958 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
959 solventMobility_ = krg * ssfnKrs.eval(Fsolgas, /*extrapolate=*/true);
960 krg *= ssfnKrg.eval(Fhydgas, /*extrapolate=*/true);
961
962 }
963
970 void solventPvtUpdate_(const ElementContext& elemCtx,
971 unsigned scvIdx,
972 unsigned timeIdx)
973 {
974 const auto& iq = asImp_();
975 unsigned pvtRegionIdx = iq.pvtRegionIndex();
976 const Evaluation& T = iq.fluidState().temperature(gasPhaseIdx);
977 const Evaluation& p = iq.fluidState().pressure(gasPhaseIdx);
978
979 const Evaluation rv = 0.0;
980 const Evaluation rvw = 0.0;
983 const auto& co2gasPvt = SolventModule::co2GasPvt();
984 solventInvFormationVolumeFactor_ = co2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
985 solventRefDensity_ = co2gasPvt.gasReferenceDensity(pvtRegionIdx);
986 solventViscosity_ = co2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
987
988 const auto& brineCo2Pvt = SolventModule::brineCo2Pvt();
989 auto& fs = asImp_().fluidState_;
990 const auto& bw = brineCo2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
991
992 const auto denw = bw*brineCo2Pvt.waterReferenceDensity(pvtRegionIdx)
993 + rsSolw()*bw*brineCo2Pvt.gasReferenceDensity(pvtRegionIdx);
994 fs.setDensity(waterPhaseIdx, denw);
995 fs.setInvB(waterPhaseIdx, bw);
996 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
997 const auto& muWat = fs.viscosity(waterPhaseIdx);
998 const auto& muWatEff = brineCo2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
999 mobw *= muWat / muWatEff;
1000 } else {
1001 const auto& h2gasPvt = SolventModule::h2GasPvt();
1002 solventInvFormationVolumeFactor_ = h2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
1003 solventRefDensity_ = h2gasPvt.gasReferenceDensity(pvtRegionIdx);
1004 solventViscosity_ = h2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
1005
1006 const auto& brineH2Pvt = SolventModule::brineH2Pvt();
1007 auto& fs = asImp_().fluidState_;
1008 const auto& bw = brineH2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
1009
1010 const auto denw = bw*brineH2Pvt.waterReferenceDensity(pvtRegionIdx)
1011 + rsSolw()*bw*brineH2Pvt.gasReferenceDensity(pvtRegionIdx);
1012 fs.setDensity(waterPhaseIdx, denw);
1013 fs.setInvB(waterPhaseIdx, bw);
1014 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
1015 const auto& muWat = fs.viscosity(waterPhaseIdx);
1016 const auto& muWatEff = brineH2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
1017 mobw *= muWat / muWatEff;
1018 }
1019 } else {
1020 const auto& solventPvt = SolventModule::solventPvt();
1021 solventInvFormationVolumeFactor_ = solventPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p);
1022 solventRefDensity_ = solventPvt.referenceDensity(pvtRegionIdx);
1023 solventViscosity_ = solventPvt.viscosity(pvtRegionIdx, T, p);
1024 }
1025
1027 effectiveProperties(elemCtx, scvIdx, timeIdx);
1028
1030
1031
1032 }
1033
1034 const Evaluation& solventSaturation() const
1035 { return solventSaturation_; }
1036
1037 const Evaluation& rsSolw() const
1038 { return rsSolw_; }
1039
1040 const Evaluation& solventDensity() const
1041 { return solventDensity_; }
1042
1043 const Evaluation& solventViscosity() const
1044 { return solventViscosity_; }
1045
1046 const Evaluation& solventMobility() const
1047 { return solventMobility_; }
1048
1049 const Evaluation& solventInverseFormationVolumeFactor() const
1051
1052 // This could be stored pr pvtRegion instead
1053 const Scalar& solventRefDensity() const
1054 { return solventRefDensity_; }
1055
1056private:
1057 // Computes the effective properties based on
1058 // Todd-Longstaff mixing model.
1059 void effectiveProperties(const ElementContext& elemCtx,
1060 unsigned scvIdx,
1061 unsigned timeIdx)
1062 {
1064 return;
1065
1066 // Don't waste calculations if no solvent
1067 // Apply a cut-off for small and negative solvent saturations
1068 if (solventSaturation() < cutOff)
1069 return;
1070
1071 // We plan to update the fluidstate with the effective
1072 // properties
1073 auto& fs = asImp_().fluidState_;
1074
1075 // Compute effective saturations
1076 Evaluation oilEffSat = 0.0;
1077 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1078 oilEffSat = fs.saturation(oilPhaseIdx);
1079 }
1080 Evaluation gasEffSat = fs.saturation(gasPhaseIdx);
1081 Evaluation solventEffSat = solventSaturation();
1082 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
1083 const auto& sorwmis = SolventModule::sorwmis(elemCtx, scvIdx, timeIdx);
1084 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, scvIdx, timeIdx);
1085 const Evaluation zero = 0.0;
1086 const Evaluation& sw = fs.saturation(waterPhaseIdx);
1087 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1088 oilEffSat = std::max(oilEffSat - sorwmis.eval(sw, /*extrapolate=*/true), zero);
1089 }
1090 gasEffSat = std::max(gasEffSat - sgcwmis.eval(sw, /*extrapolate=*/true), zero);
1091 solventEffSat = std::max(solventEffSat - sgcwmis.eval(sw, /*extrapolate=*/true), zero);
1092 }
1093 const Evaluation oilGasSolventEffSat = oilEffSat + gasEffSat + solventEffSat;
1094 const Evaluation oilSolventEffSat = oilEffSat + solventEffSat;
1095 const Evaluation solventGasEffSat = solventEffSat + gasEffSat;
1096
1097 // Compute effective viscosities
1098
1099 // Mixing parameter for viscosity
1100 // The pressureMixingParameter represent the miscibility of the solvent while the mixingParameterViscosity the effect of the porous media.
1101 // The pressureMixingParameter is not implemented in ecl100.
1102 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
1103 // account for pressure effects
1104 const auto& pmiscTable = SolventModule::pmisc(elemCtx, scvIdx, timeIdx);
1105 const Evaluation pmisc = pmiscTable.eval(p, /*extrapolate=*/true);
1106 const auto& tlPMixTable = SolventModule::tlPMixTable(elemCtx, scvIdx, timeIdx);
1107 const Evaluation tlMixParamMu = SolventModule::tlMixParamViscosity(elemCtx, scvIdx, timeIdx) * tlPMixTable.eval(p, /*extrapolate=*/true);
1108
1109 const Evaluation& muGas = fs.viscosity(gasPhaseIdx);
1110 const Evaluation& muSolvent = solventViscosity_;
1111
1112 assert(muGas.value() > 0);
1113 assert(muSolvent.value() > 0);
1114 const Evaluation muGasPow = pow(muGas, 0.25);
1115 const Evaluation muSolventPow = pow(muSolvent, 0.25);
1116
1117 Evaluation muMixSolventGas = muGas;
1118 if (solventGasEffSat > cutOff)
1119 muMixSolventGas *= muSolvent / pow(((gasEffSat / solventGasEffSat) * muSolventPow) + ((solventEffSat / solventGasEffSat) * muGasPow) , 4.0);
1120
1121 Evaluation muOil = 1.0;
1122 Evaluation muOilPow = 1.0;
1123 Evaluation muMixOilSolvent = 1.0;
1124 Evaluation muOilEff = 1.0;
1125 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1126 muOil = fs.viscosity(oilPhaseIdx);
1127 assert(muOil.value() > 0);
1128 muOilPow = pow(muOil, 0.25);
1129 muMixOilSolvent = muOil;
1130 if (oilSolventEffSat > cutOff)
1131 muMixOilSolvent *= muSolvent / pow(((oilEffSat / oilSolventEffSat) * muSolventPow) + ((solventEffSat / oilSolventEffSat) * muOilPow) , 4.0);
1132
1133 muOilEff = pow(muOil,1.0 - tlMixParamMu) * pow(muMixOilSolvent, tlMixParamMu);
1134 }
1135 Evaluation muMixSolventGasOil = muOil;
1136 if (oilGasSolventEffSat > cutOff)
1137 muMixSolventGasOil *= muSolvent * muGas / pow(((oilEffSat / oilGasSolventEffSat) * muSolventPow * muGasPow)
1138 + ((solventEffSat / oilGasSolventEffSat) * muOilPow * muGasPow) + ((gasEffSat / oilGasSolventEffSat) * muSolventPow * muOilPow), 4.0);
1139
1140 Evaluation muGasEff = pow(muGas,1.0 - tlMixParamMu) * pow(muMixSolventGas, tlMixParamMu);
1141 Evaluation muSolventEff = pow(muSolvent,1.0 - tlMixParamMu) * pow(muMixSolventGasOil, tlMixParamMu);
1142
1143 // Compute effective densities
1144 const Evaluation& rhoGas = fs.density(gasPhaseIdx);
1145 Evaluation rhoOil = 0.0;
1146 if (FluidSystem::phaseIsActive(oilPhaseIdx))
1147 rhoOil = fs.density(oilPhaseIdx);
1148
1149 const Evaluation& rhoSolvent = solventDensity_;
1150
1151 // Mixing parameter for density
1152 // The pressureMixingParameter represent the miscibility of the solvent while the mixingParameterDenisty the effect of the porous media.
1153 // The pressureMixingParameter is not implemented in ecl100.
1154 const Evaluation tlMixParamRho = SolventModule::tlMixParamDensity(elemCtx, scvIdx, timeIdx) * tlPMixTable.eval(p, /*extrapolate=*/true);
1155
1156 // compute effective viscosities for density calculations. These have to
1157 // be recomputed as a different mixing parameter may be used.
1158 const Evaluation muOilEffPow = pow(pow(muOil, 1.0 - tlMixParamRho) * pow(muMixOilSolvent, tlMixParamRho), 0.25);
1159 const Evaluation muGasEffPow = pow(pow(muGas, 1.0 - tlMixParamRho) * pow(muMixSolventGas, tlMixParamRho), 0.25);
1160 const Evaluation muSolventEffPow = pow(pow(muSolvent, 1.0 - tlMixParamRho) * pow(muMixSolventGasOil, tlMixParamRho), 0.25);
1161
1162 const Evaluation oilGasEffSaturation = oilEffSat + gasEffSat;
1163 Evaluation sof = 0.0;
1164 Evaluation sgf = 0.0;
1165 if (oilGasEffSaturation.value() > cutOff) {
1166 sof = oilEffSat / oilGasEffSaturation;
1167 sgf = gasEffSat / oilGasEffSaturation;
1168 }
1169
1170 const Evaluation muSolventOilGasPow = muSolventPow * ((sgf * muOilPow) + (sof * muGasPow));
1171
1172 Evaluation rhoMixSolventGasOil = 0.0;
1173 if (oilGasSolventEffSat.value() > cutOff)
1174 rhoMixSolventGasOil = (rhoOil * oilEffSat / oilGasSolventEffSat) + (rhoGas * gasEffSat / oilGasSolventEffSat) + (rhoSolvent * solventEffSat / oilGasSolventEffSat);
1175
1176 Evaluation rhoGasEff = 0.0;
1177 if (std::abs(muSolventPow.value() - muGasPow.value()) < cutOff)
1178 rhoGasEff = ((1.0 - tlMixParamRho) * rhoGas) + (tlMixParamRho * rhoMixSolventGasOil);
1179 else {
1180 const Evaluation solventGasEffFraction = (muGasPow * (muSolventPow - muGasEffPow)) / (muGasEffPow * (muSolventPow - muGasPow));
1181 rhoGasEff = (rhoGas * solventGasEffFraction) + (rhoSolvent * (1.0 - solventGasEffFraction));
1182 }
1183
1184 Evaluation rhoSolventEff = 0.0;
1185 if (std::abs((muSolventOilGasPow.value() - (muOilPow.value() * muGasPow.value()))) < cutOff)
1186 rhoSolventEff = ((1.0 - tlMixParamRho) * rhoSolvent) + (tlMixParamRho * rhoMixSolventGasOil);
1187 else {
1188 const Evaluation sfraction_se = (muSolventOilGasPow - (muOilPow * muGasPow * muSolventPow / muSolventEffPow)) / (muSolventOilGasPow - (muOilPow * muGasPow));
1189 rhoSolventEff = (rhoSolvent * sfraction_se) + (rhoGas * sgf * (1.0 - sfraction_se)) + (rhoOil * sof * (1.0 - sfraction_se));
1190 }
1191
1192 // compute invB from densities.
1193 unsigned pvtRegionIdx = asImp_().pvtRegionIndex();
1194 Evaluation bGasEff = rhoGasEff;
1195 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1196 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) * fs.Rv());
1197 } else {
1198 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1199 }
1200 const Evaluation bSolventEff = rhoSolventEff / solventRefDensity();
1201
1202 // copy the unmodified invB factors
1203 const Evaluation bg = fs.invB(gasPhaseIdx);
1204 const Evaluation bs = solventInverseFormationVolumeFactor();
1205 const Evaluation bg_eff = pmisc * bGasEff + (1.0 - pmisc) * bg;
1206 const Evaluation bs_eff = pmisc * bSolventEff + (1.0 - pmisc) * bs;
1207
1208 // set the densities
1209 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1210 fs.setDensity(gasPhaseIdx,
1211 bg_eff
1212 *(FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)
1213 + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)*fs.Rv()));
1214 } else {
1215 fs.setDensity(gasPhaseIdx,
1216 bg_eff
1217 *FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
1218 }
1220
1221 // set the mobility
1222 Evaluation& mobg = asImp_().mobility_[gasPhaseIdx];
1223 muGasEff = bg_eff / (pmisc * bGasEff / muGasEff + (1.0 - pmisc) * bg / muGas);
1224 mobg *= muGas / muGasEff;
1225
1226 // Update viscosity of solvent
1227 solventViscosity_ = bs_eff / (pmisc * bSolventEff / muSolventEff + (1.0 - pmisc) * bs / muSolvent);
1228
1229 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
1230 Evaluation rhoOilEff = 0.0;
1231 if (std::abs(muOilPow.value() - muSolventPow.value()) < cutOff) {
1232 rhoOilEff = ((1.0 - tlMixParamRho) * rhoOil) + (tlMixParamRho * rhoMixSolventGasOil);
1233 }
1234 else {
1235 const Evaluation solventOilEffFraction = (muOilPow * (muOilEffPow - muSolventPow)) / (muOilEffPow * (muOilPow - muSolventPow));
1236 rhoOilEff = (rhoOil * solventOilEffFraction) + (rhoSolvent * (1.0 - solventOilEffFraction));
1237 }
1238 const Evaluation bOilEff = rhoOilEff / (FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) * fs.Rs());
1239 const Evaluation bo = fs.invB(oilPhaseIdx);
1240 const Evaluation bo_eff = pmisc * bOilEff + (1.0 - pmisc) * bo;
1241 fs.setDensity(oilPhaseIdx,
1242 bo_eff
1243 *(FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)
1244 + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)*fs.Rs()));
1245
1246 // keep the mu*b interpolation
1247 Evaluation& mobo = asImp_().mobility_[oilPhaseIdx];
1248 muOilEff = bo_eff / (pmisc * bOilEff / muOilEff + (1.0 - pmisc) * bo / muOil);
1249 mobo *= muOil / muOilEff;
1250 }
1251 }
1252
1253protected:
1254 Implementation& asImp_()
1255 { return *static_cast<Implementation*>(this); }
1256
1259 Evaluation rsSolw_;
1264
1266};
1267
1268template <class TypeTag>
1270{
1274
1275
1276public:
1277 void solventPreSatFuncUpdate_(const ElementContext&,
1278 unsigned,
1279 unsigned)
1280 { }
1281
1282 void solventPostSatFuncUpdate_(const ElementContext&,
1283 unsigned,
1284 unsigned)
1285 { }
1286
1287 void solventPvtUpdate_(const ElementContext&,
1288 unsigned,
1289 unsigned)
1290 { }
1291
1292 const Evaluation& solventSaturation() const
1293 { throw std::runtime_error("solventSaturation() called but solvents are disabled"); }
1294
1295 const Evaluation& rsSolw() const
1296 { throw std::runtime_error("rsSolw() called but solvents are disabled"); }
1297
1298 const Evaluation& solventDensity() const
1299 { throw std::runtime_error("solventDensity() called but solvents are disabled"); }
1300
1301 const Evaluation& solventViscosity() const
1302 { throw std::runtime_error("solventViscosity() called but solvents are disabled"); }
1303
1304 const Evaluation& solventMobility() const
1305 { throw std::runtime_error("solventMobility() called but solvents are disabled"); }
1306
1307 const Evaluation& solventInverseFormationVolumeFactor() const
1308 { throw std::runtime_error("solventInverseFormationVolumeFactor() called but solvents are disabled"); }
1309
1310 const Scalar& solventRefDensity() const
1311 { throw std::runtime_error("solventRefDensity() called but solvents are disabled"); }
1312};
1313
1321template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
1323{
1325
1333
1334 using Toolbox = MathToolbox<Evaluation>;
1335
1336 static constexpr unsigned gasPhaseIdx = FluidSystem::gasPhaseIdx;
1337 static constexpr int dimWorld = GridView::dimensionworld;
1338
1339 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
1340 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
1341
1342public:
1347 template <class Dummy = bool> // we need to make this method a template to avoid
1348 // compiler errors if it is not instantiated!
1349 void updateVolumeFluxPerm(const ElementContext& elemCtx,
1350 unsigned scvfIdx,
1351 unsigned timeIdx)
1352 {
1353 const auto& gradCalc = elemCtx.gradientCalculator();
1354 PressureCallback<TypeTag> pressureCallback(elemCtx);
1355
1356 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
1357 const auto& faceNormal = scvf.normal();
1358
1359 unsigned i = scvf.interiorIndex();
1360 unsigned j = scvf.exteriorIndex();
1361
1362 // calculate the "raw" pressure gradient
1363 DimEvalVector solventPGrad;
1364 pressureCallback.setPhaseIndex(gasPhaseIdx);
1365 gradCalc.calculateGradient(solventPGrad,
1366 elemCtx,
1367 scvfIdx,
1368 pressureCallback);
1369 Valgrind::CheckDefined(solventPGrad);
1370
1371 // correct the pressure gradients by the gravitational acceleration
1372 if (Parameters::Get<Parameters::EnableGravity>()) {
1373 // estimate the gravitational acceleration at a given SCV face
1374 // using the arithmetic mean
1375 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
1376 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
1377
1378 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
1379 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
1380
1381 const auto& posIn = elemCtx.pos(i, timeIdx);
1382 const auto& posEx = elemCtx.pos(j, timeIdx);
1383 const auto& posFace = scvf.integrationPos();
1384
1385 // the distance between the centers of the control volumes
1386 DimVector distVecIn(posIn);
1387 DimVector distVecEx(posEx);
1388 DimVector distVecTotal(posEx);
1389
1390 distVecIn -= posFace;
1391 distVecEx -= posFace;
1392 distVecTotal -= posIn;
1393 Scalar absDistTotalSquared = distVecTotal.two_norm2();
1394
1395 // calculate the hydrostatic pressure at the integration point of the face
1396 auto rhoIn = intQuantsIn.solventDensity();
1397 auto pStatIn = - rhoIn*(gIn*distVecIn);
1398
1399 // the quantities on the exterior side of the face do not influence the
1400 // result for the TPFA scheme, so they can be treated as scalar values.
1401 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1402 Scalar pStatEx = - rhoEx*(gEx*distVecEx);
1403
1404 // compute the hydrostatic gradient between the two control volumes (this
1405 // gradient exhibitis the same direction as the vector between the two
1406 // control volume centers and the length (pStaticExterior -
1407 // pStaticInterior)/distanceInteriorToExterior
1408 DimEvalVector f(distVecTotal);
1409 f *= (pStatEx - pStatIn)/absDistTotalSquared;
1410
1411 // calculate the final potential gradient
1412 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1413 solventPGrad[dimIdx] += f[dimIdx];
1414
1415 if (!isfinite(solventPGrad[dimIdx]))
1416 throw NumericalProblem("Non-finite potential gradient for solvent 'phase'");
1417 }
1418 }
1419
1420 // determine the upstream and downstream DOFs
1421 Evaluation solventPGradNormal = 0.0;
1422 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
1423 solventPGradNormal += solventPGrad[dimIdx]*faceNormal[dimIdx];
1424
1425 if (solventPGradNormal > 0) {
1426 solventUpstreamDofIdx_ = j;
1427 solventDownstreamDofIdx_ = i;
1428 }
1429 else {
1430 solventUpstreamDofIdx_ = i;
1431 solventDownstreamDofIdx_ = j;
1432 }
1433
1434 const auto& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1435
1436 // this is also slightly hacky because it assumes that the derivative of the
1437 // flux between two DOFs only depends on the primary variables in the
1438 // upstream direction. For non-TPFA flux approximation schemes, this is not
1439 // true...
1440 if (solventUpstreamDofIdx_ == i)
1441 solventVolumeFlux_ = solventPGradNormal*up.solventMobility();
1442 else
1443 solventVolumeFlux_ = solventPGradNormal*scalarValue(up.solventMobility());
1444 }
1445
1450 template <class Dummy = bool> // we need to make this method a template to avoid
1451 // compiler errors if it is not instantiated!
1452 void updateVolumeFluxTrans(const ElementContext& elemCtx,
1453 unsigned scvfIdx,
1454 unsigned timeIdx)
1455 {
1456 const ExtensiveQuantities& extQuants = asImp_();
1457
1458 unsigned interiorDofIdx = extQuants.interiorIndex();
1459 unsigned exteriorDofIdx = extQuants.exteriorIndex();
1460 assert(interiorDofIdx != exteriorDofIdx);
1461
1462 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
1463 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
1464
1465 unsigned I = elemCtx.globalSpaceIndex(interiorDofIdx, timeIdx);
1466 unsigned J = elemCtx.globalSpaceIndex(exteriorDofIdx, timeIdx);
1467
1468 Scalar thpres = elemCtx.problem().thresholdPressure(I, J);
1469 Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
1470 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
1471
1472 Scalar zIn = elemCtx.problem().dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
1473 Scalar zEx = elemCtx.problem().dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
1474 Scalar distZ = zIn - zEx;
1475
1476 const Evaluation& rhoIn = intQuantsIn.solventDensity();
1477 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1478 const Evaluation& rhoAvg = rhoIn*0.5 + rhoEx*0.5;
1479
1480 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(gasPhaseIdx);
1481 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(gasPhaseIdx));
1482 pressureExterior += distZ*g*rhoAvg;
1483
1484 Evaluation pressureDiffSolvent = pressureExterior - pressureInterior;
1485 if (std::abs(scalarValue(pressureDiffSolvent)) > thpres) {
1486 if (pressureDiffSolvent < 0.0)
1487 pressureDiffSolvent += thpres;
1488 else
1489 pressureDiffSolvent -= thpres;
1490 }
1491 else
1492 pressureDiffSolvent = 0.0;
1493
1494 if (pressureDiffSolvent > 0.0) {
1495 solventUpstreamDofIdx_ = exteriorDofIdx;
1496 solventDownstreamDofIdx_ = interiorDofIdx;
1497 }
1498 else if (pressureDiffSolvent < 0.0) {
1499 solventUpstreamDofIdx_ = interiorDofIdx;
1500 solventDownstreamDofIdx_ = exteriorDofIdx;
1501 }
1502 else {
1503 // pressure potential gradient is zero; force consistent upstream and
1504 // downstream indices over the intersection regardless of the side which it
1505 // is looked at.
1506 solventUpstreamDofIdx_ = std::min(interiorDofIdx, exteriorDofIdx);
1507 solventDownstreamDofIdx_ = std::max(interiorDofIdx, exteriorDofIdx);
1508 solventVolumeFlux_ = 0.0;
1509 return;
1510 }
1511
1512 Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
1513 const IntensiveQuantities& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1514 if (solventUpstreamDofIdx_ == interiorDofIdx)
1515 solventVolumeFlux_ =
1516 up.solventMobility()
1517 *(-trans/faceArea)
1518 *pressureDiffSolvent;
1519 else
1520 solventVolumeFlux_ =
1521 scalarValue(up.solventMobility())
1522 *(-trans/faceArea)
1523 *pressureDiffSolvent;
1524 }
1525
1526 unsigned solventUpstreamIndex() const
1527 { return solventUpstreamDofIdx_; }
1528
1529 unsigned solventDownstreamIndex() const
1530 { return solventDownstreamDofIdx_; }
1531
1532 const Evaluation& solventVolumeFlux() const
1533 { return solventVolumeFlux_; }
1534
1536 { solventVolumeFlux_ = solventVolumeFlux; }
1537
1538private:
1539 Implementation& asImp_()
1540 { return *static_cast<Implementation*>(this); }
1541
1542 Evaluation solventVolumeFlux_;
1543 unsigned solventUpstreamDofIdx_;
1544 unsigned solventDownstreamDofIdx_;
1545};
1546
1547template <class TypeTag>
1549{
1552
1553public:
1554 void updateVolumeFluxPerm(const ElementContext&,
1555 unsigned,
1556 unsigned)
1557 { }
1558
1559 void updateVolumeFluxTrans(const ElementContext&,
1560 unsigned,
1561 unsigned)
1562 { }
1563
1564 unsigned solventUpstreamIndex() const
1565 { throw std::runtime_error("solventUpstreamIndex() called but solvents are disabled"); }
1566
1567 unsigned solventDownstreamIndex() const
1568 { throw std::runtime_error("solventDownstreamIndex() called but solvents are disabled"); }
1569
1570 const Evaluation& solventVolumeFlux() const
1571 { throw std::runtime_error("solventVolumeFlux() called but solvents are disabled"); }
1572
1573 void setSolventVolumeFlux(const Evaluation& /* solventVolumeFlux */)
1574 { throw std::runtime_error("setSolventVolumeFlux() called but solvents are disabled"); }
1575};
1576
1577} // namespace Opm
1578
1579#endif
Declares the properties required by the black oil model.
Contains the parameters required to extend the black-oil model by solvents.
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1564
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1567
void setSolventVolumeFlux(const Evaluation &)
Definition: blackoilsolventmodules.hh:1573
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1570
void updateVolumeFluxPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1554
void updateVolumeFluxTrans(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1559
Provides the solvent specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilsolventmodules.hh:1323
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1526
void updateVolumeFluxPerm(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the volume flux of the polymer "phase" using the pressure potential gradient ...
Definition: blackoilsolventmodules.hh:1349
void updateVolumeFluxTrans(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Method which calculates the volume flux of the polymer "phase" using the gas pressure potential diffe...
Definition: blackoilsolventmodules.hh:1452
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1529
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1532
void setSolventVolumeFlux(const Evaluation &solventVolumeFlux)
Definition: blackoilsolventmodules.hh:1535
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:1310
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1298
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1307
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1295
void solventPreSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1277
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1301
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1304
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1292
void solventPvtUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1287
void solventPostSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1282
Provides the volumetric quantities required for the equations needed by the solvents extension of the...
Definition: blackoilsolventmodules.hh:780
void solventPreSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called before the saturation functions are doing their magic.
Definition: blackoilsolventmodules.hh:808
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1046
Evaluation solventViscosity_
Definition: blackoilsolventmodules.hh:1261
Evaluation solventMobility_
Definition: blackoilsolventmodules.hh:1262
Scalar solventRefDensity_
Definition: blackoilsolventmodules.hh:1265
Evaluation solventDensity_
Definition: blackoilsolventmodules.hh:1260
Evaluation rsSolw_
Definition: blackoilsolventmodules.hh:1259
Implementation & asImp_()
Definition: blackoilsolventmodules.hh:1254
Evaluation solventSaturation_
Definition: blackoilsolventmodules.hh:1258
void solventPvtUpdate_(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Update the intensive PVT properties needed to handle solvents from the primary variables.
Definition: blackoilsolventmodules.hh:970
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1040
Evaluation solventInvFormationVolumeFactor_
Definition: blackoilsolventmodules.hh:1263
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1037
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1043
void solventPostSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called after the saturation functions have been doing their magic.
Definition: blackoilsolventmodules.hh:838
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1049
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1034
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:1053
Evaluation hydrocarbonSaturation_
Definition: blackoilsolventmodules.hh:1257
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:69
static void setSolventPvt(const SolventPvt &value)
Specify the solvent PVT of a all PVT regions.
Definition: blackoilsolventmodules.hh:358
static const TabulatedFunction & sof2Krn(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:650
static const TabulatedFunction & ssfnKrs(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:642
static const TabulatedFunction & pmisc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:666
static const TabulatedFunction & msfnKrsg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:674
static const TabulatedFunction & sgcwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:698
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition: blackoilsolventmodules.hh:575
static const BrineH2Pvt & brineH2Pvt()
Definition: blackoilsolventmodules.hh:629
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:399
static std::string eqName(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:415
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:422
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:582
static bool isSolubleInWater()
Definition: blackoilsolventmodules.hh:748
static bool isCO2Sol()
Definition: blackoilsolventmodules.hh:753
static const TabulatedFunction & msfnKro(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:682
static void registerOutputModules(Model &model, Simulator &simulator)
Register all solvent specific VTK and ECL output modules.
Definition: blackoilsolventmodules.hh:377
static const TabulatedFunction & ssfnKrg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:634
static void setIsMiscible(const bool isMiscible)
Definition: blackoilsolventmodules.hh:362
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilsolventmodules.hh:368
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilsolventmodules.hh:563
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:608
static bool isH2Sol()
Definition: blackoilsolventmodules.hh:758
static bool isMiscible()
Definition: blackoilsolventmodules.hh:730
static const H2GasPvt & h2GasPvt()
Definition: blackoilsolventmodules.hh:619
static const Scalar & tlMixParamDensity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:722
static void updatePrimaryVars(PrimaryVariables &newPv, const PrimaryVariables &oldPv, const EqVector &delta)
Do a Newton-Raphson update the primary variables of the solvents.
Definition: blackoilsolventmodules.hh:551
static const TabulatedFunction & tlPMixTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:706
static const Scalar & tlMixParamViscosity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:714
static void assignPrimaryVars(PrimaryVariables &priVars, Scalar solventSaturation, Scalar solventRsw)
Assign the solvent specific primary variables to a PrimaryVariables object.
Definition: blackoilsolventmodules.hh:530
static const TabulatedFunction & sorwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:690
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:392
static bool eqApplies(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:407
static const Value solubilityLimit(unsigned pvtIdx, const Value &temperature, const Value &pressure, const Value &saltConcentration)
Definition: blackoilsolventmodules.hh:736
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:431
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:593
static const Co2GasPvt & co2GasPvt()
Definition: blackoilsolventmodules.hh:614
static const BrineCo2Pvt & brineCo2Pvt()
Definition: blackoilsolventmodules.hh:624
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:463
static const TabulatedFunction & misc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:658
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:384
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:84
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:108
VTK output module for the black oil model's solvent related quantities.
Definition: vtkblackoilsolventmodule.hh:63
Defines the common parameters for the porous medium multi-phase models.
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
This method contains all callback classes for quantities that are required by some extensive quantiti...
Struct holding the parameters for the BlackOilSolventModule class.
Definition: blackoilsolventparams.hh:43
::Opm::SolventPvt< Scalar > SolventPvt
Definition: blackoilsolventparams.hh:46
::Opm::H2GasPvt< Scalar > H2GasPvt
Definition: blackoilsolventparams.hh:52
::Opm::BrineCo2Pvt< Scalar > BrineCo2Pvt
Definition: blackoilsolventparams.hh:55
::Opm::Co2GasPvt< Scalar > Co2GasPvt
Definition: blackoilsolventparams.hh:49
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilsolventparams.hh:44
::Opm::BrineH2Pvt< Scalar > BrineH2Pvt
Definition: blackoilsolventparams.hh:58