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