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/constraintsolvers/NcpFlash.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
50
54
55#include <cassert>
56#include <memory>
57#include <sstream>
58#include <string>
59#include <tuple>
60
61namespace Opm {
62
63template <class TypeTag>
64class FlashModel;
65
66}
67
68namespace Opm::Properties {
69
70namespace TTag {
71
74{ using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
75
76} // namespace TTag
77
79template<class TypeTag>
80struct LocalResidual<TypeTag, TTag::FlashModel>
82
84template<class TypeTag>
85struct FlashSolver<TypeTag, TTag::FlashModel>
86{
87 using type = NcpFlash<GetPropType<TypeTag, Properties::Scalar>,
89};
90
92template<class TypeTag>
93struct Model<TypeTag, TTag::FlashModel>
95
97template<class TypeTag>
98struct PrimaryVariables<TypeTag, TTag::FlashModel>
100
102template<class TypeTag>
103struct RateVector<TypeTag, TTag::FlashModel>
105
107template<class TypeTag>
108struct BoundaryRateVector<TypeTag, TTag::FlashModel>
110
112template<class TypeTag>
113struct IntensiveQuantities<TypeTag, TTag::FlashModel>
115
117template<class TypeTag>
118struct ExtensiveQuantities<TypeTag, TTag::FlashModel>
120
122template<class TypeTag>
123struct Indices<TypeTag, TTag::FlashModel>
124{ using type = FlashIndices<TypeTag, /*PVIdx=*/0>; };
125
126// disable molecular diffusion by default
127template<class TypeTag>
128struct EnableDiffusion<TypeTag, TTag::FlashModel>
129{ static constexpr bool value = false; };
130
132template<class TypeTag>
133struct EnableEnergy<TypeTag, TTag::FlashModel>
134{ static constexpr bool value = false; };
135
136} // namespace Opm::Properties
137
138namespace Opm {
139
199template <class TypeTag>
200class FlashModel
201 : public MultiPhaseBaseModel<TypeTag>
202{
203 using ParentType = MultiPhaseBaseModel<TypeTag>;
204
208
210
211 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
212 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
213 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
214
216
217public:
218 explicit FlashModel(Simulator& simulator)
219 : ParentType(simulator)
220 {}
221
225 static void registerParameters()
226 {
228
229 // register runtime parameters of the VTK output modules
231
232 if constexpr (enableDiffusion) {
234 }
235
236 if constexpr (enableEnergy) {
238 }
239
240 Parameters::Register<Parameters::FlashTolerance<Scalar>>
241 ("The maximum tolerance for the flash solver to "
242 "consider the solution converged");
243
244 // The updates of intensive quantities tend to be _very_ expensive for this
245 // model, so let's try to minimize the number of required ones
246 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(true);
247
248 // since thermodynamic hints are basically free if the cache for intensive quantities is
249 // enabled, and this model usually shows quite a performance improvment if they are
250 // enabled, let's enable them by default.
251 Parameters::SetDefault<Parameters::EnableThermodynamicHints>(true);
252 }
253
257 static std::string name()
258 { return "flash"; }
259
263 std::string primaryVarName(unsigned pvIdx) const
264 {
265 const std::string& tmp = EnergyModule::primaryVarName(pvIdx);
266 if (!tmp.empty()) {
267 return tmp;
268 }
269
270 if (Indices::cTot0Idx <= pvIdx && pvIdx < Indices::cTot0Idx + numComponents) {
271 std::ostringstream oss;
272 oss << "c_tot," << FluidSystem::componentName(/*compIdx=*/pvIdx - Indices::cTot0Idx);
273 return oss.str();
274 }
275 else {
276 assert(false);
277 return "";
278 }
279 }
280
284 std::string eqName(unsigned eqIdx) const
285 {
286 const std::string& tmp = EnergyModule::eqName(eqIdx);
287 if (!tmp.empty()) {
288 return tmp;
289 }
290
291 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx + numComponents) {
292 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
293 std::ostringstream oss;
294 oss << "continuity^" << FluidSystem::componentName(compIdx);
295 return oss.str();
296 }
297 else {
298 assert(false);
299 return "";
300 }
301 }
302
306 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
307 {
308 const Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
309 if (tmp > 0) {
310 return tmp;
311 }
312
313 const unsigned compIdx = pvIdx - Indices::cTot0Idx;
314
315 // make all kg equal. also, divide the weight of all total
316 // compositions by 100 to make the relative errors more
317 // comparable to the ones of the other models (at 10% porosity
318 // the medium is fully saturated with water at atmospheric
319 // conditions if 100 kg/m^3 are present!)
320 return FluidSystem::molarMass(compIdx) / 100.0;
321 }
322
326 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
327 {
328 const Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
329 if (tmp > 0) {
330 return tmp;
331 }
332
333 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
334
335 // make all kg equal
336 return FluidSystem::molarMass(compIdx);
337 }
338
340 {
342
343 // add the VTK output modules which are meaningful for the model
344 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
345 if constexpr (enableDiffusion) {
346 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
347 }
348 if constexpr (enableEnergy) {
349 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
350 }
351 }
352};
353
354} // namespace Opm
355
356#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: flash/flashmodel.hh:218
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: flash/flashmodel.hh:225
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: flash/flashmodel.hh:306
static std::string name()
Definition: flash/flashmodel.hh:257
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: flash/flashmodel.hh:326
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: flash/flashmodel.hh:263
void registerOutputModules_()
Definition: flash/flashmodel.hh:339
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: flash/flashmodel.hh:284
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
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:84
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
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declares the parameters for the compositional multi-phase model based on flash calculations.
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
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
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
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