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