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 <opm/common/Exceptions.hpp>
32
33#include <opm/material/fluidsystems/blackoilpvt/SolventPvt.hpp>
34
37
40
42
43#include <opm/material/common/Valgrind.hpp>
44
45#include <dune/common/fvector.hh>
46
47#include <string>
48
49namespace Opm {
50
56template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
58{
71 using Toolbox = MathToolbox<Evaluation>;
72 using SolventPvt = typename BlackOilSolventParams<Scalar>::SolventPvt;
73 using Co2GasPvt = typename BlackOilSolventParams<Scalar>::Co2GasPvt;
74 using H2GasPvt = typename BlackOilSolventParams<Scalar>::H2GasPvt;
75 using BrineCo2Pvt = typename BlackOilSolventParams<Scalar>::BrineCo2Pvt;
76 using BrineH2Pvt = typename BlackOilSolventParams<Scalar>::BrineH2Pvt;
77
78 using TabulatedFunction = typename BlackOilSolventParams<Scalar>::TabulatedFunction;
79
80 static constexpr unsigned solventSaturationIdx = Indices::solventSaturationIdx;
81 static constexpr unsigned contiSolventEqIdx = Indices::contiSolventEqIdx;
82 static constexpr unsigned enableSolvent = enableSolventV;
83 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
84 static constexpr unsigned numPhases = FluidSystem::numPhases;
85 static constexpr bool blackoilConserveSurfaceVolume = getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>();
86 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
87
88public:
91 {
92 params_ = params;
93 }
94
98 static void setSolventPvt(const SolventPvt& value)
99 { params_.solventPvt_ = value; }
100
101
102 static void setIsMiscible(const bool isMiscible)
103 { params_.isMiscible_ = isMiscible; }
104
108 static void registerParameters()
109 {
110 if constexpr (enableSolvent)
112 }
113
117 static void registerOutputModules(Model& model,
118 Simulator& simulator)
119 {
120 if constexpr (enableSolvent)
121 model.addOutputModule(new VtkBlackOilSolventModule<TypeTag>(simulator));
122 }
123
124 static bool primaryVarApplies(unsigned pvIdx)
125 {
126 if constexpr (enableSolvent)
127 return pvIdx == solventSaturationIdx;
128 else
129 return false;
130 }
131
132 static std::string primaryVarName([[maybe_unused]] unsigned pvIdx)
133 {
134 assert(primaryVarApplies(pvIdx));
135
136 return "saturation_solvent";
137 }
138
139 static Scalar primaryVarWeight([[maybe_unused]] unsigned pvIdx)
140 {
141 assert(primaryVarApplies(pvIdx));
142
143 // TODO: it may be beneficial to chose this differently.
144 return static_cast<Scalar>(1.0);
145 }
146
147 static bool eqApplies(unsigned eqIdx)
148 {
149 if constexpr (enableSolvent)
150 return eqIdx == contiSolventEqIdx;
151 else
152 return false;
153 }
154
155 static std::string eqName([[maybe_unused]] unsigned eqIdx)
156 {
157 assert(eqApplies(eqIdx));
158
159 return "conti^solvent";
160 }
161
162 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
163 {
164 assert(eqApplies(eqIdx));
165
166 // TODO: it may be beneficial to chose this differently.
167 return static_cast<Scalar>(1.0);
168 }
169
170 template <class LhsEval>
171 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
172 const IntensiveQuantities& intQuants)
173 {
174 if constexpr (enableSolvent) {
175 if constexpr (blackoilConserveSurfaceVolume) {
176 storage[contiSolventEqIdx] +=
177 Toolbox::template decay<LhsEval>(intQuants.porosity())
178 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
179 * Toolbox::template decay<LhsEval>(intQuants.solventInverseFormationVolumeFactor());
180 if (isSolubleInWater()) {
181 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
182 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
183 * Toolbox::template decay<LhsEval>(intQuants.fluidState().invB(waterPhaseIdx))
184 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
185 }
186 }
187 else {
188 storage[contiSolventEqIdx] +=
189 Toolbox::template decay<LhsEval>(intQuants.porosity())
190 * Toolbox::template decay<LhsEval>(intQuants.solventSaturation())
191 * Toolbox::template decay<LhsEval>(intQuants.solventDensity());
192 if (isSolubleInWater()) {
193 storage[contiSolventEqIdx] += Toolbox::template decay<LhsEval>(intQuants.porosity())
194 * Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(waterPhaseIdx))
195 * Toolbox::template decay<LhsEval>(intQuants.fluidState().density(waterPhaseIdx))
196 * Toolbox::template decay<LhsEval>(intQuants.rsSolw());
197 }
198
199 }
200 }
201 }
202
203 static void computeFlux([[maybe_unused]] RateVector& flux,
204 [[maybe_unused]] const ElementContext& elemCtx,
205 [[maybe_unused]] unsigned scvfIdx,
206 [[maybe_unused]] unsigned timeIdx)
207
208 {
209 if constexpr (enableSolvent) {
210 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
211
212 unsigned upIdx = extQuants.solventUpstreamIndex();
213 unsigned inIdx = extQuants.interiorIndex();
214 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
215
216 if constexpr (blackoilConserveSurfaceVolume) {
217 if (upIdx == inIdx)
218 flux[contiSolventEqIdx] =
219 extQuants.solventVolumeFlux()
220 *up.solventInverseFormationVolumeFactor();
221 else
222 flux[contiSolventEqIdx] =
223 extQuants.solventVolumeFlux()
224 *decay<Scalar>(up.solventInverseFormationVolumeFactor());
225
226
227 if (isSolubleInWater()) {
228 if (upIdx == inIdx)
229 flux[contiSolventEqIdx] +=
230 extQuants.volumeFlux(waterPhaseIdx)
231 * up.fluidState().invB(waterPhaseIdx)
232 * up.rsSolw();
233 else
234 flux[contiSolventEqIdx] +=
235 extQuants.volumeFlux(waterPhaseIdx)
236 *decay<Scalar>(up.fluidState().invB(waterPhaseIdx))
237 *decay<Scalar>(up.rsSolw());
238 }
239 }
240 else {
241 if (upIdx == inIdx)
242 flux[contiSolventEqIdx] =
243 extQuants.solventVolumeFlux()
244 *up.solventDensity();
245 else
246 flux[contiSolventEqIdx] =
247 extQuants.solventVolumeFlux()
248 *decay<Scalar>(up.solventDensity());
249
250
251 if (isSolubleInWater()) {
252 if (upIdx == inIdx)
253 flux[contiSolventEqIdx] +=
254 extQuants.volumeFlux(waterPhaseIdx)
255 * up.fluidState().density(waterPhaseIdx)
256 * up.rsSolw();
257 else
258 flux[contiSolventEqIdx] +=
259 extQuants.volumeFlux(waterPhaseIdx)
260 *decay<Scalar>(up.fluidState().density(waterPhaseIdx))
261 *decay<Scalar>(up.rsSolw());
262 }
263 }
264 }
265 }
266
270 static void assignPrimaryVars(PrimaryVariables& priVars,
271 Scalar solventSaturation,
272 Scalar solventRsw)
273 {
274 if constexpr (!enableSolvent) {
275 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Disabled);
276 return;
277 }
278 // Determine the meaning of the solvent primary variables
279 if (solventSaturation > 0 || !isSolubleInWater()) {
280 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Ss);
281 priVars[solventSaturationIdx] = solventSaturation;
282 } else {
283 priVars.setPrimaryVarsMeaningSolvent(PrimaryVariables::SolventMeaning::Rsolw);
284 priVars[solventSaturationIdx] = solventRsw;
285 }
286 }
287
291 static void updatePrimaryVars(PrimaryVariables& newPv,
292 const PrimaryVariables& oldPv,
293 const EqVector& delta)
294 {
295 if constexpr (enableSolvent)
296 // do a plain unchopped Newton update
297 newPv[solventSaturationIdx] = oldPv[solventSaturationIdx] - delta[solventSaturationIdx];
298 }
299
303 static Scalar computeUpdateError(const PrimaryVariables&,
304 const EqVector&)
305 {
306 // do not consider consider the cange of solvent primary variables for
307 // convergence
308 // TODO: maybe this should be changed
309 return static_cast<Scalar>(0.0);
310 }
311
315 static Scalar computeResidualError(const EqVector& resid)
316 {
317 // do not weight the residual of solvents when it comes to convergence
318 return std::abs(Toolbox::scalarValue(resid[contiSolventEqIdx]));
319 }
320
321 template <class DofEntity>
322 static void serializeEntity(const Model& model, std::ostream& outstream, const DofEntity& dof)
323 {
324 if constexpr (enableSolvent) {
325 unsigned dofIdx = model.dofMapper().index(dof);
326
327 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/0)[dofIdx];
328 outstream << priVars[solventSaturationIdx];
329 }
330 }
331
332 template <class DofEntity>
333 static void deserializeEntity(Model& model, std::istream& instream, const DofEntity& dof)
334 {
335 if constexpr (enableSolvent) {
336 unsigned dofIdx = model.dofMapper().index(dof);
337
338 PrimaryVariables& priVars0 = model.solution(/*timeIdx=*/0)[dofIdx];
339 PrimaryVariables& priVars1 = model.solution(/*timeIdx=*/1)[dofIdx];
340
341 instream >> priVars0[solventSaturationIdx];
342
343 // set the primary variables for the beginning of the current time step.
344 priVars1 = priVars0[solventSaturationIdx];
345 }
346 }
347
348 static const SolventPvt& solventPvt()
349 {
350 return params_.solventPvt_;
351 }
352
353
354 static const Co2GasPvt& co2GasPvt()
355 {
356 return params_.co2GasPvt_;
357 }
358
359 static const H2GasPvt& h2GasPvt()
360 {
361 return params_.h2GasPvt_;
362 }
363
364 static const BrineCo2Pvt& brineCo2Pvt()
365 {
366 return params_.brineCo2Pvt_;
367 }
368
369 static const BrineH2Pvt& brineH2Pvt()
370 {
371 return params_.brineH2Pvt_;
372 }
373
374 static const TabulatedFunction& ssfnKrg(const ElementContext& elemCtx,
375 unsigned scvIdx,
376 unsigned timeIdx)
377 {
378 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
379 return params_.ssfnKrg_[satnumRegionIdx];
380 }
381
382 static const TabulatedFunction& ssfnKrs(const ElementContext& elemCtx,
383 unsigned scvIdx,
384 unsigned timeIdx)
385 {
386 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
387 return params_.ssfnKrs_[satnumRegionIdx];
388 }
389
390 static const TabulatedFunction& sof2Krn(const ElementContext& elemCtx,
391 unsigned scvIdx,
392 unsigned timeIdx)
393 {
394 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
395 return params_.sof2Krn_[satnumRegionIdx];
396 }
397
398 static const TabulatedFunction& misc(const ElementContext& elemCtx,
399 unsigned scvIdx,
400 unsigned timeIdx)
401 {
402 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
403 return params_.misc_[miscnumRegionIdx];
404 }
405
406 static const TabulatedFunction& pmisc(const ElementContext& elemCtx,
407 unsigned scvIdx,
408 unsigned timeIdx)
409 {
410 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
411 return params_.pmisc_[miscnumRegionIdx];
412 }
413
414 static const TabulatedFunction& msfnKrsg(const ElementContext& elemCtx,
415 unsigned scvIdx,
416 unsigned timeIdx)
417 {
418 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
419 return params_.msfnKrsg_[satnumRegionIdx];
420 }
421
422 static const TabulatedFunction& msfnKro(const ElementContext& elemCtx,
423 unsigned scvIdx,
424 unsigned timeIdx)
425 {
426 unsigned satnumRegionIdx = elemCtx.problem().satnumRegionIndex(elemCtx, scvIdx, timeIdx);
427 return params_.msfnKro_[satnumRegionIdx];
428 }
429
430 static const TabulatedFunction& sorwmis(const ElementContext& elemCtx,
431 unsigned scvIdx,
432 unsigned timeIdx)
433 {
434 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
435 return params_.sorwmis_[miscnumRegionIdx];
436 }
437
438 static const TabulatedFunction& sgcwmis(const ElementContext& elemCtx,
439 unsigned scvIdx,
440 unsigned timeIdx)
441 {
442 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
443 return params_.sgcwmis_[miscnumRegionIdx];
444 }
445
446 static const TabulatedFunction& tlPMixTable(const ElementContext& elemCtx,
447 unsigned scvIdx,
448 unsigned timeIdx)
449 {
450 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
451 return params_.tlPMixTable_[miscnumRegionIdx];
452 }
453
454 static const Scalar& tlMixParamViscosity(const ElementContext& elemCtx,
455 unsigned scvIdx,
456 unsigned timeIdx)
457 {
458 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
459 return params_.tlMixParamViscosity_[miscnumRegionIdx];
460 }
461
462 static const Scalar& tlMixParamDensity(const ElementContext& elemCtx,
463 unsigned scvIdx,
464 unsigned timeIdx)
465 {
466 unsigned miscnumRegionIdx = elemCtx.problem().miscnumRegionIndex(elemCtx, scvIdx, timeIdx);
467 return params_.tlMixParamDensity_[miscnumRegionIdx];
468 }
469
470 static bool isMiscible()
471 {
472 return params_.isMiscible_;
473 }
474
475 template <class Value>
476 static const Value solubilityLimit(unsigned pvtIdx, const Value& temperature, const Value& pressure, const Value& saltConcentration)
477 {
478 if (!isSolubleInWater())
479 return 0.0;
480
481 assert(isCO2Sol() || isH2Sol());
482 if (isCO2Sol())
483 return brineCo2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
484 else
485 return brineH2Pvt().saturatedGasDissolutionFactor(pvtIdx, temperature, pressure, saltConcentration);
486 }
487
488 static bool isSolubleInWater()
489 {
490 return params_.rsSolw_active_;
491 }
492
493 static bool isCO2Sol()
494 {
495 return params_.co2sol_;
496 }
497
498 static bool isH2Sol()
499 {
500 return params_.h2sol_;
501 }
502
503private:
504 static BlackOilSolventParams<Scalar> params_; // the krg(Fs) column of the SSFN table
505};
506
507template <class TypeTag, bool enableSolventV>
508BlackOilSolventParams<typename BlackOilSolventModule<TypeTag, enableSolventV>::Scalar>
509BlackOilSolventModule<TypeTag, enableSolventV>::params_;
510
518template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
520{
522
530
532
533 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
534 static constexpr int solventSaturationIdx = Indices::solventSaturationIdx;
535 static constexpr int oilPhaseIdx = FluidSystem::oilPhaseIdx;
536 static constexpr int gasPhaseIdx = FluidSystem::gasPhaseIdx;
537 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
538 static constexpr double cutOff = 1e-12;
539
540
541public:
548 void solventPreSatFuncUpdate_(const ElementContext& elemCtx,
549 unsigned dofIdx,
550 unsigned timeIdx)
551 {
552 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
553
554 auto& fs = asImp_().fluidState_;
555 solventSaturation_ = 0.0;
556 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
557 solventSaturation_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
558 }
559
560 hydrocarbonSaturation_ = fs.saturation(gasPhaseIdx);
561
562 // apply a cut-off. Don't waste calculations if no solvent
563 if (solventSaturation().value() < cutOff)
564 return;
565
566 // make the saturation of the gas phase which is used by the saturation functions
567 // the sum of the solvent "saturation" and the saturation the hydrocarbon gas.
568 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_ + solventSaturation_);
569 }
570
578 void solventPostSatFuncUpdate_(const ElementContext& elemCtx,
579 unsigned dofIdx,
580 unsigned timeIdx)
581 {
582 // revert the gas "saturation" of the fluid state back to the saturation of the
583 // hydrocarbon gas.
584 auto& fs = asImp_().fluidState_;
585 fs.setSaturation(gasPhaseIdx, hydrocarbonSaturation_);
586
587
588 // update rsSolw. This needs to be done after the pressure is defined in the fluid state.
589 rsSolw_ = 0.0;
590 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
591 if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
592 rsSolw_ = SolventModule::solubilityLimit(asImp_().pvtRegionIndex(), fs.temperature(waterPhaseIdx), fs.pressure(waterPhaseIdx), fs.saltConcentration());
593 } else if (priVars.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Rsolw) {
594 rsSolw_ = priVars.makeEvaluation(solventSaturationIdx, timeIdx, elemCtx.linearizationType());
595 }
596
597 solventMobility_ = 0.0;
598
599 // apply a cut-off. Don't waste calculations if no solvent
600 if (solventSaturation().value() < cutOff)
601 return;
602
603 // Pressure effects on capillary pressure miscibility
605 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
606 const Evaluation pmisc = SolventModule::pmisc(elemCtx, dofIdx, timeIdx).eval(p, /*extrapolate=*/true);
607 const Evaluation& pgImisc = fs.pressure(gasPhaseIdx);
608
609 // compute capillary pressure for miscible fluid
610 const auto& problem = elemCtx.problem();
611 Evaluation pgMisc = 0.0;
612 std::array<Evaluation, numPhases> pC;
613 const auto& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
614 MaterialLaw::capillaryPressures(pC, materialParams, fs);
615
616 //oil is the reference phase for pressure
617 const auto linearizationType = elemCtx.linearizationType();
618 if (priVars.primaryVarsMeaningPressure() == PrimaryVariables::PressureMeaning::Pg)
619 pgMisc = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
620 else {
621 const Evaluation& po = priVars.makeEvaluation(Indices::pressureSwitchIdx, timeIdx, linearizationType);
622 pgMisc = po + (pC[gasPhaseIdx] - pC[oilPhaseIdx]);
623 }
624
625 fs.setPressure(gasPhaseIdx, pmisc * pgMisc + (1.0 - pmisc) * pgImisc);
626 }
627
628
629 Evaluation gasSolventSat = hydrocarbonSaturation_ + solventSaturation_;
630
631 if (gasSolventSat.value() < cutOff) // avoid division by zero
632 return;
633
634 Evaluation Fhydgas = hydrocarbonSaturation_/gasSolventSat;
635 Evaluation Fsolgas = solventSaturation_/gasSolventSat;
636
637 // account for miscibility of oil and solvent
638 if (SolventModule::isMiscible() && FluidSystem::phaseIsActive(oilPhaseIdx)) {
639 const auto& misc = SolventModule::misc(elemCtx, dofIdx, timeIdx);
640 const auto& pmisc = SolventModule::pmisc(elemCtx, dofIdx, timeIdx);
641 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
642 const Evaluation miscibility = misc.eval(Fsolgas, /*extrapolate=*/true) * pmisc.eval(p, /*extrapolate=*/true);
643
644 // TODO adjust endpoints of sn and ssg
645 unsigned cellIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
646 const auto& materialLawManager = elemCtx.problem().materialLawManager();
647 const auto& scaledDrainageInfo =
648 materialLawManager->oilWaterScaledEpsInfoDrainage(cellIdx);
649
650 const Scalar& sogcr = scaledDrainageInfo.Sogcr;
651 Evaluation sor = sogcr;
652 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
653 const Evaluation& sw = fs.saturation(waterPhaseIdx);
654 const auto& sorwmis = SolventModule::sorwmis(elemCtx, dofIdx, timeIdx);
655 sor = miscibility * sorwmis.eval(sw, /*extrapolate=*/true) + (1.0 - miscibility) * sogcr;
656 }
657 const Scalar& sgcr = scaledDrainageInfo.Sgcr;
658 Evaluation sgc = sgcr;
659 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
660 const Evaluation& sw = fs.saturation(waterPhaseIdx);
661 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, dofIdx, timeIdx);
662 sgc = miscibility * sgcwmis.eval(sw, /*extrapolate=*/true) + (1.0 - miscibility) * sgcr;
663 }
664
665 Evaluation oilGasSolventSat = gasSolventSat;
666 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
667 oilGasSolventSat += fs.saturation(oilPhaseIdx);
668 }
669 const Evaluation zero = 0.0;
670 const Evaluation oilGasSolventEffSat = std::max(oilGasSolventSat - sor - sgc, zero);
671
672 Evaluation F_totalGas = 0.0;
673 if (oilGasSolventEffSat.value() > cutOff) {
674 const Evaluation gasSolventEffSat = std::max(gasSolventSat - sgc, zero);
675 F_totalGas = gasSolventEffSat / oilGasSolventEffSat;
676 }
677 const auto& msfnKro = SolventModule::msfnKro(elemCtx, dofIdx, timeIdx);
678 const auto& msfnKrsg = SolventModule::msfnKrsg(elemCtx, dofIdx, timeIdx);
679 const auto& sof2Krn = SolventModule::sof2Krn(elemCtx, dofIdx, timeIdx);
680
681 const Evaluation mkrgt = msfnKrsg.eval(F_totalGas, /*extrapolate=*/true) * sof2Krn.eval(oilGasSolventSat, /*extrapolate=*/true);
682 const Evaluation mkro = msfnKro.eval(F_totalGas, /*extrapolate=*/true) * sof2Krn.eval(oilGasSolventSat, /*extrapolate=*/true);
683
684 Evaluation& kro = asImp_().mobility_[oilPhaseIdx];
685 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
686
687 // combine immiscible and miscible part of the relperm
688 krg *= (1.0 - miscibility);
689 krg += miscibility * mkrgt;
690 kro *= (1.0 - miscibility);
691 kro += miscibility * mkro;
692 }
693
694 // compute the mobility of the solvent "phase" and modify the gas phase
695 const auto& ssfnKrg = SolventModule::ssfnKrg(elemCtx, dofIdx, timeIdx);
696 const auto& ssfnKrs = SolventModule::ssfnKrs(elemCtx, dofIdx, timeIdx);
697
698 Evaluation& krg = asImp_().mobility_[gasPhaseIdx];
699 solventMobility_ = krg * ssfnKrs.eval(Fsolgas, /*extrapolate=*/true);
700 krg *= ssfnKrg.eval(Fhydgas, /*extrapolate=*/true);
701
702 }
703
710 void solventPvtUpdate_(const ElementContext& elemCtx,
711 unsigned scvIdx,
712 unsigned timeIdx)
713 {
714 const auto& iq = asImp_();
715 unsigned pvtRegionIdx = iq.pvtRegionIndex();
716 const Evaluation& T = iq.fluidState().temperature(gasPhaseIdx);
717 const Evaluation& p = iq.fluidState().pressure(gasPhaseIdx);
718
719 const Evaluation rv = 0.0;
720 const Evaluation rvw = 0.0;
723 const auto& co2gasPvt = SolventModule::co2GasPvt();
724 solventInvFormationVolumeFactor_ = co2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
725 solventRefDensity_ = co2gasPvt.gasReferenceDensity(pvtRegionIdx);
726 solventViscosity_ = co2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
727
728 const auto& brineCo2Pvt = SolventModule::brineCo2Pvt();
729 auto& fs = asImp_().fluidState_;
730 const auto& bw = brineCo2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
731
732 const auto denw = bw*brineCo2Pvt.waterReferenceDensity(pvtRegionIdx)
733 + rsSolw()*bw*brineCo2Pvt.gasReferenceDensity(pvtRegionIdx);
734 fs.setDensity(waterPhaseIdx, denw);
735 fs.setInvB(waterPhaseIdx, bw);
736 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
737 const auto& muWat = fs.viscosity(waterPhaseIdx);
738 const auto& muWatEff = brineCo2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
739 mobw *= muWat / muWatEff;
740 } else {
741 const auto& h2gasPvt = SolventModule::h2GasPvt();
742 solventInvFormationVolumeFactor_ = h2gasPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rv, rvw);
743 solventRefDensity_ = h2gasPvt.gasReferenceDensity(pvtRegionIdx);
744 solventViscosity_ = h2gasPvt.viscosity(pvtRegionIdx, T, p, rv, rvw);
745
746 const auto& brineH2Pvt = SolventModule::brineH2Pvt();
747 auto& fs = asImp_().fluidState_;
748 const auto& bw = brineH2Pvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p, rsSolw());
749
750 const auto denw = bw*brineH2Pvt.waterReferenceDensity(pvtRegionIdx)
751 + rsSolw()*bw*brineH2Pvt.gasReferenceDensity(pvtRegionIdx);
752 fs.setDensity(waterPhaseIdx, denw);
753 fs.setInvB(waterPhaseIdx, bw);
754 Evaluation& mobw = asImp_().mobility_[waterPhaseIdx];
755 const auto& muWat = fs.viscosity(waterPhaseIdx);
756 const auto& muWatEff = brineH2Pvt.viscosity(pvtRegionIdx, T, p, rsSolw());
757 mobw *= muWat / muWatEff;
758 }
759 } else {
760 const auto& solventPvt = SolventModule::solventPvt();
761 solventInvFormationVolumeFactor_ = solventPvt.inverseFormationVolumeFactor(pvtRegionIdx, T, p);
762 solventRefDensity_ = solventPvt.referenceDensity(pvtRegionIdx);
763 solventViscosity_ = solventPvt.viscosity(pvtRegionIdx, T, p);
764 }
765
767 effectiveProperties(elemCtx, scvIdx, timeIdx);
768
770
771
772 }
773
774 const Evaluation& solventSaturation() const
775 { return solventSaturation_; }
776
777 const Evaluation& rsSolw() const
778 { return rsSolw_; }
779
780 const Evaluation& solventDensity() const
781 { return solventDensity_; }
782
783 const Evaluation& solventViscosity() const
784 { return solventViscosity_; }
785
786 const Evaluation& solventMobility() const
787 { return solventMobility_; }
788
789 const Evaluation& solventInverseFormationVolumeFactor() const
791
792 // This could be stored pr pvtRegion instead
793 const Scalar& solventRefDensity() const
794 { return solventRefDensity_; }
795
796private:
797 // Computes the effective properties based on
798 // Todd-Longstaff mixing model.
799 void effectiveProperties(const ElementContext& elemCtx,
800 unsigned scvIdx,
801 unsigned timeIdx)
802 {
804 return;
805
806 // Don't waste calculations if no solvent
807 // Apply a cut-off for small and negative solvent saturations
808 if (solventSaturation() < cutOff)
809 return;
810
811 // We plan to update the fluidstate with the effective
812 // properties
813 auto& fs = asImp_().fluidState_;
814
815 // Compute effective saturations
816 Evaluation oilEffSat = 0.0;
817 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
818 oilEffSat = fs.saturation(oilPhaseIdx);
819 }
820 Evaluation gasEffSat = fs.saturation(gasPhaseIdx);
821 Evaluation solventEffSat = solventSaturation();
822 if (FluidSystem::phaseIsActive(waterPhaseIdx)) {
823 const auto& sorwmis = SolventModule::sorwmis(elemCtx, scvIdx, timeIdx);
824 const auto& sgcwmis = SolventModule::sgcwmis(elemCtx, scvIdx, timeIdx);
825 const Evaluation zero = 0.0;
826 const Evaluation& sw = fs.saturation(waterPhaseIdx);
827 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
828 oilEffSat = std::max(oilEffSat - sorwmis.eval(sw, /*extrapolate=*/true), zero);
829 }
830 gasEffSat = std::max(gasEffSat - sgcwmis.eval(sw, /*extrapolate=*/true), zero);
831 solventEffSat = std::max(solventEffSat - sgcwmis.eval(sw, /*extrapolate=*/true), zero);
832 }
833 const Evaluation oilGasSolventEffSat = oilEffSat + gasEffSat + solventEffSat;
834 const Evaluation oilSolventEffSat = oilEffSat + solventEffSat;
835 const Evaluation solventGasEffSat = solventEffSat + gasEffSat;
836
837 // Compute effective viscosities
838
839 // Mixing parameter for viscosity
840 // The pressureMixingParameter represent the miscibility of the solvent while the mixingParameterViscosity the effect of the porous media.
841 // The pressureMixingParameter is not implemented in ecl100.
842 const Evaluation& p = FluidSystem::phaseIsActive(oilPhaseIdx)? fs.pressure(oilPhaseIdx) : fs.pressure(gasPhaseIdx);
843 // account for pressure effects
844 const auto& pmiscTable = SolventModule::pmisc(elemCtx, scvIdx, timeIdx);
845 const Evaluation pmisc = pmiscTable.eval(p, /*extrapolate=*/true);
846 const auto& tlPMixTable = SolventModule::tlPMixTable(elemCtx, scvIdx, timeIdx);
847 const Evaluation tlMixParamMu = SolventModule::tlMixParamViscosity(elemCtx, scvIdx, timeIdx) * tlPMixTable.eval(p, /*extrapolate=*/true);
848
849 const Evaluation& muGas = fs.viscosity(gasPhaseIdx);
850 const Evaluation& muSolvent = solventViscosity_;
851
852 assert(muGas.value() > 0);
853 assert(muSolvent.value() > 0);
854 const Evaluation muGasPow = pow(muGas, 0.25);
855 const Evaluation muSolventPow = pow(muSolvent, 0.25);
856
857 Evaluation muMixSolventGas = muGas;
858 if (solventGasEffSat > cutOff)
859 muMixSolventGas *= muSolvent / pow(((gasEffSat / solventGasEffSat) * muSolventPow) + ((solventEffSat / solventGasEffSat) * muGasPow) , 4.0);
860
861 Evaluation muOil = 1.0;
862 Evaluation muOilPow = 1.0;
863 Evaluation muMixOilSolvent = 1.0;
864 Evaluation muOilEff = 1.0;
865 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
866 muOil = fs.viscosity(oilPhaseIdx);
867 assert(muOil.value() > 0);
868 muOilPow = pow(muOil, 0.25);
869 muMixOilSolvent = muOil;
870 if (oilSolventEffSat > cutOff)
871 muMixOilSolvent *= muSolvent / pow(((oilEffSat / oilSolventEffSat) * muSolventPow) + ((solventEffSat / oilSolventEffSat) * muOilPow) , 4.0);
872
873 muOilEff = pow(muOil,1.0 - tlMixParamMu) * pow(muMixOilSolvent, tlMixParamMu);
874 }
875 Evaluation muMixSolventGasOil = muOil;
876 if (oilGasSolventEffSat > cutOff)
877 muMixSolventGasOil *= muSolvent * muGas / pow(((oilEffSat / oilGasSolventEffSat) * muSolventPow * muGasPow)
878 + ((solventEffSat / oilGasSolventEffSat) * muOilPow * muGasPow) + ((gasEffSat / oilGasSolventEffSat) * muSolventPow * muOilPow), 4.0);
879
880 Evaluation muGasEff = pow(muGas,1.0 - tlMixParamMu) * pow(muMixSolventGas, tlMixParamMu);
881 Evaluation muSolventEff = pow(muSolvent,1.0 - tlMixParamMu) * pow(muMixSolventGasOil, tlMixParamMu);
882
883 // Compute effective densities
884 const Evaluation& rhoGas = fs.density(gasPhaseIdx);
885 Evaluation rhoOil = 0.0;
886 if (FluidSystem::phaseIsActive(oilPhaseIdx))
887 rhoOil = fs.density(oilPhaseIdx);
888
889 const Evaluation& rhoSolvent = solventDensity_;
890
891 // Mixing parameter for density
892 // The pressureMixingParameter represent the miscibility of the solvent while the mixingParameterDenisty the effect of the porous media.
893 // The pressureMixingParameter is not implemented in ecl100.
894 const Evaluation tlMixParamRho = SolventModule::tlMixParamDensity(elemCtx, scvIdx, timeIdx) * tlPMixTable.eval(p, /*extrapolate=*/true);
895
896 // compute effective viscosities for density calculations. These have to
897 // be recomputed as a different mixing parameter may be used.
898 const Evaluation muOilEffPow = pow(pow(muOil, 1.0 - tlMixParamRho) * pow(muMixOilSolvent, tlMixParamRho), 0.25);
899 const Evaluation muGasEffPow = pow(pow(muGas, 1.0 - tlMixParamRho) * pow(muMixSolventGas, tlMixParamRho), 0.25);
900 const Evaluation muSolventEffPow = pow(pow(muSolvent, 1.0 - tlMixParamRho) * pow(muMixSolventGasOil, tlMixParamRho), 0.25);
901
902 const Evaluation oilGasEffSaturation = oilEffSat + gasEffSat;
903 Evaluation sof = 0.0;
904 Evaluation sgf = 0.0;
905 if (oilGasEffSaturation.value() > cutOff) {
906 sof = oilEffSat / oilGasEffSaturation;
907 sgf = gasEffSat / oilGasEffSaturation;
908 }
909
910 const Evaluation muSolventOilGasPow = muSolventPow * ((sgf * muOilPow) + (sof * muGasPow));
911
912 Evaluation rhoMixSolventGasOil = 0.0;
913 if (oilGasSolventEffSat.value() > cutOff)
914 rhoMixSolventGasOil = (rhoOil * oilEffSat / oilGasSolventEffSat) + (rhoGas * gasEffSat / oilGasSolventEffSat) + (rhoSolvent * solventEffSat / oilGasSolventEffSat);
915
916 Evaluation rhoGasEff = 0.0;
917 if (std::abs(muSolventPow.value() - muGasPow.value()) < cutOff)
918 rhoGasEff = ((1.0 - tlMixParamRho) * rhoGas) + (tlMixParamRho * rhoMixSolventGasOil);
919 else {
920 const Evaluation solventGasEffFraction = (muGasPow * (muSolventPow - muGasEffPow)) / (muGasEffPow * (muSolventPow - muGasPow));
921 rhoGasEff = (rhoGas * solventGasEffFraction) + (rhoSolvent * (1.0 - solventGasEffFraction));
922 }
923
924 Evaluation rhoSolventEff = 0.0;
925 if (std::abs((muSolventOilGasPow.value() - (muOilPow.value() * muGasPow.value()))) < cutOff)
926 rhoSolventEff = ((1.0 - tlMixParamRho) * rhoSolvent) + (tlMixParamRho * rhoMixSolventGasOil);
927 else {
928 const Evaluation sfraction_se = (muSolventOilGasPow - (muOilPow * muGasPow * muSolventPow / muSolventEffPow)) / (muSolventOilGasPow - (muOilPow * muGasPow));
929 rhoSolventEff = (rhoSolvent * sfraction_se) + (rhoGas * sgf * (1.0 - sfraction_se)) + (rhoOil * sof * (1.0 - sfraction_se));
930 }
931
932 // compute invB from densities.
933 unsigned pvtRegionIdx = asImp_().pvtRegionIndex();
934 Evaluation bGasEff = rhoGasEff;
935 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
936 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) * fs.Rv());
937 } else {
938 bGasEff /= (FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
939 }
940 const Evaluation bSolventEff = rhoSolventEff / solventRefDensity();
941
942 // copy the unmodified invB factors
943 const Evaluation bg = fs.invB(gasPhaseIdx);
944 const Evaluation bs = solventInverseFormationVolumeFactor();
945 const Evaluation bg_eff = pmisc * bGasEff + (1.0 - pmisc) * bg;
946 const Evaluation bs_eff = pmisc * bSolventEff + (1.0 - pmisc) * bs;
947
948 // set the densities
949 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
950 fs.setDensity(gasPhaseIdx,
951 bg_eff
952 *(FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)
953 + FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)*fs.Rv()));
954 } else {
955 fs.setDensity(gasPhaseIdx,
956 bg_eff
957 *FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx));
958 }
960
961 // set the mobility
962 Evaluation& mobg = asImp_().mobility_[gasPhaseIdx];
963 muGasEff = bg_eff / (pmisc * bGasEff / muGasEff + (1.0 - pmisc) * bg / muGas);
964 mobg *= muGas / muGasEff;
965
966 // Update viscosity of solvent
967 solventViscosity_ = bs_eff / (pmisc * bSolventEff / muSolventEff + (1.0 - pmisc) * bs / muSolvent);
968
969 if (FluidSystem::phaseIsActive(oilPhaseIdx)) {
970 Evaluation rhoOilEff = 0.0;
971 if (std::abs(muOilPow.value() - muSolventPow.value()) < cutOff) {
972 rhoOilEff = ((1.0 - tlMixParamRho) * rhoOil) + (tlMixParamRho * rhoMixSolventGasOil);
973 }
974 else {
975 const Evaluation solventOilEffFraction = (muOilPow * (muOilEffPow - muSolventPow)) / (muOilEffPow * (muOilPow - muSolventPow));
976 rhoOilEff = (rhoOil * solventOilEffFraction) + (rhoSolvent * (1.0 - solventOilEffFraction));
977 }
978 const Evaluation bOilEff = rhoOilEff / (FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx) + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx) * fs.Rs());
979 const Evaluation bo = fs.invB(oilPhaseIdx);
980 const Evaluation bo_eff = pmisc * bOilEff + (1.0 - pmisc) * bo;
981 fs.setDensity(oilPhaseIdx,
982 bo_eff
983 *(FluidSystem::referenceDensity(oilPhaseIdx, pvtRegionIdx)
984 + FluidSystem::referenceDensity(gasPhaseIdx, pvtRegionIdx)*fs.Rs()));
985
986 // keep the mu*b interpolation
987 Evaluation& mobo = asImp_().mobility_[oilPhaseIdx];
988 muOilEff = bo_eff / (pmisc * bOilEff / muOilEff + (1.0 - pmisc) * bo / muOil);
989 mobo *= muOil / muOilEff;
990 }
991 }
992
993protected:
994 Implementation& asImp_()
995 { return *static_cast<Implementation*>(this); }
996
999 Evaluation rsSolw_;
1004
1006};
1007
1008template <class TypeTag>
1010{
1014
1015
1016public:
1017 void solventPreSatFuncUpdate_(const ElementContext&,
1018 unsigned,
1019 unsigned)
1020 { }
1021
1022 void solventPostSatFuncUpdate_(const ElementContext&,
1023 unsigned,
1024 unsigned)
1025 { }
1026
1027 void solventPvtUpdate_(const ElementContext&,
1028 unsigned,
1029 unsigned)
1030 { }
1031
1032 const Evaluation& solventSaturation() const
1033 { throw std::runtime_error("solventSaturation() called but solvents are disabled"); }
1034
1035 const Evaluation& rsSolw() const
1036 { throw std::runtime_error("rsSolw() called but solvents are disabled"); }
1037
1038 const Evaluation& solventDensity() const
1039 { throw std::runtime_error("solventDensity() called but solvents are disabled"); }
1040
1041 const Evaluation& solventViscosity() const
1042 { throw std::runtime_error("solventViscosity() called but solvents are disabled"); }
1043
1044 const Evaluation& solventMobility() const
1045 { throw std::runtime_error("solventMobility() called but solvents are disabled"); }
1046
1047 const Evaluation& solventInverseFormationVolumeFactor() const
1048 { throw std::runtime_error("solventInverseFormationVolumeFactor() called but solvents are disabled"); }
1049
1050 const Scalar& solventRefDensity() const
1051 { throw std::runtime_error("solventRefDensity() called but solvents are disabled"); }
1052};
1053
1061template <class TypeTag, bool enableSolventV = getPropValue<TypeTag, Properties::EnableSolvent>()>
1063{
1065
1073
1074 using Toolbox = MathToolbox<Evaluation>;
1075
1076 static constexpr unsigned gasPhaseIdx = FluidSystem::gasPhaseIdx;
1077 static constexpr int dimWorld = GridView::dimensionworld;
1078
1079 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
1080 using DimEvalVector = Dune::FieldVector<Evaluation, dimWorld>;
1081
1082public:
1087 template <class Dummy = bool> // we need to make this method a template to avoid
1088 // compiler errors if it is not instantiated!
1089 void updateVolumeFluxPerm(const ElementContext& elemCtx,
1090 unsigned scvfIdx,
1091 unsigned timeIdx)
1092 {
1093 const auto& gradCalc = elemCtx.gradientCalculator();
1094 PressureCallback<TypeTag> pressureCallback(elemCtx);
1095
1096 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
1097 const auto& faceNormal = scvf.normal();
1098
1099 unsigned i = scvf.interiorIndex();
1100 unsigned j = scvf.exteriorIndex();
1101
1102 // calculate the "raw" pressure gradient
1103 DimEvalVector solventPGrad;
1104 pressureCallback.setPhaseIndex(gasPhaseIdx);
1105 gradCalc.calculateGradient(solventPGrad,
1106 elemCtx,
1107 scvfIdx,
1108 pressureCallback);
1109 Valgrind::CheckDefined(solventPGrad);
1110
1111 // correct the pressure gradients by the gravitational acceleration
1112 if (Parameters::Get<Parameters::EnableGravity>()) {
1113 // estimate the gravitational acceleration at a given SCV face
1114 // using the arithmetic mean
1115 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
1116 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
1117
1118 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
1119 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
1120
1121 const auto& posIn = elemCtx.pos(i, timeIdx);
1122 const auto& posEx = elemCtx.pos(j, timeIdx);
1123 const auto& posFace = scvf.integrationPos();
1124
1125 // the distance between the centers of the control volumes
1126 DimVector distVecIn(posIn);
1127 DimVector distVecEx(posEx);
1128 DimVector distVecTotal(posEx);
1129
1130 distVecIn -= posFace;
1131 distVecEx -= posFace;
1132 distVecTotal -= posIn;
1133 Scalar absDistTotalSquared = distVecTotal.two_norm2();
1134
1135 // calculate the hydrostatic pressure at the integration point of the face
1136 auto rhoIn = intQuantsIn.solventDensity();
1137 auto pStatIn = - rhoIn*(gIn*distVecIn);
1138
1139 // the quantities on the exterior side of the face do not influence the
1140 // result for the TPFA scheme, so they can be treated as scalar values.
1141 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1142 Scalar pStatEx = - rhoEx*(gEx*distVecEx);
1143
1144 // compute the hydrostatic gradient between the two control volumes (this
1145 // gradient exhibitis the same direction as the vector between the two
1146 // control volume centers and the length (pStaticExterior -
1147 // pStaticInterior)/distanceInteriorToExterior
1148 DimEvalVector f(distVecTotal);
1149 f *= (pStatEx - pStatIn)/absDistTotalSquared;
1150
1151 // calculate the final potential gradient
1152 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1153 solventPGrad[dimIdx] += f[dimIdx];
1154
1155 if (!isfinite(solventPGrad[dimIdx]))
1156 throw NumericalProblem("Non-finite potential gradient for solvent 'phase'");
1157 }
1158 }
1159
1160 // determine the upstream and downstream DOFs
1161 Evaluation solventPGradNormal = 0.0;
1162 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
1163 solventPGradNormal += solventPGrad[dimIdx]*faceNormal[dimIdx];
1164
1165 if (solventPGradNormal > 0) {
1166 solventUpstreamDofIdx_ = j;
1167 solventDownstreamDofIdx_ = i;
1168 }
1169 else {
1170 solventUpstreamDofIdx_ = i;
1171 solventDownstreamDofIdx_ = j;
1172 }
1173
1174 const auto& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1175
1176 // this is also slightly hacky because it assumes that the derivative of the
1177 // flux between two DOFs only depends on the primary variables in the
1178 // upstream direction. For non-TPFA flux approximation schemes, this is not
1179 // true...
1180 if (solventUpstreamDofIdx_ == i)
1181 solventVolumeFlux_ = solventPGradNormal*up.solventMobility();
1182 else
1183 solventVolumeFlux_ = solventPGradNormal*scalarValue(up.solventMobility());
1184 }
1185
1190 template <class Dummy = bool> // we need to make this method a template to avoid
1191 // compiler errors if it is not instantiated!
1192 void updateVolumeFluxTrans(const ElementContext& elemCtx,
1193 unsigned scvfIdx,
1194 unsigned timeIdx)
1195 {
1196 const ExtensiveQuantities& extQuants = asImp_();
1197
1198 unsigned interiorDofIdx = extQuants.interiorIndex();
1199 unsigned exteriorDofIdx = extQuants.exteriorIndex();
1200 assert(interiorDofIdx != exteriorDofIdx);
1201
1202 const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx, timeIdx);
1203 const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx, timeIdx);
1204
1205 unsigned I = elemCtx.globalSpaceIndex(interiorDofIdx, timeIdx);
1206 unsigned J = elemCtx.globalSpaceIndex(exteriorDofIdx, timeIdx);
1207
1208 Scalar thpres = elemCtx.problem().thresholdPressure(I, J);
1209 Scalar trans = elemCtx.problem().transmissibility(elemCtx, interiorDofIdx, exteriorDofIdx);
1210 Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
1211
1212 Scalar zIn = elemCtx.problem().dofCenterDepth(elemCtx, interiorDofIdx, timeIdx);
1213 Scalar zEx = elemCtx.problem().dofCenterDepth(elemCtx, exteriorDofIdx, timeIdx);
1214 Scalar distZ = zIn - zEx;
1215
1216 const Evaluation& rhoIn = intQuantsIn.solventDensity();
1217 Scalar rhoEx = Toolbox::value(intQuantsEx.solventDensity());
1218 const Evaluation& rhoAvg = rhoIn*0.5 + rhoEx*0.5;
1219
1220 const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(gasPhaseIdx);
1221 Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(gasPhaseIdx));
1222 pressureExterior += distZ*g*rhoAvg;
1223
1224 Evaluation pressureDiffSolvent = pressureExterior - pressureInterior;
1225 if (std::abs(scalarValue(pressureDiffSolvent)) > thpres) {
1226 if (pressureDiffSolvent < 0.0)
1227 pressureDiffSolvent += thpres;
1228 else
1229 pressureDiffSolvent -= thpres;
1230 }
1231 else
1232 pressureDiffSolvent = 0.0;
1233
1234 if (pressureDiffSolvent > 0.0) {
1235 solventUpstreamDofIdx_ = exteriorDofIdx;
1236 solventDownstreamDofIdx_ = interiorDofIdx;
1237 }
1238 else if (pressureDiffSolvent < 0.0) {
1239 solventUpstreamDofIdx_ = interiorDofIdx;
1240 solventDownstreamDofIdx_ = exteriorDofIdx;
1241 }
1242 else {
1243 // pressure potential gradient is zero; force consistent upstream and
1244 // downstream indices over the intersection regardless of the side which it
1245 // is looked at.
1246 solventUpstreamDofIdx_ = std::min(interiorDofIdx, exteriorDofIdx);
1247 solventDownstreamDofIdx_ = std::max(interiorDofIdx, exteriorDofIdx);
1248 solventVolumeFlux_ = 0.0;
1249 return;
1250 }
1251
1252 Scalar faceArea = elemCtx.stencil(timeIdx).interiorFace(scvfIdx).area();
1253 const IntensiveQuantities& up = elemCtx.intensiveQuantities(solventUpstreamDofIdx_, timeIdx);
1254 if (solventUpstreamDofIdx_ == interiorDofIdx)
1255 solventVolumeFlux_ =
1256 up.solventMobility()
1257 *(-trans/faceArea)
1258 *pressureDiffSolvent;
1259 else
1260 solventVolumeFlux_ =
1261 scalarValue(up.solventMobility())
1262 *(-trans/faceArea)
1263 *pressureDiffSolvent;
1264 }
1265
1266 unsigned solventUpstreamIndex() const
1267 { return solventUpstreamDofIdx_; }
1268
1269 unsigned solventDownstreamIndex() const
1270 { return solventDownstreamDofIdx_; }
1271
1272 const Evaluation& solventVolumeFlux() const
1273 { return solventVolumeFlux_; }
1274
1276 { solventVolumeFlux_ = solventVolumeFlux; }
1277
1278private:
1279 Implementation& asImp_()
1280 { return *static_cast<Implementation*>(this); }
1281
1282 Evaluation solventVolumeFlux_;
1283 unsigned solventUpstreamDofIdx_;
1284 unsigned solventDownstreamDofIdx_;
1285};
1286
1287template <class TypeTag>
1289{
1292
1293public:
1294 void updateVolumeFluxPerm(const ElementContext&,
1295 unsigned,
1296 unsigned)
1297 { }
1298
1299 void updateVolumeFluxTrans(const ElementContext&,
1300 unsigned,
1301 unsigned)
1302 { }
1303
1304 unsigned solventUpstreamIndex() const
1305 { throw std::runtime_error("solventUpstreamIndex() called but solvents are disabled"); }
1306
1307 unsigned solventDownstreamIndex() const
1308 { throw std::runtime_error("solventDownstreamIndex() called but solvents are disabled"); }
1309
1310 const Evaluation& solventVolumeFlux() const
1311 { throw std::runtime_error("solventVolumeFlux() called but solvents are disabled"); }
1312
1313 void setSolventVolumeFlux(const Evaluation& /* solventVolumeFlux */)
1314 { throw std::runtime_error("setSolventVolumeFlux() called but solvents are disabled"); }
1315};
1316
1317} // namespace Opm
1318
1319#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:1304
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1307
void setSolventVolumeFlux(const Evaluation &)
Definition: blackoilsolventmodules.hh:1313
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1310
void updateVolumeFluxPerm(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1294
void updateVolumeFluxTrans(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1299
Provides the solvent specific extensive quantities to the generic black-oil module's extensive quanti...
Definition: blackoilsolventmodules.hh:1063
unsigned solventUpstreamIndex() const
Definition: blackoilsolventmodules.hh:1266
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:1089
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:1192
unsigned solventDownstreamIndex() const
Definition: blackoilsolventmodules.hh:1269
const Evaluation & solventVolumeFlux() const
Definition: blackoilsolventmodules.hh:1272
void setSolventVolumeFlux(const Evaluation &solventVolumeFlux)
Definition: blackoilsolventmodules.hh:1275
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:1050
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:1038
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:1047
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:1035
void solventPreSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1017
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:1041
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:1044
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:1032
void solventPvtUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1027
void solventPostSatFuncUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilsolventmodules.hh:1022
Provides the volumetric quantities required for the equations needed by the solvents extension of the...
Definition: blackoilsolventmodules.hh:520
void solventPreSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called before the saturation functions are doing their magic.
Definition: blackoilsolventmodules.hh:548
const Evaluation & solventMobility() const
Definition: blackoilsolventmodules.hh:786
Evaluation solventViscosity_
Definition: blackoilsolventmodules.hh:1001
Evaluation solventMobility_
Definition: blackoilsolventmodules.hh:1002
Scalar solventRefDensity_
Definition: blackoilsolventmodules.hh:1005
Evaluation solventDensity_
Definition: blackoilsolventmodules.hh:1000
Evaluation rsSolw_
Definition: blackoilsolventmodules.hh:999
Implementation & asImp_()
Definition: blackoilsolventmodules.hh:994
Evaluation solventSaturation_
Definition: blackoilsolventmodules.hh:998
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:710
const Evaluation & solventDensity() const
Definition: blackoilsolventmodules.hh:780
Evaluation solventInvFormationVolumeFactor_
Definition: blackoilsolventmodules.hh:1003
const Evaluation & rsSolw() const
Definition: blackoilsolventmodules.hh:777
const Evaluation & solventViscosity() const
Definition: blackoilsolventmodules.hh:783
void solventPostSatFuncUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Called after the saturation functions have been doing their magic.
Definition: blackoilsolventmodules.hh:578
const Evaluation & solventInverseFormationVolumeFactor() const
Definition: blackoilsolventmodules.hh:789
const Evaluation & solventSaturation() const
Definition: blackoilsolventmodules.hh:774
const Scalar & solventRefDensity() const
Definition: blackoilsolventmodules.hh:793
Evaluation hydrocarbonSaturation_
Definition: blackoilsolventmodules.hh:997
Contains the high level supplements required to extend the black oil model by solvents.
Definition: blackoilsolventmodules.hh:58
static void setSolventPvt(const SolventPvt &value)
Specify the solvent PVT of a all PVT regions.
Definition: blackoilsolventmodules.hh:98
static const TabulatedFunction & sof2Krn(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:390
static const TabulatedFunction & ssfnKrs(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:382
static const TabulatedFunction & pmisc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:406
static const TabulatedFunction & msfnKrsg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:414
static const TabulatedFunction & sgcwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:438
static Scalar computeResidualError(const EqVector &resid)
Return how much a residual is considered an error.
Definition: blackoilsolventmodules.hh:315
static const BrineH2Pvt & brineH2Pvt()
Definition: blackoilsolventmodules.hh:369
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:139
static void setParams(BlackOilSolventParams< Scalar > &&params)
Set parameters.
Definition: blackoilsolventmodules.hh:90
static std::string eqName(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:155
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:162
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:322
static bool isSolubleInWater()
Definition: blackoilsolventmodules.hh:488
static bool isCO2Sol()
Definition: blackoilsolventmodules.hh:493
static const TabulatedFunction & msfnKro(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:422
static void registerOutputModules(Model &model, Simulator &simulator)
Register all solvent specific VTK and ECL output modules.
Definition: blackoilsolventmodules.hh:117
static const TabulatedFunction & ssfnKrg(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:374
static void setIsMiscible(const bool isMiscible)
Definition: blackoilsolventmodules.hh:102
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilsolventmodules.hh:108
static Scalar computeUpdateError(const PrimaryVariables &, const EqVector &)
Return how much a Newton-Raphson update is considered an error.
Definition: blackoilsolventmodules.hh:303
static const SolventPvt & solventPvt()
Definition: blackoilsolventmodules.hh:348
static bool isH2Sol()
Definition: blackoilsolventmodules.hh:498
static bool isMiscible()
Definition: blackoilsolventmodules.hh:470
static const H2GasPvt & h2GasPvt()
Definition: blackoilsolventmodules.hh:359
static const Scalar & tlMixParamDensity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:462
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:291
static const TabulatedFunction & tlPMixTable(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:446
static const Scalar & tlMixParamViscosity(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:454
static void assignPrimaryVars(PrimaryVariables &priVars, Scalar solventSaturation, Scalar solventRsw)
Assign the solvent specific primary variables to a PrimaryVariables object.
Definition: blackoilsolventmodules.hh:270
static const TabulatedFunction & sorwmis(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:430
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:132
static bool eqApplies(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:147
static const Value solubilityLimit(unsigned pvtIdx, const Value &temperature, const Value &pressure, const Value &saltConcentration)
Definition: blackoilsolventmodules.hh:476
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilsolventmodules.hh:171
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:333
static const Co2GasPvt & co2GasPvt()
Definition: blackoilsolventmodules.hh:354
static const BrineCo2Pvt & brineCo2Pvt()
Definition: blackoilsolventmodules.hh:364
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:203
static const TabulatedFunction & misc(const ElementContext &elemCtx, unsigned scvIdx, unsigned timeIdx)
Definition: blackoilsolventmodules.hh:398
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:124
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.hpp:54
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.hpp:49
::Opm::SolventPvt< Scalar > SolventPvt
Definition: blackoilsolventparams.hpp:54
::Opm::H2GasPvt< Scalar > H2GasPvt
Definition: blackoilsolventparams.hpp:53
::Opm::BrineCo2Pvt< Scalar > BrineCo2Pvt
Definition: blackoilsolventparams.hpp:50
::Opm::Co2GasPvt< Scalar > Co2GasPvt
Definition: blackoilsolventparams.hpp:52
Tabulated1DFunction< Scalar > TabulatedFunction
Definition: blackoilsolventparams.hpp:55
::Opm::BrineH2Pvt< Scalar > BrineH2Pvt
Definition: blackoilsolventparams.hpp:51