flash/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 EWOMS_FLASH_MODEL_HH
29#define EWOMS_FLASH_MODEL_HH
30
31#include <opm/material/densead/Math.hpp>
32
33#include "flashproperties.hh"
35#include "flashlocalresidual.hh"
36#include "flashratevector.hh"
40#include "flashindices.hh"
41
47#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
48#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
49#include <opm/material/constraintsolvers/NcpFlash.hpp>
50
51#include <sstream>
52#include <string>
53
54namespace Opm {
55template <class TypeTag>
56class FlashModel;
57}
58
59namespace Opm::Properties {
60
61namespace TTag {
63struct FlashModel { using InheritsFrom = std::tuple<VtkDiffusion,
67} // namespace TTag
68
70template<class TypeTag>
72
74template<class TypeTag>
75struct FlashSolver<TypeTag, TTag::FlashModel>
76{ using type = Opm::NcpFlash<GetPropType<TypeTag, Properties::Scalar>,
78
80template<class TypeTag>
81struct FlashTolerance<TypeTag, TTag::FlashModel>
82{
84 static constexpr type value = -1.0;
85};
86
88template<class TypeTag>
89struct Model<TypeTag, TTag::FlashModel> { using type = Opm::FlashModel<TypeTag>; };
90
92template<class TypeTag>
94
96template<class TypeTag>
97struct RateVector<TypeTag, TTag::FlashModel> { using type = Opm::FlashRateVector<TypeTag>; };
98
100template<class TypeTag>
102
104template<class TypeTag>
106
108template<class TypeTag>
110
112template<class TypeTag>
113struct Indices<TypeTag, TTag::FlashModel> { using type = Opm::FlashIndices<TypeTag, /*PVIdx=*/0>; };
114
115// The updates of intensive quantities tend to be _very_ expensive for this
116// model, so let's try to minimize the number of required ones
117template<class TypeTag>
118struct EnableIntensiveQuantityCache<TypeTag, TTag::FlashModel> { static constexpr bool value = true; };
119
120// since thermodynamic hints are basically free if the cache for intensive quantities is
121// enabled, and this model usually shows quite a performance improvment if they are
122// enabled, let's enable them by default.
123template<class TypeTag>
124struct EnableThermodynamicHints<TypeTag, TTag::FlashModel> { static constexpr bool value = true; };
125
126// disable molecular diffusion by default
127template<class TypeTag>
128struct EnableDiffusion<TypeTag, TTag::FlashModel> { static constexpr bool value = false; };
129
131template<class TypeTag>
132struct EnableEnergy<TypeTag, TTag::FlashModel> { static constexpr bool value = false; };
133
134} // namespace Opm::Properties
135
136namespace Opm {
137
197template <class TypeTag>
198class FlashModel
199 : public MultiPhaseBaseModel<TypeTag>
200{
201 using ParentType = MultiPhaseBaseModel<TypeTag>;
202
206
208
209 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
210 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
211 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
212
213
214 using EnergyModule = Opm::EnergyModule<TypeTag, enableEnergy>;
215
216public:
217 FlashModel(Simulator& simulator)
218 : ParentType(simulator)
219 {}
220
224 static void registerParameters()
225 {
227
228 // register runtime parameters of the VTK output modules
230
231 if (enableDiffusion)
233
234 if (enableEnergy)
236
237 Parameters::registerParam<TypeTag, Properties::FlashTolerance>
238 ("The maximum tolerance for the flash solver to "
239 "consider the solution converged");
240 }
241
245 static std::string name()
246 { return "flash"; }
247
251 std::string primaryVarName(unsigned pvIdx) const
252 {
253 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
254 if (tmp != "")
255 return tmp;
256
257 std::ostringstream oss;
258 if (Indices::cTot0Idx <= pvIdx && pvIdx < Indices::cTot0Idx
259 + numComponents)
260 oss << "c_tot," << FluidSystem::componentName(/*compIdx=*/pvIdx
261 - Indices::cTot0Idx);
262 else
263 assert(false);
264
265 return oss.str();
266 }
267
271 std::string eqName(unsigned eqIdx) const
272 {
273 const std::string& tmp = EnergyModule::eqName(eqIdx);
274 if (tmp != "")
275 return tmp;
276
277 std::ostringstream oss;
278 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
279 + numComponents) {
280 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
281 oss << "continuity^" << FluidSystem::componentName(compIdx);
282 }
283 else
284 assert(false);
285
286 return oss.str();
287 }
288
292 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
293 {
294 Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
295 if (tmp > 0)
296 return tmp;
297
298 unsigned compIdx = pvIdx - Indices::cTot0Idx;
299
300 // make all kg equal. also, divide the weight of all total
301 // compositions by 100 to make the relative errors more
302 // comparable to the ones of the other models (at 10% porosity
303 // the medium is fully saturated with water at atmospheric
304 // conditions if 100 kg/m^3 are present!)
305 return FluidSystem::molarMass(compIdx) / 100.0;
306 }
307
311 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
312 {
313 Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
314 if (tmp > 0)
315 return tmp;
316
317 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
318
319 // make all kg equal
320 return FluidSystem::molarMass(compIdx);
321 }
322
324 {
326
327 // add the VTK output modules which are meaningful for the model
328 this->addOutputModule(new Opm::VtkCompositionModule<TypeTag>(this->simulator_));
329 if (enableDiffusion)
330 this->addOutputModule(new Opm::VtkDiffusionModule<TypeTag>(this->simulator_));
331 if (enableEnergy)
332 this->addOutputModule(new Opm::VtkEnergyModule<TypeTag>(this->simulator_));
333 }
334};
335
336} // namespace Opm
337
338#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: flash/flashmodel.hh:217
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: flash/flashmodel.hh:224
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: flash/flashmodel.hh:292
static std::string name()
Definition: flash/flashmodel.hh:245
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: flash/flashmodel.hh:311
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: flash/flashmodel.hh:251
void registerOutputModules_()
Definition: flash/flashmodel.hh:323
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: flash/flashmodel.hh:271
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
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:102
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
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
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:82
Specify whether energy should be considered as a conservation quantity or not.
Definition: multiphasebaseproperties.hh:76
Specify whether all intensive quantities for the grid should be cached in the discretization.
Definition: fvbaseproperties.hh:297
Specify whether to use the already calculated solutions as starting values of the intensive quantitie...
Definition: fvbaseproperties.hh:317
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
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
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