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>
63template <
class TypeTag>
68namespace Opm::Properties {
72struct FlashModel {
using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
76template<
class TypeTag>
77struct LocalResidual<TypeTag, TTag::FlashModel>
81template<
class TypeTag>
86template<
class TypeTag>
94template<
class TypeTag>
99template<
class TypeTag>
100struct PrimaryVariables<TypeTag, TTag::FlashModel>
104template<
class TypeTag>
105struct RateVector<TypeTag, TTag::FlashModel>
109template<
class TypeTag>
110struct BoundaryRateVector<TypeTag, TTag::FlashModel>
114template<
class TypeTag>
115struct IntensiveQuantities<TypeTag, TTag::FlashModel>
119template<
class TypeTag>
120struct ExtensiveQuantities<TypeTag, TTag::FlashModel>
124template<
class TypeTag>
125struct Indices<TypeTag, TTag::FlashModel>
129template<
class TypeTag>
130struct EnableDiffusion<TypeTag, TTag::FlashModel>
131{
static constexpr bool value =
false; };
134template<
class TypeTag>
135struct EnableEnergy<TypeTag, TTag::FlashModel>
136{
static constexpr bool value =
false; };
138template<
class TypeTag>
188template <
class TypeTag>
229 Parameters::Register<Parameters::FlashTolerance<Scalar>>
230 (
"The maximum tolerance for the flash solver to "
231 "consider the solution converged");
232 Parameters::Register<Parameters::FlashVerbosity>
233 (
"Flash solver verbosity level");
234 Parameters::Register<Parameters::FlashTwoPhaseMethod>
235 (
"Method for solving vapor-liquid composition. Available options include: "
236 "ssi, newton, ssi+newton");
238 Parameters::SetDefault<Parameters::FlashTolerance<Scalar>>(1.e-8);
239 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(
true);
244 Parameters::SetDefault<Parameters::EnableThermodynamicHints>(
true);
252 const std::string& tmp = EnergyModule::primaryVarName(
pvIdx);
256 std::ostringstream
oss;
257 if (Indices::z0Idx <=
pvIdx &&
pvIdx < Indices::z0Idx + numComponents - 1)
258 oss <<
"z_," << FluidSystem::componentName(
pvIdx - Indices::z0Idx);
259 else if (
pvIdx==Indices::pressure0Idx)
260 oss <<
"pressure_" << FluidSystem::phaseName(0);
272 const std::string& tmp = EnergyModule::eqName(
eqIdx);
276 std::ostringstream
oss;
277 if (Indices::conti0EqIdx <=
eqIdx &&
eqIdx < Indices::conti0EqIdx
280 oss <<
"continuity^" << FluidSystem::componentName(
compIdx);
288 void registerOutputModules_()
290 ParentType::registerOutputModules_();
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 flashindices.hh:47
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition flashintensivequantities.hh:58
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition flashlocalresidual.hh:46
A compositional multi-phase model based on flash-calculations.
Definition flashmodel.hh:191
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition flashmodel.hh:215
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition flashmodel.hh:250
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition flashmodel.hh:270
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 flashprimaryvariables.hh:51
Implements a vector representing rates of conserved quantities.
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
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:86
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:87
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:86
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:82
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...
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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
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.
Definition multiphasebaseproperties.hh:87
The type of the flash constraint solver.
Definition flashproperties.hh:39
Specifies the type of the actual Newton method.
Definition newtonmethodproperties.hh:32
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,...