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