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 <cassert>
59#include <memory>
60#include <sstream>
61#include <string>
62#include <tuple>
63
64namespace Opm {
65
66template <class TypeTag>
67class FlashModel;
68
69}
70
71namespace Opm::Properties {
72
73namespace TTag {
75struct FlashModel { using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
76} // namespace TTag
77
79template<class TypeTag>
80struct LocalResidual<TypeTag, TTag::FlashModel>
82
84template<class TypeTag>
85struct NewtonMethod<TypeTag, TTag::FlashModel>
87
89template<class TypeTag>
90struct FlashSolver<TypeTag, TTag::FlashModel>
91{
92 using type = PTFlash<GetPropType<TypeTag, Properties::Scalar>,
94};
95
97template<class TypeTag>
98struct Model<TypeTag, TTag::FlashModel>
100
102template<class TypeTag>
103struct PrimaryVariables<TypeTag, TTag::FlashModel>
105
107template<class TypeTag>
108struct RateVector<TypeTag, TTag::FlashModel>
110
112template<class TypeTag>
113struct BoundaryRateVector<TypeTag, TTag::FlashModel>
115
117template<class TypeTag>
118struct IntensiveQuantities<TypeTag, TTag::FlashModel>
120
122template<class TypeTag>
123struct ExtensiveQuantities<TypeTag, TTag::FlashModel>
125
127template<class TypeTag>
128struct Indices<TypeTag, TTag::FlashModel>
129{ using type = FlashIndices<TypeTag, /*PVIdx=*/0>; };
130
131// disable molecular diffusion by default
132template<class TypeTag>
133struct EnableDiffusion<TypeTag, TTag::FlashModel>
134{ static constexpr bool value = false; };
135
137template<class TypeTag>
138struct EnableEnergy<TypeTag, TTag::FlashModel>
139{ static constexpr bool value = false; };
140
141template<class TypeTag>
142struct EnableWater<TypeTag, TTag::MultiPhaseBaseModel>
144
145} // namespace Opm::Properties
146
147namespace Opm {
148
191template <class TypeTag>
193 : public MultiPhaseBaseModel<TypeTag>
194{
195 using ParentType = MultiPhaseBaseModel<TypeTag>;
196
200
202
203 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
204 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
205 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
206
207
209
210public:
211 explicit FlashModel(Simulator& simulator)
212 : ParentType(simulator)
213 {}
214
218 static void registerParameters()
219 {
221
222 // register runtime parameters of the VTK output modules
225
226 if constexpr (enableDiffusion) {
228 }
229
230 if constexpr (enableEnergy) {
232 }
233
234 Parameters::Register<Parameters::FlashTolerance<Scalar>>
235 ("The maximum tolerance for the flash solver to "
236 "consider the solution converged");
237 Parameters::Register<Parameters::FlashVerbosity>
238 ("Flash solver verbosity level");
239 Parameters::Register<Parameters::FlashTwoPhaseMethod>
240 ("Method for solving vapor-liquid composition. Available options include: "
241 "ssi, newton, ssi+newton");
242
243 Parameters::SetDefault<Parameters::FlashTolerance<Scalar>>(1.e-8);
244 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(true);
245
246 // since thermodynamic hints are basically free if the cache for intensive quantities is
247 // enabled, and this model usually shows quite a performance improvment if they are
248 // enabled, let's enable them by default.
249 Parameters::SetDefault<Parameters::EnableThermodynamicHints>(true);
250 }
251
255 std::string primaryVarName(unsigned pvIdx) const
256 {
257 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
258 if (!tmp.empty()) {
259 return tmp;
260 }
261
262 std::ostringstream oss;
263 if (Indices::z0Idx <= pvIdx && pvIdx < Indices::z0Idx + numComponents - 1) {
264 oss << "z_," << FluidSystem::componentName(/*compIdx=*/pvIdx - Indices::z0Idx);
265 }
266 else if (pvIdx == Indices::pressure0Idx) {
267 oss << "pressure_" << FluidSystem::phaseName(0);
268 }
269 else {
270 assert(false);
271 }
272
273 return oss.str();
274 }
275
279 std::string eqName(unsigned eqIdx) const
280 {
281 const std::string& tmp = EnergyModule::eqName(eqIdx);
282 if (!tmp.empty()) {
283 return tmp;
284 }
285
286 std::ostringstream oss;
287 if (Indices::conti0EqIdx <= eqIdx &&
288 eqIdx < Indices::conti0EqIdx + numComponents)
289 {
290 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
291 oss << "continuity^" << FluidSystem::componentName(compIdx);
292 }
293 else {
294 assert(false);
295 }
296
297 return oss.str();
298 }
299
301 {
303
304 // add the VTK output modules which are meaningful for the model
305 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
306 this->addOutputModule(std::make_unique<VtkPTFlashModule<TypeTag>>(this->simulator_));
307 if constexpr (enableDiffusion) {
308 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
309 }
310 if constexpr (enableEnergy) {
311 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
312 }
313 }
314};
315
316} // namespace Opm
317
318#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:49
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:60
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition: flash/flashlocalresidual.hh:47
A compositional multi-phase model based on flash-calculations.
Definition: ptflash/flashmodel.hh:194
FlashModel(Simulator &simulator)
Definition: ptflash/flashmodel.hh:211
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: ptflash/flashmodel.hh:218
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: ptflash/flashmodel.hh:255
void registerOutputModules_()
Definition: ptflash/flashmodel.hh:300
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: ptflash/flashmodel.hh:279
A Newton solver specific to the PTFlash model.
Definition: flashnewtonmethod.hh:55
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:52
Definition: flashratevector.hh:48
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:168
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:194
void registerOutputModules_()
Definition: multiphasebasemodel.hh:270
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hpp:57
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hpp:87
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hpp:58
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hpp:88
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hpp:58
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hpp:87
VTK output module for the PT Flash calculation This module deals with the following quantities: K,...
Definition: vtkptflashmodule.hpp:53
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkptflashmodule.hpp:83
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:79
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
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:91
Enable water?
Definition: multiphasebaseproperties.hh:103
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:149
NcpFlash< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem > > type
Definition: flash/flashmodel.hh:88
The type of the flash constraint solver.
Definition: flashproperties.hh:39
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:51
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:74
std::tuple< MultiPhaseBaseModel > InheritsFrom
Definition: flash/flashmodel.hh:74