flash/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 EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
29#define EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <opm/material/common/Valgrind.hpp>
35#include <opm/material/fluidstates/CompositionalFluidState.hpp>
36
39
43
44namespace Opm {
45
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
69 // primary variable indices
70 enum { cTot0Idx = Indices::cTot0Idx };
71 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
72 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
73 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
74 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
75 enum { dimWorld = GridView::dimensionworld };
76
81
82 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
83 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
84
85 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
88
89public:
91 using FluidState = Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
92
94 { }
95
97
99
103 void update(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
104 {
105 ParentType::update(elemCtx, dofIdx, timeIdx);
106 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
107
108 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
109 const auto& problem = elemCtx.problem();
110 Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
111
112 // extract the total molar densities of the components
113 ComponentVector cTotal;
114 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
115 cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
116
117 const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
118 if (hint) {
119 // use the same fluid state as the one of the hint, but
120 // make sure that we don't overwrite the temperature
121 // specified by the primary variables
122 Evaluation T = fluidState_.temperature(/*phaseIdx=*/0);
123 fluidState_.assign(hint->fluidState());
124 fluidState_.setTemperature(T);
125 }
126 else
127 FlashSolver::guessInitial(fluidState_, cTotal);
128
129 // compute the phase compositions, densities and pressures
130 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
131 const MaterialLawParams& materialParams =
132 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
133 FlashSolver::template solve<MaterialLaw>(fluidState_,
134 materialParams,
135 paramCache,
136 cTotal,
137 flashTolerance);
138
139 // calculate relative permeabilities
140 MaterialLaw::relativePermeabilities(relativePermeability_,
141 materialParams, fluidState_);
142 Opm::Valgrind::CheckDefined(relativePermeability_);
143
144 // set the phase viscosities
145 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
146 paramCache.updatePhase(fluidState_, phaseIdx);
147
148 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
149 fluidState_.setViscosity(phaseIdx, mu);
150
151 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
152 Opm::Valgrind::CheckDefined(mobility_[phaseIdx]);
153 }
154
156 // calculate the remaining quantities
158
159 // porosity
160 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
161 Opm::Valgrind::CheckDefined(porosity_);
162
163 // intrinsic permeability
164 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
165
166 // update the quantities specific for the velocity model
167 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
168
169 // energy related quantities
170 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
171
172 // update the diffusion specific quantities of the intensive quantities
173 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
174 }
175
179 const FluidState& fluidState() const
180 { return fluidState_; }
181
185 const DimMatrix& intrinsicPermeability() const
186 { return intrinsicPerm_; }
187
191 const Evaluation& relativePermeability(unsigned phaseIdx) const
192 { return relativePermeability_[phaseIdx]; }
193
197 const Evaluation& mobility(unsigned phaseIdx) const
198 {
199 return mobility_[phaseIdx];
200 }
201
205 const Evaluation& porosity() const
206 { return porosity_; }
207
208private:
209 DimMatrix intrinsicPerm_;
210 FluidState fluidState_;
211 Evaluation porosity_;
212 std::array<Evaluation,numPhases> relativePermeability_;
213 std::array<Evaluation,numPhases> mobility_;
214};
215
216} // namespace Opm
217
218#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: flash/flashintensivequantities.hh:197
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flash/flashintensivequantities.hh:185
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flash/flashintensivequantities.hh:179
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flash/flashintensivequantities.hh:205
FlashIntensiveQuantities()
Definition: flash/flashintensivequantities.hh:93
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:191
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flash/flashintensivequantities.hh:103
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 parameters for the compositional multi-phase model based on flash calculations.
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