ncpintensivequantities.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_NCP_INTENSIVE_QUANTITIES_HH
29#define EWOMS_NCP_INTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fvector.hh>
32#include <dune/common/fmatrix.hh>
33
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/constraintsolvers/CompositionFromFugacities.hpp>
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38
43
44#include <array>
45
46namespace Opm {
47
55template <class TypeTag>
57 : public GetPropType<TypeTag, Properties::DiscIntensiveQuantities>
58 , public DiffusionIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableDiffusion>() >
59 , public EnergyIntensiveQuantities<TypeTag, getPropValue<TypeTag, Properties::EnableEnergy>() >
60 , public GetPropType<TypeTag, Properties::FluxModule>::FluxIntensiveQuantities
61{
63
71
74
75 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
76 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
77 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
78 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
79 enum { fugacity0Idx = Indices::fugacity0Idx };
80 enum { saturation0Idx = Indices::saturation0Idx };
81 enum { pressure0Idx = Indices::pressure0Idx };
82 enum { dimWorld = GridView::dimensionworld };
83
84 using CompositionFromFugacitiesSolver = ::Opm::CompositionFromFugacities<Scalar, FluidSystem, Evaluation>;
85 using FluidState = CompositionalFluidState<Evaluation, FluidSystem, /*storeEnthalpy=*/enableEnergy>;
86 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
87 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
90 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
91
92public:
94
96
98
102 void update(const ElementContext& elemCtx,
103 unsigned dofIdx,
104 unsigned timeIdx)
105 {
106 ParentType::update(elemCtx, dofIdx, timeIdx);
107 ParentType::checkDefined();
108
109 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
110 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
111
112 // set the phase saturations
113 Evaluation sumSat = 0;
114 for (unsigned phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
115 const Evaluation& val = priVars.makeEvaluation(saturation0Idx + phaseIdx, timeIdx);
116 fluidState_.setSaturation(phaseIdx, val);
117 sumSat += val;
118 }
119 fluidState_.setSaturation(numPhases - 1, 1.0 - sumSat);
120 Opm::Valgrind::CheckDefined(sumSat);
121
122 // set the fluid phase temperature
123 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
124
125 // retrieve capillary pressure parameters
126 const auto& problem = elemCtx.problem();
127 const MaterialLawParams& materialParams =
128 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
129
130 // calculate capillary pressures
131 std::array<Evaluation, numPhases> capPress;
132 MaterialLaw::capillaryPressures(capPress, materialParams, fluidState_);
133
134 // add to the pressure of the first fluid phase
135 const Evaluation& pressure0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
136 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
137 fluidState_.setPressure(phaseIdx, pressure0 + (capPress[phaseIdx] - capPress[0]));
138 }
139
140 ComponentVector fug;
141 // retrieve component fugacities
142 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
143 fug[compIdx] = priVars.makeEvaluation(fugacity0Idx + compIdx, timeIdx);
144 }
145
146 // calculate phase compositions
147 const auto* hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
148 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
149 // initial guess
150 if (hint) {
151 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
152 // use the hint for the initial mole fraction!
153 const Evaluation& moleFracIJ = hint->fluidState().moleFraction(phaseIdx, compIdx);
154 fluidState_.setMoleFraction(phaseIdx, compIdx, moleFracIJ);
155 }
156 }
157 else { // !hint
158 CompositionFromFugacitiesSolver::guessInitial(fluidState_, phaseIdx, fug);
159 }
160
161 // calculate the phase composition from the component
162 // fugacities
163 CompositionFromFugacitiesSolver::solve(fluidState_, paramCache, phaseIdx, fug);
164 }
165
166 // porosity
167 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
168 Opm::Valgrind::CheckDefined(porosity_);
169
170 // relative permeabilities
171 MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
172
173 // dynamic viscosities
174 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
175 // viscosities
176 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
177 fluidState_.setViscosity(phaseIdx, mu);
178
179 mobility_[phaseIdx] = relativePermeability_[phaseIdx]/mu;
180 }
181
182 // intrinsic permeability
183 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
184
185 // update the quantities specific for the velocity model
186 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
187
188 // energy related quantities
189 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
190
191 // update the diffusion specific quantities of the intensive quantities
192 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
193
194 checkDefined();
195 }
196
200 const FluidState& fluidState() const
201 { return fluidState_; }
202
206 const DimMatrix& intrinsicPermeability() const
207 { return intrinsicPerm_; }
208
212 const Evaluation& relativePermeability(unsigned phaseIdx) const
213 { return relativePermeability_[phaseIdx]; }
214
218 const Evaluation& mobility(unsigned phaseIdx) const
219 { return mobility_[phaseIdx]; }
220
224 const Evaluation& porosity() const
225 { return porosity_; }
226
230 void checkDefined() const
231 {
232#if !defined NDEBUG && HAVE_VALGRIND
233 ParentType::checkDefined();
234
235 Valgrind::CheckDefined(porosity_);
236 Valgrind::CheckDefined(relativePermeability_);
237
238 fluidState_.checkDefined();
239#endif
240 }
241
242private:
243 DimMatrix intrinsicPerm_;
244 FluidState fluidState_;
245 Evaluation porosity_;
246 std::array<Evaluation, numPhases> relativePermeability_{};
247 std::array<Evaluation, numPhases> mobility_{};
248};
249
250} // namespace Opm
251
252#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: ncpintensivequantities.hh:61
void checkDefined() const
IntensiveQuantities::checkDefined.
Definition: ncpintensivequantities.hh:230
const FluidState & fluidState() const
ImmiscibleIntensiveQuantities::fluidState.
Definition: ncpintensivequantities.hh:200
const Evaluation & mobility(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::mobility.
Definition: ncpintensivequantities.hh:218
NcpIntensiveQuantities & operator=(const NcpIntensiveQuantities &other)=default
NcpIntensiveQuantities(const NcpIntensiveQuantities &other)=default
const Evaluation & porosity() const
ImmiscibleIntensiveQuantities::porosity.
Definition: ncpintensivequantities.hh:224
const Evaluation & relativePermeability(unsigned phaseIdx) const
ImmiscibleIntensiveQuantities::relativePermeability.
Definition: ncpintensivequantities.hh:212
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
IntensiveQuantities::update.
Definition: ncpintensivequantities.hh:102
const DimMatrix & intrinsicPermeability() const
ImmiscibleIntensiveQuantities::intrinsicPermeability.
Definition: ncpintensivequantities.hh:206
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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