28#ifndef OPM_PTFLASH_MODEL_HH
29#define OPM_PTFLASH_MODEL_HH
31#include <opm/material/constraintsolvers/PTFlash.hpp>
33#include <opm/material/densead/Math.hpp>
35#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
36#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
66template <
class TypeTag>
71namespace Opm::Properties {
75struct FlashModel {
using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
79template<
class TypeTag>
81{
using type = FlashLocalResidual<TypeTag>; };
84template<
class TypeTag>
89template<
class TypeTag>
92 using type = PTFlash<GetPropType<TypeTag, Properties::Scalar>,
97template<
class TypeTag>
102template<
class TypeTag>
103struct PrimaryVariables<TypeTag, TTag::FlashModel>
104{
using type = FlashPrimaryVariables<TypeTag>; };
107template<
class TypeTag>
109{
using type = FlashRateVector<TypeTag>; };
112template<
class TypeTag>
114{
using type = FlashBoundaryRateVector<TypeTag>; };
117template<
class TypeTag>
119{
using type = FlashIntensiveQuantities<TypeTag>; };
122template<
class TypeTag>
124{
using type = FlashExtensiveQuantities<TypeTag>; };
127template<
class TypeTag>
129{
using type = FlashIndices<TypeTag, 0>; };
132template<
class TypeTag>
134{
static constexpr bool value =
false; };
137template<
class TypeTag>
139{
static constexpr bool value =
false; };
141template<
class TypeTag>
191template <
class TypeTag>
193 :
public MultiPhaseBaseModel<TypeTag>
195 using ParentType = MultiPhaseBaseModel<TypeTag>;
211 explicit FlashModel(Simulator& simulator)
212 : ParentType(simulator)
226 if constexpr (enableDiffusion) {
230 if constexpr (enableEnergy) {
235 (
"The maximum tolerance for the flash solver to "
236 "consider the solution converged");
238 (
"Flash solver verbosity level");
240 (
"Method for solving vapor-liquid composition. Available options include: "
241 "ssi, newton, ssi+newton");
257 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
262 std::ostringstream oss;
263 if (Indices::z0Idx <= pvIdx && pvIdx < Indices::z0Idx + numComponents - 1) {
264 oss <<
"z_," << FluidSystem::componentName(pvIdx - Indices::z0Idx);
266 else if (pvIdx == Indices::pressure0Idx) {
267 oss <<
"pressure_" << FluidSystem::phaseName(0);
281 const std::string& tmp = EnergyModule::eqName(eqIdx);
286 std::ostringstream oss;
287 if (Indices::conti0EqIdx <= eqIdx &&
288 eqIdx < Indices::conti0EqIdx + numComponents)
290 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
291 oss <<
"continuity^" << FluidSystem::componentName(compIdx);
300 void registerOutputModules_()
302 ParentType::registerOutputModules_();
307 if constexpr (enableDiffusion) {
310 if constexpr (enableEnergy) {
311 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:54
A compositional multi-phase model based on flash-calculations.
Definition flashmodel.hh:194
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition flashmodel.hh:218
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition flashmodel.hh:255
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition flashmodel.hh:279
A Newton solver specific to the PTFlash model.
Definition flashnewtonmethod.hh:55
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition multiphasebasemodel.hh:175
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition multiphasebasemodel.hh:197
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
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.
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
This template class contains the data which is required to calculate all fluxes of components over a ...
A Newton solver specific to the PTFlash model.
Declares the properties required by the compositional multi-phase model based on flash calculations.
Implements a vector representing rates of conserved quantities.
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:82
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void SetDefault(decltype(Param::value) new_value)
Set a runtime parameter.
Definition parametersystem.hpp:212
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
Contains the intensive quantities of the flash-based compositional multi-phase model.
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Declares the parameters for the compositional multi-phase model based on flash calculations.
Represents the primary variables used by the compositional flow model based on flash calculations.
Type of object for specifying boundary conditions.
Definition fvbaseproperties.hh:119
Enable diffusive fluxes?
Definition multiphasebaseproperties.hh:91
Specify whether energy should be considered as a conservation quantity or not.
Definition multiphasebaseproperties.hh:87
Enable water?
Definition multiphasebaseproperties.hh:103
Data required to calculate a flux over a face.
Definition fvbaseproperties.hh:153
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
Specifies the type of the actual Newton method.
Definition newtonmethodproperties.hh:32
Vector containing volumetric or areal rates of quantities.
Definition fvbaseproperties.hh:116
The type tag for the isothermal single phase problems.
Definition flashmodel.hh:74
VTK output module for the fluid composition.
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
VTK output module for quantities which make sense for models which assume thermal equilibrium.
VTK output module for the PT Flash calculation This module deals with the following quantities: K,...