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/densead/Math.hpp>
32
33#include "flashproperties.hh"
35#include "flashlocalresidual.hh"
40#include "flashindices.hh"
41#include "flashnewtonmethod.hh"
42
49#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
50#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
51#include <opm/material/constraintsolvers/PTFlash.hpp>
52
53#include <sstream>
54#include <string>
55
56namespace Opm {
57template <class TypeTag>
58class FlashModel;
59}
60
61namespace Opm::Properties {
62
63namespace TTag {
65struct FlashModel { using InheritsFrom = std::tuple<VtkDiffusion,
70} // namespace TTag
71
73template<class TypeTag>
75
77template<class TypeTag>
78struct NewtonMethod<TypeTag, TTag::FlashModel> { using type = Opm::FlashNewtonMethod<TypeTag>; };
79
81template<class TypeTag>
82struct FlashSolver<TypeTag, TTag::FlashModel>
83{ using type = Opm::PTFlash<GetPropType<TypeTag, Properties::Scalar>,
85
87template<class TypeTag>
88struct FlashTolerance<TypeTag, TTag::FlashModel>
89{
91 static constexpr type value = 1.e-12;
92};
93
94// Flash solver verbosity
95template<class TypeTag>
96struct FlashVerbosity<TypeTag, TTag::FlashModel> { static constexpr int value = 0; };
97
98// Flash two-phase method
99template<class TypeTag>
100struct FlashTwoPhaseMethod<TypeTag, TTag::FlashModel> { static constexpr auto value = "ssi"; };
101
103template<class TypeTag>
104struct Model<TypeTag, TTag::FlashModel> { using type = Opm::FlashModel<TypeTag>; };
105
107template<class TypeTag>
109
111template<class TypeTag>
112struct RateVector<TypeTag, TTag::FlashModel> { using type = Opm::FlashRateVector<TypeTag>; };
113
115template<class TypeTag>
117
119template<class TypeTag>
121
123template<class TypeTag>
125
127template<class TypeTag>
128struct Indices<TypeTag, TTag::FlashModel> { using type = Opm::FlashIndices<TypeTag, /*PVIdx=*/0>; };
129
130// The updates of intensive quantities tend to be _very_ expensive for this
131// model, so let's try to minimize the number of required ones
132template<class TypeTag>
133struct EnableIntensiveQuantityCache<TypeTag, TTag::FlashModel> { static constexpr bool value = true; };
134
135// since thermodynamic hints are basically free if the cache for intensive quantities is
136// enabled, and this model usually shows quite a performance improvment if they are
137// enabled, let's enable them by default.
138template<class TypeTag>
139struct EnableThermodynamicHints<TypeTag, TTag::FlashModel> { static constexpr bool value = true; };
140
141// disable molecular diffusion by default
142template<class TypeTag>
143struct EnableDiffusion<TypeTag, TTag::FlashModel> { static constexpr bool value = false; };
144
146template<class TypeTag>
147struct EnableEnergy<TypeTag, TTag::FlashModel> { static constexpr bool value = false; };
148
149} // namespace Opm::Properties
150
151namespace Opm {
152
195template <class TypeTag>
197 : public MultiPhaseBaseModel<TypeTag>
198{
199 using ParentType = MultiPhaseBaseModel<TypeTag>;
200
204
206
207 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
208 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
209 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
210
211
213
214public:
215 explicit FlashModel(Simulator& simulator)
216 : ParentType(simulator)
217 {}
218
222 static void registerParameters()
223 {
225
226 // register runtime parameters of the VTK output modules
229
230 if (enableDiffusion)
232
233 if (enableEnergy)
235
236 Parameters::registerParam<TypeTag, Properties::FlashTolerance>
237 ("The maximum tolerance for the flash solver to "
238 "consider the solution converged");
239 Parameters::registerParam<TypeTag, Properties::FlashVerbosity>
240 ("Flash solver verbosity level");
241 Parameters::registerParam<TypeTag, Properties::FlashTwoPhaseMethod>
242 ("Method for solving vapor-liquid composition. Available options include: "
243 "ssi, newton, ssi+newton");
244 }
245
249 std::string primaryVarName(unsigned pvIdx) const
250 {
251 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
252 if (!tmp.empty())
253 return tmp;
254
255 std::ostringstream oss;
256 if (Indices::z0Idx <= pvIdx && pvIdx < Indices::z0Idx + numComponents - 1)
257 oss << "z_," << FluidSystem::componentName(/*compIdx=*/pvIdx - Indices::z0Idx);
258 else if (pvIdx==Indices::pressure0Idx)
259 oss << "pressure_" << FluidSystem::phaseName(0);
260 else
261 assert(false);
262
263 return oss.str();
264 }
265
269 std::string eqName(unsigned eqIdx) const
270 {
271 const std::string& tmp = EnergyModule::eqName(eqIdx);
272 if (!tmp.empty())
273 return tmp;
274
275 std::ostringstream oss;
276 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
277 + numComponents) {
278 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
279 oss << "continuity^" << FluidSystem::componentName(compIdx);
280 }
281 else
282 assert(false);
283
284 return oss.str();
285 }
286
288 {
290
291 // add the VTK output modules which are meaningful for the model
292 this->addOutputModule(new Opm::VtkCompositionModule<TypeTag>(this->simulator_));
293 this->addOutputModule(new Opm::VtkPTFlashModule<TypeTag>(this->simulator_));
294 if (enableDiffusion)
295 this->addOutputModule(new Opm::VtkDiffusionModule<TypeTag>(this->simulator_));
296 if (enableEnergy)
297 this->addOutputModule(new Opm::VtkEnergyModule<TypeTag>(this->simulator_));
298 }
299};
300
301} // namespace Opm
302
303#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:48
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:57
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:198
FlashModel(Simulator &simulator)
Definition: ptflash/flashmodel.hh:215
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: ptflash/flashmodel.hh:222
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: ptflash/flashmodel.hh:249
void registerOutputModules_()
Definition: ptflash/flashmodel.hh:287
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: ptflash/flashmodel.hh:269
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:97
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:124
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hh:82
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:109
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:85
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:112
VTK output module for the PT Flash calculation This module deals with the following quantities: K,...
Definition: vtkptflashmodule.hh:73
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkptflashmodule.hh:100
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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:242
Declares the properties required by the compositional multi-phase model based on flash calculations.
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:136
Specify whether all intensive quantities for the grid should be cached in the discretization.
Definition: fvbaseproperties.hh:297
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:166
Opm::NcpFlash< GetPropType< TypeTag, Properties::Scalar >, GetPropType< TypeTag, Properties::FluidSystem > > type
Definition: flash/flashmodel.hh:77
The type of the flash constraint solver.
Definition: flash/flashproperties.hh:42
GetPropType< TypeTag, Scalar > type
Definition: flash/flashmodel.hh:83
The maximum accepted error of the flash solver.
Definition: flash/flashproperties.hh:45
Two-phase flash method.
Definition: ptflash/flashproperties.hh:51
The verbosity level of the flash solver.
Definition: ptflash/flashproperties.hh:48
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:48
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:150
The type of the local residual function.
Definition: fvbaseproperties.hh:111
The type of the model.
Definition: basicproperties.hh:81
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:147
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:133
The type tag for the isothermal single phase problems.
Definition: flash/flashmodel.hh:63
std::tuple< VtkDiffusion, VtkEnergy, VtkComposition, MultiPhaseBaseModel > InheritsFrom
Definition: flash/flashmodel.hh:66
Definition: multiphasebasemodel.hh:57
Definition: vtkcompositionmodule.hh:43
Definition: vtkdiffusionmodule.hh:46
Definition: vtkenergymodule.hh:43
Definition: vtkptflashmodule.hh:43