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
46namespace Opm {
47
54template <class TypeTag>
55class FlashIntensiveQuantities
56 : public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
57 , public DiffusionIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableDiffusion>() >
58 , public EnergyIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableEnergy>() >
59 , public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
60{
61 using ParentType = GetPropType<TypeTag, Properties::DiscIntensiveQuantities>;
62
63 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
64 using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
65 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
66 using Indices = GetPropType<TypeTag, Properties::Indices>;
67 using FluxModule = GetPropType<TypeTag, Properties::FluxModule>;
68 using GridView = GetPropType<TypeTag, Properties::GridView>;
69 using ThreadManager = GetPropType<TypeTag, Properties::ThreadManager>;
70
71 // primary variable indices
72 enum { z0Idx = Indices::z0Idx };
73 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
74 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
75 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
76 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
77 enum { dimWorld = GridView::dimensionworld };
78 enum { pressure0Idx = Indices::pressure0Idx };
79
80 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
81 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
82 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
83 using FlashSolver = GetPropType<TypeTag, Properties::FlashSolver>;
84
85 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
86 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
87
88 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
89 using DiffusionIntensiveQuantities = Opm::DiffusionIntensiveQuantities<TypeTag, enableDiffusion>;
90 using EnergyIntensiveQuantities = Opm::EnergyIntensiveQuantities<TypeTag, enableEnergy>;
91
92public:
94 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
95
97
99
101
105 void update(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
106 {
107 ParentType::update(elemCtx, dofIdx, timeIdx);
108 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
109
110 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
111 const auto& problem = elemCtx.problem();
112
113 const Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
114 const int flashVerbosity = Parameters::Get<Parameters::FlashVerbosity>();
115 const std::string flashTwoPhaseMethod = Parameters::Get<Parameters::FlashTwoPhaseMethod>();
116
117 // extract the total molar densities of the components
118 ComponentVector z(0.);
119 {
120 Evaluation lastZ = 1.0;
121 for (unsigned compIdx = 0; compIdx < numComponents - 1; ++compIdx) {
122 z[compIdx] = priVars.makeEvaluation(z0Idx + compIdx, timeIdx);
123 lastZ -= z[compIdx];
124 }
125 z[numComponents - 1] = lastZ;
126
127 Evaluation sumz = 0.0;
128 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
129 z[compIdx] = Opm::max(z[compIdx], 1e-8);
130 sumz +=z[compIdx];
131 }
132 z /= sumz;
133 }
134
135 Evaluation p = priVars.makeEvaluation(pressure0Idx, timeIdx);
136 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
137 fluidState_.setPressure(phaseIdx, p);
138
139 // Get initial K and L from storage initially (if enabled)
140 const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
141 if (hint) {
142 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
143 const Evaluation& Ktmp = hint->fluidState().K(compIdx);
144 fluidState_.setKvalue(compIdx, Ktmp);
145 }
146 const Evaluation& Ltmp = hint->fluidState().L();
147 fluidState_.setLvalue(Ltmp);
148 }
149 else if (timeIdx == 0 && elemCtx.thermodynamicHint(dofIdx, 1)) {
150 // checking the storage cache
151 const auto& hint2 = elemCtx.thermodynamicHint(dofIdx, 1);
152 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
153 const Evaluation& Ktmp = hint2->fluidState().K(compIdx);
154 fluidState_.setKvalue(compIdx, Ktmp);
155 }
156 const Evaluation& Ltmp = hint2->fluidState().L();
157 fluidState_.setLvalue(Ltmp);
158 }
159 else {
160 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
161 const Evaluation Ktmp = fluidState_.wilsonK_(compIdx);
162 fluidState_.setKvalue(compIdx, Ktmp);
163 }
164 const Evaluation& Ltmp = -1.0;
165 fluidState_.setLvalue(Ltmp);
166 }
167
169 // Compute the phase compositions and densities
171 if (flashVerbosity >= 1) {
172 const int spatialIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
173 std::cout << " updating the intensive quantities for Cell " << spatialIdx << std::endl;
174 }
175 FlashSolver::solve(fluidState_, z, flashTwoPhaseMethod, flashTolerance, flashVerbosity);
176
177 if (flashVerbosity >= 5) {
178 // printing of flash result after solve
179 const int spatialIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
180 std::cout << " \n After flash solve for cell " << spatialIdx << std::endl;
181 ComponentVector x, y;
182 for (unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
183 x[comp_idx] = fluidState_.moleFraction(FluidSystem::oilPhaseIdx, comp_idx);
184 y[comp_idx] = fluidState_.moleFraction(FluidSystem::gasPhaseIdx, comp_idx);
185 }
186 for (unsigned comp_idx = 0; comp_idx < numComponents; ++comp_idx) {
187 std::cout << " x for component: " << comp_idx << " is:" << std::endl;
188 std::cout << x[comp_idx] << std::endl;
189
190 std::cout << " y for component: " << comp_idx << "is:" << std::endl;
191 std::cout << y[comp_idx] << std::endl;
192 }
193 const Evaluation& L = fluidState_.L();
194 std::cout << " L is:" << std::endl;
195 std::cout << L << std::endl;
196 }
197
198
199 // Update phases
200 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
201 paramCache.updatePhase(fluidState_, FluidSystem::oilPhaseIdx);
202
203 const Scalar R = Opm::Constants<Scalar>::R;
204 Evaluation Z_L = (paramCache.molarVolume(FluidSystem::oilPhaseIdx) * fluidState_.pressure(FluidSystem::oilPhaseIdx) )/
205 (R * fluidState_.temperature(FluidSystem::oilPhaseIdx));
206 paramCache.updatePhase(fluidState_, FluidSystem::gasPhaseIdx);
207 Evaluation Z_V = (paramCache.molarVolume(FluidSystem::gasPhaseIdx) * fluidState_.pressure(FluidSystem::gasPhaseIdx) )/
208 (R * fluidState_.temperature(FluidSystem::gasPhaseIdx));
209
210
211 // Update saturation
212 // \Note: the current implementation assume oil-gas system.
213 Evaluation L = fluidState_.L();
214 Evaluation So = Opm::max((L * Z_L / ( L * Z_L + (1 - L) * Z_V)), 0.0);
215 Evaluation Sg = Opm::max(1 - So, 0.0);
216 Scalar sumS = Opm::getValue(So) + Opm::getValue(Sg);
217 So /= sumS;
218 Sg /= sumS;
219
220 fluidState_.setSaturation(0, So);
221 fluidState_.setSaturation(1, Sg);
222
223 fluidState_.setCompressFactor(0, Z_L);
224 fluidState_.setCompressFactor(1, Z_V);
225
226 // Print saturation
227 if (flashVerbosity >= 5) {
228 std::cout << "So = " << So <<std::endl;
229 std::cout << "Sg = " << Sg <<std::endl;
230 }
231
232 // Print saturation
233 if (flashVerbosity >= 5) {
234 std::cout << "So = " << So <<std::endl;
235 std::cout << "Sg = " << Sg <<std::endl;
236 std::cout << "Z_L = " << Z_L <<std::endl;
237 std::cout << "Z_V = " << Z_V <<std::endl;
238 }
239
241 // Compute rel. perm and viscosity and densities
243 const MaterialLawParams& materialParams = problem.materialLawParams(elemCtx, dofIdx, timeIdx);
244
245 // calculate relative permeability
246 MaterialLaw::relativePermeabilities(relativePermeability_,
247 materialParams, fluidState_);
248 Opm::Valgrind::CheckDefined(relativePermeability_);
249
250 // set the phase viscosity and density
251 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
252 paramCache.updatePhase(fluidState_, phaseIdx);
253
254 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
255
256 fluidState_.setViscosity(phaseIdx, mu);
257
258 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
259 Opm::Valgrind::CheckDefined(mobility_[phaseIdx]);
260
261 const Evaluation& rho = FluidSystem::density(fluidState_, paramCache, phaseIdx);
262 fluidState_.setDensity(phaseIdx, rho);
263 }
264
266 // Compute the remaining quantities
268
269 // porosity
270 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
271 Opm::Valgrind::CheckDefined(porosity_);
272
273 // intrinsic permeability
274 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
275
276 // update the quantities specific for the velocity model
277 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
278
279 // energy related quantities
280 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
281
282 // update the diffusion specific quantities of the intensive quantities
283 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
284 }
285
289 const FluidState& fluidState() const
290 { return fluidState_; }
291
295 const DimMatrix& intrinsicPermeability() const
296 { return intrinsicPerm_; }
297
301 const Evaluation& relativePermeability(unsigned phaseIdx) const
302 { return relativePermeability_[phaseIdx]; }
303
307 const Evaluation& mobility(unsigned phaseIdx) const
308 {
309 return mobility_[phaseIdx];
310 }
311
315 const Evaluation& porosity() const
316 { return porosity_; }
317
318private:
319 DimMatrix intrinsicPerm_;
320 FluidState fluidState_;
321 Evaluation porosity_;
322 std::array<Evaluation,numPhases> relativePermeability_;
323 std::array<Evaluation,numPhases> mobility_;
324};
325
326} // namespace Opm
327
328#endif
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes.
Definition: diffusionmodule.hh:141
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:532
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:58
Opm::CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:91
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: ptflash/flashintensivequantities.hh:307
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: ptflash/flashintensivequantities.hh:295
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: ptflash/flashintensivequantities.hh:289
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: ptflash/flashintensivequantities.hh:315
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: ptflash/flashintensivequantities.hh:301
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: ptflash/flashintensivequantities.hh:105
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:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Declares the parameters for the compositional multi-phase model based on flash calculations.