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  Copyright (C) 2011-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
26 #ifndef EWOMS_FLASH_MODEL_HH
27 #define EWOMS_FLASH_MODEL_HH
28 
29 #include <opm/material/localad/Math.hpp>
30 
31 #include "flashproperties.hh"
32 #include "flashprimaryvariables.hh"
33 #include "flashlocalresidual.hh"
34 #include "flashratevector.hh"
38 #include "flashindices.hh"
39 
45 #include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
46 #include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
47 #include <opm/material/constraintsolvers/NcpFlash.hpp>
48 
49 #include <sstream>
50 #include <string>
51 
52 namespace Ewoms {
53 template <class TypeTag>
54 class FlashModel;
55 }
56 
57 namespace Ewoms {
58 namespace Properties {
61  VtkComposition,
62  VtkEnergy,
63  VtkDiffusion));
64 
66 SET_TYPE_PROP(FlashModel, LocalResidual,
68 
70 SET_TYPE_PROP(FlashModel, FlashSolver,
71  Opm::NcpFlash<typename GET_PROP_TYPE(TypeTag, Scalar),
72  typename GET_PROP_TYPE(TypeTag, FluidSystem)>);
73 
75 SET_SCALAR_PROP(FlashModel, FlashTolerance, 0.0);
76 
79 
82 
85 
88 
91 
94 
97 
98 // The updates of intensive quantities tend to be _very_ expensive for this
99 // model, so let's try to minimize the number of required ones
100 SET_BOOL_PROP(FlashModel, EnableIntensiveQuantityCache, true);
101 
102 // since thermodynamic hints are basically free if the cache for intensive quantities is
103 // enabled, and this model usually shows quite a performance improvment if they are
104 // enabled, let's enable them by default.
105 SET_BOOL_PROP(FlashModel, EnableThermodynamicHints, true);
106 
107 // disable molecular diffusion by default
108 SET_BOOL_PROP(FlashModel, EnableDiffusion, false);
109 
111 SET_BOOL_PROP(FlashModel, EnableEnergy, false);
112 } // namespace Properties
113 
172 template <class TypeTag>
173 class FlashModel
174  : public MultiPhaseBaseModel<TypeTag>
175 {
176  typedef MultiPhaseBaseModel<TypeTag> ParentType;
177 
178  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
179  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
180  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
181 
182  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
183 
184  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
185  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
186  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
187 
188 
189  typedef Ewoms::EnergyModule<TypeTag, enableEnergy> EnergyModule;
190 
191 public:
192  FlashModel(Simulator &simulator)
193  : ParentType(simulator)
194  {}
195 
199  static void registerParameters()
200  {
202 
203  // register runtime parameters of the VTK output modules
205 
206  if (enableDiffusion)
208 
209  if (enableEnergy)
211 
212  EWOMS_REGISTER_PARAM(TypeTag, Scalar, FlashTolerance,
213  "The maximum tolerance for the flash solver to "
214  "consider the solution converged");
215  }
216 
220  static std::string name()
221  { return "flash"; }
222 
226  std::string primaryVarName(int pvIdx) const
227  {
228  const std::string &tmp = EnergyModule::primaryVarName(pvIdx);
229  if (tmp != "")
230  return tmp;
231 
232  std::ostringstream oss;
233  if (Indices::cTot0Idx <= pvIdx && pvIdx < Indices::cTot0Idx
234  + numComponents)
235  oss << "c_tot," << FluidSystem::componentName(/*compIdx=*/pvIdx
236  - Indices::cTot0Idx);
237  else
238  assert(false);
239 
240  return oss.str();
241  }
242 
246  std::string eqName(int eqIdx) const
247  {
248  const std::string &tmp = EnergyModule::eqName(eqIdx);
249  if (tmp != "")
250  return tmp;
251 
252  std::ostringstream oss;
253  if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
254  + numComponents) {
255  int compIdx = eqIdx - Indices::conti0EqIdx;
256  oss << "continuity^" << FluidSystem::componentName(compIdx);
257  }
258  else
259  assert(false);
260 
261  return oss.str();
262  }
263 
267  Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
268  {
269  Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
270  if (tmp > 0)
271  return tmp;
272 
273  int compIdx = pvIdx - Indices::cTot0Idx;
274 
275  // make all kg equal. also, divide the weight of all total
276  // compositions by 100 to make the relative errors more
277  // comparable to the ones of the other models (at 10% porosity
278  // the medium is fully saturated with water at atmospheric
279  // conditions if 100 kg/m^3 are present!)
280  return FluidSystem::molarMass(compIdx) / 100.0;
281  }
282 
286  Scalar eqWeight(int globalDofIdx, int eqIdx) const
287  {
288  Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
289  if (tmp > 0)
290  return tmp;
291 
292  int compIdx = eqIdx - Indices::conti0EqIdx;
293 
294  // make all kg equal
295  return FluidSystem::molarMass(compIdx);
296  }
297 
299  {
301 
302  // add the VTK output modules which are meaningful for the model
303  this->addOutputModule(new Ewoms::VtkCompositionModule<TypeTag>(this->simulator_));
304  if (enableDiffusion)
305  this->addOutputModule(new Ewoms::VtkDiffusionModule<TypeTag>(this->simulator_));
306  if (enableEnergy)
307  this->addOutputModule(new Ewoms::VtkEnergyModule<TypeTag>(this->simulator_));
308  }
309 };
310 
311 } // namespace Ewoms
312 
313 #endif
This template class contains the data which is required to calculate all fluxes of components over a ...
Definition: flashextensivequantities.hh:50
void registerOutputModules_()
Definition: multiphasebasemodel.hh:244
VTK output module for quantities which make sense for models which assume thermal equilibrium...
Implements a vector representing rates of conserved quantities.
Definition: flashratevector.hh:45
std::string eqName(int eqIdx) const
Given an equation index, return a human readable name.
Definition: flashmodel.hh:246
Implements a vector representing rates of conserved quantities.
SET_SCALAR_PROP(NumericModel, EndTime,-1e100)
The default value for the simulation's end time.
A compositional multi-phase model based on flash-calculations.
Definition: flashmodel.hh:54
VTK output module for quantities which make sense for models which assume thermal equilibrium...
Definition: vtkenergymodule.hh:70
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
Definition: flashindices.hh:43
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
Contains the intensive quantities of the flash-based compositional multi-phase model.
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:45
This template class contains the data which is required to calculate all fluxes of components over a ...
Calculates the local residual of the compositional multi-phase model based on flash calculations...
Definition: flashlocalresidual.hh:42
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flashintensivequantities.hh:48
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:44
FlashModel(Simulator &simulator)
Definition: flashmodel.hh:192
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:93
Calculates the local residual of the compositional multi-phase model based on flash calculations...
SET_TYPE_PROP(NumericModel, Scalar, double)
Set the default type of scalar values to double.
std::string primaryVarName(int pvIdx) const
Given an primary variable index, return a human readable name.
Definition: flashmodel.hh:226
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: flashmodel.hh:267
Definition: baseauxiliarymodule.hh:35
VTK output module for the fluid composition.
VTK output module for quantities which make sense for models which incorperate molecular diffusion...
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:145
#define EWOMS_REGISTER_PARAM(TypeTag, ParamType, ParamName, Description)
Register a run-time parameter.
Definition: parametersystem.hh:64
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:102
static std::string name()
Definition: flashmodel.hh:220
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
void registerOutputModules_()
Definition: flashmodel.hh:298
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
NEW_TYPE_TAG(AuxModule)
Declares the properties required by the compositional multi-phase model based on flash calculations...
Represents the primary variables used by the compositional flow model based on flash calculations...
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: flashmodel.hh:199
SET_BOOL_PROP(FvBaseDiscretization, EnableVtkOutput, true)
Enable the VTK output by default.
Scalar eqWeight(int globalDofIdx, int eqIdx) const
Returns the relative weight of an equation.
Definition: flashmodel.hh:286
VTK output module for quantities which make sense for models which incorperate molecular diffusion...
Definition: vtkdiffusionmodule.hh:67
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
#define INHERITS_FROM(...)
Syntactic sugar for NEW_TYPE_TAG.
Definition: propertysystem.hh:229
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:98
Represents the primary variables used by the compositional flow model based on flash calculations...
Definition: flashprimaryvariables.hh:54
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:74