ptflash/flashmodel.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 OPM_PTFLASH_MODEL_HH
29#define OPM_PTFLASH_MODEL_HH
30
31#include <opm/material/constraintsolvers/PTFlash.hpp>
32
33#include <opm/material/densead/Math.hpp>
34
35#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
36#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
37
40
45
50
57
58#include <sstream>
59#include <string>
60
61namespace Opm {
62
63template <class TypeTag>
64class FlashModel;
65
66}
67
68namespace Opm::Properties {
69
70namespace TTag {
72struct FlashModel { using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
73} // namespace TTag
74
76template<class TypeTag>
77struct LocalResidual<TypeTag, TTag::FlashModel>
79
81template<class TypeTag>
82struct NewtonMethod<TypeTag, TTag::FlashModel>
84
86template<class TypeTag>
87struct FlashSolver<TypeTag, TTag::FlashModel>
88{
89 using type = Opm::PTFlash<GetPropType<TypeTag, Properties::Scalar>,
91};
92
94template<class TypeTag>
95struct Model<TypeTag, TTag::FlashModel>
97
99template<class TypeTag>
100struct PrimaryVariables<TypeTag, TTag::FlashModel>
102
104template<class TypeTag>
105struct RateVector<TypeTag, TTag::FlashModel>
107
109template<class TypeTag>
110struct BoundaryRateVector<TypeTag, TTag::FlashModel>
112
114template<class TypeTag>
115struct IntensiveQuantities<TypeTag, TTag::FlashModel>
117
119template<class TypeTag>
120struct ExtensiveQuantities<TypeTag, TTag::FlashModel>
122
124template<class TypeTag>
125struct Indices<TypeTag, TTag::FlashModel>
126{ using type = Opm::FlashIndices<TypeTag, /*PVIdx=*/0>; };
127
128// disable molecular diffusion by default
129template<class TypeTag>
130struct EnableDiffusion<TypeTag, TTag::FlashModel>
131{ static constexpr bool value = false; };
132
134template<class TypeTag>
135struct EnableEnergy<TypeTag, TTag::FlashModel>
136{ static constexpr bool value = false; };
137
138} // namespace Opm::Properties
139
140namespace Opm {
141
184template <class TypeTag>
186 : public MultiPhaseBaseModel<TypeTag>
187{
188 using ParentType = MultiPhaseBaseModel<TypeTag>;
189
193
195
196 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
197 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
198 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
199
200
202
203public:
204 explicit FlashModel(Simulator& simulator)
205 : ParentType(simulator)
206 {}
207
211 static void registerParameters()
212 {
214
215 // register runtime parameters of the VTK output modules
218
219 if (enableDiffusion)
221
222 if (enableEnergy)
224
225 Parameters::Register<Parameters::FlashTolerance<Scalar>>
226 ("The maximum tolerance for the flash solver to "
227 "consider the solution converged");
228 Parameters::Register<Parameters::FlashVerbosity>
229 ("Flash solver verbosity level");
230 Parameters::Register<Parameters::FlashTwoPhaseMethod>
231 ("Method for solving vapor-liquid composition. Available options include: "
232 "ssi, newton, ssi+newton");
233
234 Parameters::SetDefault<Parameters::FlashTolerance<Scalar>>(1e-12);
235 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(true);
236
237 // since thermodynamic hints are basically free if the cache for intensive quantities is
238 // enabled, and this model usually shows quite a performance improvment if they are
239 // enabled, let's enable them by default.
240 Parameters::SetDefault<Parameters::EnableThermodynamicHints>(true);
241 }
242
246 std::string primaryVarName(unsigned pvIdx) const
247 {
248 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
249 if (!tmp.empty())
250 return tmp;
251
252 std::ostringstream oss;
253 if (Indices::z0Idx <= pvIdx && pvIdx < Indices::z0Idx + numComponents - 1)
254 oss << "z_," << FluidSystem::componentName(/*compIdx=*/pvIdx - Indices::z0Idx);
255 else if (pvIdx==Indices::pressure0Idx)
256 oss << "pressure_" << FluidSystem::phaseName(0);
257 else
258 assert(false);
259
260 return oss.str();
261 }
262
266 std::string eqName(unsigned eqIdx) const
267 {
268 const std::string& tmp = EnergyModule::eqName(eqIdx);
269 if (!tmp.empty())
270 return tmp;
271
272 std::ostringstream oss;
273 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
274 + numComponents) {
275 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
276 oss << "continuity^" << FluidSystem::componentName(compIdx);
277 }
278 else
279 assert(false);
280
281 return oss.str();
282 }
283
285 {
287
288 // add the VTK output modules which are meaningful for the model
289 this->addOutputModule(new Opm::VtkCompositionModule<TypeTag>(this->simulator_));
290 this->addOutputModule(new Opm::VtkPTFlashModule<TypeTag>(this->simulator_));
291 if (enableDiffusion)
292 this->addOutputModule(new Opm::VtkDiffusionModule<TypeTag>(this->simulator_));
293 if (enableEnergy)
294 this->addOutputModule(new Opm::VtkEnergyModule<TypeTag>(this->simulator_));
295 }
296};
297
298} // namespace Opm
299
300#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:45
This template class contains the data which is required to calculate all fluxes of components over a ...
Definition: flashextensivequantities.hh:54
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
Definition: flash/flashindices.hh:47
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flash/flashintensivequantities.hh:58
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition: flash/flashlocalresidual.hh:46
A compositional multi-phase model based on flash-calculations.
Definition: ptflash/flashmodel.hh:187
FlashModel(Simulator &simulator)
Definition: ptflash/flashmodel.hh:204
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: ptflash/flashmodel.hh:211
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: ptflash/flashmodel.hh:246
void registerOutputModules_()
Definition: ptflash/flashmodel.hh:284
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: ptflash/flashmodel.hh:266
A Newton solver specific to the PTFlash model.
Definition: flashnewtonmethod.hh:53
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:57
Definition: flashratevector.hh:48
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:153
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:179
void registerOutputModules_()
Definition: multiphasebasemodel.hh:254
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:69
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:96
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hh:65
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:92
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:66
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:93
VTK output module for the PT Flash calculation This module deals with the following quantities: K,...
Definition: vtkptflashmodule.hh:60
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkptflashmodule.hh:87
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declares the properties required by the compositional multi-phase model based on flash calculations.
Definition: blackoilmodel.hh:72
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
Declares the parameters for the compositional multi-phase model based on flash calculations.
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:119
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:79
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:149
Opm::NcpFlash< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem > > type
Definition: flash/flashmodel.hh:78
The type of the flash constraint solver.
Definition: flashproperties.hh:39
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:48
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:133
The type of the local residual function.
Definition: fvbaseproperties.hh:94
The type of the model.
Definition: basicproperties.hh:88
Specifies the type of the actual Newton method.
Definition: newtonmethodproperties.hh:32
A vector of primary variables within a sub-control volume.
Definition: fvbaseproperties.hh:130
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:116
The type tag for the isothermal single phase problems.
Definition: flash/flashmodel.hh:67
std::tuple< MultiPhaseBaseModel > InheritsFrom
Definition: flash/flashmodel.hh:67