pvsintensivequantities.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_PVS_INTENSIVE_QUANTITIES_HH
29#define EWOMS_PVS_INTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <opm/material/common/MathToolbox.hpp>
35#include <opm/material/common/Valgrind.hpp>
36#include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
37#include <opm/material/constraintsolvers/MiscibleMultiPhaseComposition.hpp>
38#include <opm/material/fluidstates/CompositionalFluidState.hpp>
39
42
43#include <array>
44#include <limits>
45
46namespace Opm {
47
56template <class TypeTag>
58 : public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
59 , public DiffusionIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableDiffusion>() >
60 , public EnergyIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableEnergy>() >
61 , public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
62{
64
74
75 enum { switch0Idx = Indices::switch0Idx };
76 enum { pressure0Idx = Indices::pressure0Idx };
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
83 using Toolbox = MathToolbox<Evaluation>;
84 using MiscibleMultiPhaseComposition = ::Opm::MiscibleMultiPhaseComposition<Scalar, FluidSystem, Evaluation>;
85 using ComputeFromReferencePhase = ::Opm::ComputeFromReferencePhase<Scalar, FluidSystem, Evaluation>;
86
87 using PhaseVector = Dune::FieldVector<Scalar, numPhases>;
88 using EvalPhaseVector = Dune::FieldVector<Evaluation, numPhases>;
89 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
90
93 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
94
95public:
97 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem>;
98
100
102
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
117 // set the saturations
119 Evaluation sumSat = 0.0;
120 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
121 fluidState_.setSaturation(phaseIdx, priVars.explicitSaturationValue(phaseIdx, timeIdx));
122 Valgrind::CheckDefined(fluidState_.saturation(phaseIdx));
123 sumSat += fluidState_.saturation(phaseIdx);
124 }
125 Valgrind::CheckDefined(priVars.implicitSaturationIdx());
126 Valgrind::CheckDefined(sumSat);
127 fluidState_.setSaturation(priVars.implicitSaturationIdx(), 1.0 - sumSat);
128
130 // set the pressures of the fluid phases
132
133 // calculate capillary pressure
134 const MaterialLawParams& materialParams =
135 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
136 EvalPhaseVector pC;
137 MaterialLaw::capillaryPressures(pC, materialParams, fluidState_);
138
139 // set the absolute phase pressures in the fluid state
140 const Evaluation& p0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
141 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
142 fluidState_.setPressure(phaseIdx, p0 + (pC[phaseIdx] - pC[0]));
143 }
144
146 // calculate the phase compositions
148
149 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
150 unsigned lowestPresentPhaseIdx = priVars.lowestPresentPhaseIdx();
151 unsigned numNonPresentPhases = 0;
152 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
153 if (!priVars.phaseIsPresent(phaseIdx)) {
154 ++numNonPresentPhases;
155 }
156 }
157
158 // now comes the tricky part: calculate phase compositions
159 if (numNonPresentPhases == numPhases - 1) {
160 // only one phase is present, i.e. the primary variables
161 // contain the complete composition of the phase
162 Evaluation sumx = 0.0;
163 for (unsigned compIdx = 1; compIdx < numComponents; ++compIdx) {
164 const Evaluation& x = priVars.makeEvaluation(switch0Idx + compIdx - 1, timeIdx);
165 fluidState_.setMoleFraction(lowestPresentPhaseIdx, compIdx, x);
166 sumx += x;
167 }
168
169 // set the mole fraction of the first component
170 fluidState_.setMoleFraction(lowestPresentPhaseIdx, 0, 1 - sumx);
171
172 // calculate the composition of the remaining phases (as
173 // well as the densities of all phases). this is the job
174 // of the "ComputeFromReferencePhase" constraint solver
175 ComputeFromReferencePhase::solve(fluidState_, paramCache,
176 lowestPresentPhaseIdx,
177 /*setViscosity=*/true,
178 /*setEnthalpy=*/false);
179 }
180 else {
181 // create the auxiliary constraints
182 unsigned numAuxConstraints = numComponents + numNonPresentPhases - numPhases;
183 Opm::MMPCAuxConstraint<Evaluation> auxConstraints[numComponents];
184
185 unsigned auxIdx = 0;
186 unsigned switchIdx = 0;
187 for (; switchIdx < numPhases - 1; ++switchIdx) {
188 const unsigned compIdx = switchIdx + 1;
189 unsigned switchPhaseIdx = switchIdx;
190 if (switchIdx >= lowestPresentPhaseIdx)
191 switchPhaseIdx += 1;
192
193 if (!priVars.phaseIsPresent(switchPhaseIdx)) {
194 auxConstraints[auxIdx].set(lowestPresentPhaseIdx, compIdx,
195 priVars.makeEvaluation(switch0Idx + switchIdx, timeIdx));
196 ++auxIdx;
197 }
198 }
199
200 for (; auxIdx < numAuxConstraints; ++auxIdx, ++switchIdx) {
201 const unsigned compIdx = numPhases - numNonPresentPhases + auxIdx;
202 auxConstraints[auxIdx].set(lowestPresentPhaseIdx, compIdx,
203 priVars.makeEvaluation(switch0Idx + switchIdx, timeIdx));
204 }
205
206 // both phases are present, i.e. phase compositions are a result of the the
207 // gas <-> liquid equilibrium. This is the job of the
208 // "MiscibleMultiPhaseComposition" constraint solver
209 MiscibleMultiPhaseComposition::solve(fluidState_, paramCache,
210 priVars.phasePresence(),
211 auxConstraints,
212 numAuxConstraints,
213 /*setViscosity=*/true,
214 /*setEnthalpy=*/false);
215 }
216
217#ifndef NDEBUG
218 // make valgrind happy and set the enthalpies to NaN
219 if constexpr (!enableEnergy) {
220 constexpr Scalar myNan = std::numeric_limits<Scalar>::quiet_NaN();
221 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
222 fluidState_.setEnthalpy(phaseIdx, myNan);
223 }
224 }
225#endif
226
228 // calculate the remaining quantities
230
231 // calculate relative permeabilities
232 MaterialLaw::relativePermeabilities(relativePermeability_,
233 materialParams, fluidState_);
234 Valgrind::CheckDefined(relativePermeability_);
235
236 // mobilities
237 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
238 mobility_[phaseIdx] =
239 relativePermeability_[phaseIdx] / fluidState().viscosity(phaseIdx);
240 }
241
242 // porosity
243 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
244 Valgrind::CheckDefined(porosity_);
245
246 // intrinsic permeability
247 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
248
249 // update the quantities specific for the velocity model
250 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
251
252 // energy related quantities
253 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
254
255 // update the diffusion specific quantities of the intensive quantities
256 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
257
258 fluidState_.checkDefined();
259 }
260
264 const FluidState& fluidState() const
265 { return fluidState_; }
266
270 const DimMatrix& intrinsicPermeability() const
271 { return intrinsicPerm_; }
272
276 const Evaluation& relativePermeability(unsigned phaseIdx) const
277 { return relativePermeability_[phaseIdx]; }
278
282 const Evaluation& mobility(unsigned phaseIdx) const
283 { return mobility_[phaseIdx]; }
284
288 const Evaluation& porosity() const
289 { return porosity_; }
290
291private:
292 FluidState fluidState_;
293 Evaluation porosity_;
294 DimMatrix intrinsicPerm_;
295 std::array<Evaluation, numPhases> relativePermeability_{};
296 std::array<Evaluation, numPhases> mobility_{};
297};
298
299} // namespace Opm
300
301#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 quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: pvsintensivequantities.hh:62
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: pvsintensivequantities.hh:288
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: pvsintensivequantities.hh:276
Opm::CompositionalFluidState< Evaluation, FluidSystem > FluidState
The type of the object returned by the fluidState() method.
Definition: pvsintensivequantities.hh:97
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: pvsintensivequantities.hh:108
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: pvsintensivequantities.hh:282
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: pvsintensivequantities.hh:270
PvsIntensiveQuantities & operator=(const PvsIntensiveQuantities &other)=default
PvsIntensiveQuantities(const PvsIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: pvsintensivequantities.hh:264
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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