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
44#include <array>
45
46namespace Opm {
47
54template <class TypeTag>
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 { cTot0Idx = Indices::cTot0Idx };
73 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
74 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
75 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
76 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
77 enum { dimWorld = GridView::dimensionworld };
78
83
84 using ComponentVector = Dune::FieldVector<Evaluation, numComponents>;
85 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
86
87 using FluxIntensiveQuantities = typename FluxModule::FluxIntensiveQuantities;
90
91public:
93 using FluidState = CompositionalFluidState<Evaluation, FluidSystem, enableEnergy>;
94
96
98
100
104 void update(const ElementContext& elemCtx, unsigned dofIdx, unsigned timeIdx)
105 {
106 ParentType::update(elemCtx, dofIdx, timeIdx);
107 EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
108
109 const auto& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
110 const auto& problem = elemCtx.problem();
111 const Scalar flashTolerance = Parameters::Get<Parameters::FlashTolerance<Scalar>>();
112
113 // extract the total molar densities of the components
114 ComponentVector cTotal;
115 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
116 cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
117 }
118
119 const auto* hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
120 if (hint) {
121 // use the same fluid state as the one of the hint, but
122 // make sure that we don't overwrite the temperature
123 // specified by the primary variables
124 const Evaluation T = fluidState_.temperature(/*phaseIdx=*/0);
125 fluidState_.assign(hint->fluidState());
126 fluidState_.setTemperature(T);
127 }
128 else {
129 FlashSolver::guessInitial(fluidState_, cTotal);
130 }
131
132 // compute the phase compositions, densities and pressures
133 typename FluidSystem::template ParameterCache<Evaluation> paramCache;
134 const MaterialLawParams& materialParams =
135 problem.materialLawParams(elemCtx, dofIdx, timeIdx);
136 FlashSolver::template solve<MaterialLaw>(fluidState_,
137 materialParams,
138 paramCache,
139 cTotal,
140 flashTolerance);
141
142 // calculate relative permeabilities
143 MaterialLaw::relativePermeabilities(relativePermeability_,
144 materialParams, fluidState_);
145 Valgrind::CheckDefined(relativePermeability_);
146
147 // set the phase viscosities
148 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
149 paramCache.updatePhase(fluidState_, phaseIdx);
150
151 const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
152 fluidState_.setViscosity(phaseIdx, mu);
153
154 mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
155 Valgrind::CheckDefined(mobility_[phaseIdx]);
156 }
157
159 // calculate the remaining quantities
161
162 // porosity
163 porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
164 Valgrind::CheckDefined(porosity_);
165
166 // intrinsic permeability
167 intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
168
169 // update the quantities specific for the velocity model
170 FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
171
172 // energy related quantities
173 EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
174
175 // update the diffusion specific quantities of the intensive quantities
176 DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
177 }
178
182 const FluidState& fluidState() const
183 { return fluidState_; }
184
188 const DimMatrix& intrinsicPermeability() const
189 { return intrinsicPerm_; }
190
194 const Evaluation& relativePermeability(unsigned phaseIdx) const
195 { return relativePermeability_[phaseIdx]; }
196
200 const Evaluation& mobility(unsigned phaseIdx) const
201 { return mobility_[phaseIdx]; }
202
206 const Evaluation& porosity() const
207 { return porosity_; }
208
209private:
210 DimMatrix intrinsicPerm_;
211 FluidState fluidState_;
212 Evaluation porosity_;
213 std::array<Evaluation,numPhases> relativePermeability_;
214 std::array<Evaluation,numPhases> mobility_;
215};
216
217} // namespace Opm
218
219#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 intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:60
const Evaluation & mobility(unsigned phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:200
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flash/flashintensivequantities.hh:188
FlashIntensiveQuantities(const FlashIntensiveQuantities &other)=default
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: flash/flashintensivequantities.hh:182
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flash/flashintensivequantities.hh:206
const Evaluation & relativePermeability(unsigned phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flash/flashintensivequantities.hh:194
CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flash/flashintensivequantities.hh:93
void update(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: flash/flashintensivequantities.hh:104
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: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