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