ptflash/flashintensivequantities.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 OPM_FLASH_INTENSIVE_QUANTITIES_HH
29#define OPM_FLASH_INTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <opm/material/Constants.hpp>
35#include <opm/material/common/Valgrind.hpp>
36#include <opm/material/fluidstates/CompositionalFluidState.hpp>
37
40
42
45
46#include <array>
47#include <iostream>
48#include <string>
49
50namespace Opm {
51
58template <class TypeTag>
59class FlashIntensiveQuantities
60 : public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
61 , public DiffusionIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableDiffusion>() >
62 , public EnergyIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableEnergy>() >
63 , public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
64{
65 using ParentType = GetPropType<TypeTag, Properties::DiscIntensiveQuantities>;
66
67 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
68 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
69 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
70 using Indices = GetPropType<TypeTag, Properties::Indices>;
71 using FluxModule = GetPropType<TypeTag, Properties::FluxModule>;
72 using GridView = GetPropType<TypeTag, Properties::GridView>;
73 using ThreadManager = GetPropType<TypeTag, Properties::ThreadManager>;
74
75 // primary variable indices
76 enum { z0Idx = Indices::z0Idx };
77 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
78 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
79 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
80 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
81 enum { dimWorld = GridView::dimensionworld };
82 enum { pressure0Idx = Indices::pressure0Idx };
83 enum { water0Idx = Indices::water0Idx};
84
85 static constexpr bool waterEnabled = Indices::waterEnabled;
86
87 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
88 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
89 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
90 using FlashSolver = GetPropType<TypeTag, Properties::FlashSolver>;
91
92 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
93 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
94
95 using DiffusionIntensiveQuantities = ::Opm::DiffusionIntensiveQuantities<TypeTag, enableDiffusion>;
96 using EnergyIntensiveQuantities = ::Opm::EnergyIntensiveQuantities<TypeTag, enableEnergy>;
97 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
98
99public:
101 using FluidState = CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
102
104
106
108
112 void update(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
113 {
114 ParentType::update(elemCtx, dofIdx, timeIdx);
115 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
116
117 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
118 const auto& problem = elemCtx.problem();
119
120 const Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
121 const int flashVerbosity = Parameters::Get<Parameters::FlashVerbosity>();
122 const std::string flashTwoPhaseMethod = Parameters::Get<Parameters::FlashTwoPhaseMethod>();
123 // TODO: the formulation here is still to begin with XMF and YMF values to derive ZMF value
124 // TODO: we should check how we update ZMF in the newton update, since it is the primary variables.
125
126 // extract the total molar densities of the components
127 ComponentVector z(0.);
128 {
129 Evaluation lastZ = 1.0;
130 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
131 z[compIdx] = priVars.makeEvaluation(z0Idx + compIdx, timeIdx);
132 lastZ -= z[compIdx];
133 }
134 z[numComponents - 1] = lastZ;
135
136 Evaluation sumz = 0.0;
137 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
138 z[compIdx] = max(z[compIdx], 1e-8);
139 sumz += z[compIdx];
140 }
141 z /= sumz;
142 }
143
144 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
145 fluidState_.setMoleFraction(compIdx, z[compIdx]);
146 }
147
148 Evaluation p = priVars.makeEvaluation(pressure0Idx, timeIdx);
149 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
150 fluidState_.setPressure(phaseIdx, p);
151 }
152
153 // Get initial K and L from storage initially (if enabled)
154 const auto* hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
155 if (hint) {
156 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
157 const Evaluation& Ktmp = hint->fluidState().K(compIdx);
158 fluidState_.setKvalue(compIdx, Ktmp);
159 }
160 const Evaluation& Ltmp = hint->fluidState().L();
161 fluidState_.setLvalue(Ltmp);
162 }
163 else if (timeIdx == 0 && elemCtx.thermodynamicHint(dofIdx, 1)) {
164 // checking the storage cache
165 const auto& hint2 = elemCtx.thermodynamicHint(dofIdx, 1);
166 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
167 const Evaluation& Ktmp = hint2->fluidState().K(compIdx);
168 fluidState_.setKvalue(compIdx, Ktmp);
169 }
170 const Evaluation& Ltmp = hint2->fluidState().L();
171 fluidState_.setLvalue(Ltmp);
172 }
173 else {
174 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
175 const Evaluation Ktmp = fluidState_.wilsonK_(compIdx);
176 fluidState_.setKvalue(compIdx, Ktmp);
177 }
178 const Evaluation& Ltmp = -1.0;
179 fluidState_.setLvalue(Ltmp);
180 }
181
183 // Compute the phase compositions and densities
185 if (flashVerbosity >= 1) {
186 const int spatialIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
187 std::cout << " updating the intensive quantities for Cell " << spatialIdx << std::endl;
188 }
189 const auto& eos_type = problem.getEosType();
190 FlashSolver::solve(fluidState_, flashTwoPhaseMethod, flashTolerance, eos_type, flashVerbosity);
191
192 if (flashVerbosity >= 5) {
193 // printing of flash result after solve
194 const int spatialIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
195 std::cout << " \n After flash solve for cell " << spatialIdx << std::endl;
196 ComponentVector x, y;
197 for (unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
198 x[comp_idx] = fluidState_.moleFraction(FluidSystem::oilPhaseIdx, comp_idx);
199 y[comp_idx] = fluidState_.moleFraction(FluidSystem::gasPhaseIdx, comp_idx);
200 }
201 for (unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
202 std::cout << " x for component: " << comp_idx << " is:" << std::endl;
203 std::cout << x[comp_idx] << std::endl;
204
205 std::cout << " y for component: " << comp_idx << "is:" << std::endl;
206 std::cout << y[comp_idx] << std::endl;
207 }
208 const Evaluation& L = fluidState_.L();
209 std::cout << " L is:" << std::endl;
210 std::cout << L << std::endl;
211 }
212
213 // Update phases
214 typename FluidSystem::template ParameterCache<Evaluation> paramCache(eos_type);
215 paramCache.updatePhase(fluidState_, FluidSystem::oilPhaseIdx);
216
217 const Scalar R = Opm::Constants<Scalar>::R;
218 const Evaluation Z_L = (paramCache.molarVolume(FluidSystem::oilPhaseIdx) *
219 fluidState_.pressure(FluidSystem::oilPhaseIdx)) /
220 (R * fluidState_.temperature(FluidSystem::oilPhaseIdx));
221 paramCache.updatePhase(fluidState_, FluidSystem::gasPhaseIdx);
222 const Evaluation Z_V = (paramCache.molarVolume(FluidSystem::gasPhaseIdx) *
223 fluidState_.pressure(FluidSystem::gasPhaseIdx)) /
224 (R * fluidState_.temperature(FluidSystem::gasPhaseIdx));
225
226 // Update saturation
227 Evaluation Sw = 0.0;
228 if constexpr (waterEnabled) {
229 Sw = priVars.makeEvaluation(water0Idx, timeIdx);
230 }
231 const Evaluation L = fluidState_.L();
232 Evaluation So = max((1 - Sw) * (L * Z_L / ( L * Z_L + (1 - L) * Z_V)), 0.0);
233 Evaluation Sg = max(1 - So - Sw, 0.0);
234 const Scalar sumS = getValue(So) + getValue(Sg) + getValue(Sw);
235 So /= sumS;
236 Sg /= sumS;
237
238 fluidState_.setSaturation(FluidSystem::oilPhaseIdx, So);
239 fluidState_.setSaturation(FluidSystem::gasPhaseIdx, Sg);
240 if constexpr (waterEnabled) {
241 Sw /= sumS;
242 fluidState_.setSaturation(FluidSystem::waterPhaseIdx, Sw);
243 }
244
245 fluidState_.setCompressFactor(FluidSystem::oilPhaseIdx, Z_L);
246 fluidState_.setCompressFactor(FluidSystem::gasPhaseIdx, Z_V);
247
248 // Print saturation
249 if (flashVerbosity >= 5) {
250 std::cout << "So = " << So << std::endl;
251 std::cout << "Sg = " << Sg << std::endl;
252 std::cout << "Z_L = " << Z_L << std::endl;
253 std::cout << "Z_V = " << Z_V << std::endl;
254 }
255
257 // Compute rel. perm and viscosity and densities
259 const MaterialLawParams& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
260
261 // calculate relative permeability
262 MaterialLaw::relativePermeabilities(relativePermeability_,
263 materialParams, fluidState_);
264 Valgrind::CheckDefined(relativePermeability_);
265
266 // set the phase viscosity and density
267 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
268 if (phaseIdx == static_cast<unsigned int>(FluidSystem::oilPhaseIdx) ||
269 phaseIdx == static_cast<unsigned int>(FluidSystem::gasPhaseIdx))
270 {
271 paramCache.updatePhase(fluidState_, phaseIdx);
272 }
273
274 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
275
276 fluidState_.setViscosity(phaseIdx, mu);
277
278 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
279 Valgrind::CheckDefined(mobility_[phaseIdx]);
280
281 const Evaluation& rho = FluidSystem::density(fluidState_, paramCache, phaseIdx);
282 fluidState_.setDensity(phaseIdx, rho);
283 }
284
286 // Compute the remaining quantities
288
289 // porosity
290 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
291 Valgrind::CheckDefined(porosity_);
292
293 // intrinsic permeability
294 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
295
296 // update the quantities specific for the velocity model
297 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
298
299 // energy related quantities
300 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
301
302 // update the diffusion specific quantities of the intensive quantities
303 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
304 }
305
309 const FluidState& fluidState() const
310 { return fluidState_; }
311
315 const DimMatrix& intrinsicPermeability() const
316 { return intrinsicPerm_; }
317
321 const Evaluation& relativePermeability(unsigned phaseIdx) const
322 { return relativePermeability_[phaseIdx]; }
323
327 const Evaluation& mobility(unsigned phaseIdx) const
328 { return mobility_[phaseIdx]; }
329
333 const Evaluation& porosity() const
334 { return porosity_; }
335
336private:
337 DimMatrix intrinsicPerm_;
338 FluidState fluidState_;
339 Evaluation porosity_;
340 std::array<Evaluation,numPhases> relativePermeability_;
341 std::array<Evaluation,numPhases> mobility_;
342};
343
344} // namespace Opm
345
346#endif
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:145
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:540
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:60
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: ptflash/flashintensivequantities.hh:327
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: ptflash/flashintensivequantities.hh:315
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: ptflash/flashintensivequantities.hh:309
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: ptflash/flashintensivequantities.hh:333
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: ptflash/flashintensivequantities.hh:321
CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:93
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: ptflash/flashintensivequantities.hh:112
FlashIntensiveQuantities & operator=(const FlashIntensiveQuantities &other)=default
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declares the properties required by the compositional multi-phase model based on flash calculations.
Definition: blackoilboundaryratevector.hh:39
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
Declares the parameters for the compositional multi-phase model based on flash calculations.