blackoilmodel.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) 2010-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_BLACK_OIL_MODEL_HH
27 #define EWOMS_BLACK_OIL_MODEL_HH
28 
29 #include <opm/material/localad/Math.hpp>
30 
31 #include "blackoilproblem.hh"
32 #include "blackoilindices.hh"
36 #include "blackoilratevector.hh"
38 #include "blackoillocalresidual.hh"
39 #include "blackoilnewtonmethod.hh"
40 #include "blackoilproperties.hh"
41 
45 
46 #include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
47 
48 #include <sstream>
49 #include <string>
50 
51 namespace Ewoms {
52 template <class TypeTag>
54 
55 template <class TypeTag>
57 }
58 
59 namespace Ewoms {
60 namespace Properties {
63  VtkBlackOil,
64  VtkComposition));
65 
67 SET_TYPE_PROP(BlackOilModel, LocalResidual,
69 
72 
75 
78 
80 SET_PROP(BlackOilModel, BlackOilFluidState)
81 { private:
82  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
83  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
84 public:
85  typedef Opm::CompositionalFluidState<Evaluation,
86  FluidSystem,
87  /*enableEnthalpy=*/false> type;
88 };
89 
92 
95 
98 
101 
104 
107 
109 SET_PROP(BlackOilModel, FluidSystem)
110 {
111 private:
112  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
113  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
114 
115 public:
116  typedef Opm::FluidSystems::BlackOil<Scalar> type;
117 };
118 
121 SET_INT_PROP(BlackOilModel, BlackoilNumChoppedIterations, 4);
122 } // namespace Properties
123 
186 template<class TypeTag >
187 class BlackOilModel
188  : public MultiPhaseBaseModel<TypeTag>
189 {
190  typedef typename GET_PROP_TYPE(TypeTag, Model) Implementation;
191  typedef MultiPhaseBaseModel<TypeTag> ParentType;
192 
193  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
194  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
195  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
196  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
197  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
198  typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
199 
200  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
201  enum { numComponents = FluidSystem::numComponents };
202  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
203 
204 public:
206  : ParentType(simulator)
207  {}
208 
212  static void registerParameters()
213  {
215 
216  // register runtime parameters of the VTK output modules
219  }
220 
224  void finishInit()
225  {
227 
228  Dune::FMatrixPrecision<Scalar>::set_singular_limit(1e-35);
229  }
230 
234  static std::string name()
235  { return "blackoil"; }
236 
240  std::string primaryVarName(int pvIdx) const
241  {
242  std::ostringstream oss;
243 
244  if (pvIdx == Indices::gasPressureIdx)
245  oss << "pressure_" << FluidSystem::phaseName(FluidSystem::gasPhaseIdx);
246  else if (pvIdx == Indices::waterSaturationIdx)
247  oss << "saturation_" << FluidSystem::phaseName(FluidSystem::waterPhaseIdx);
248  else if (pvIdx == Indices::switchIdx)
249  oss << "switching";
250  else
251  assert(false);
252 
253  return oss.str();
254  }
255 
259  std::string eqName(int eqIdx) const
260  {
261  std::ostringstream oss;
262 
263  if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx + numComponents)
264  oss << "conti_" << FluidSystem::phaseName(eqIdx - Indices::conti0EqIdx);
265  else
266  assert(false);
267 
268  return oss.str();
269  }
270 
274  Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
275  {
276  // do not care about the auxiliary equations as they are supposed to scale
277  // themselves
278  if (globalDofIdx >= (int) this->numGridDof())
279  return 1;
280 
281  if (Indices::gasPressureIdx == pvIdx)
282  return 10/referencePressure_;
283 
284  return 1;
285  }
286 
290  Scalar eqWeight(int globalDofIdx, int eqIdx) const
291  {
292  // do not care about the auxiliary equations as they are supposed to scale
293  // themselves
294  if (globalDofIdx >= (int) this->numGridDof())
295  return 1;
296 
297  int compIdx = eqIdx - Indices::conti0EqIdx;
298  assert(0 <= compIdx && compIdx <= numPhases);
299 
300  // make all kg equal
301  return FluidSystem::molarMass(compIdx);
302  }
303 
307  void updateBegin()
308  {
309  ParentType::updateBegin();
310 
311  // find the a reference pressure. The first degree of freedom
312  // might correspond to non-interior entities which would lead
313  // to an undefined value, so we have to iterate...
314  for (size_t dofIdx = 0; dofIdx < this->numGridDof(); ++ dofIdx) {
315  if (this->isLocalDof(dofIdx)) {
316  referencePressure_ = this->solution(/*timeIdx=*/0)[dofIdx][Indices::gasPressureIdx];
317  break;
318  }
319  }
320  }
321 
330  {
331  numSwitched_ = 0;
332 
333  int numDof = this->numGridDof();
334  for (int globalDofIdx = 0; globalDofIdx < numDof; ++globalDofIdx) {
335  auto &priVars = this->solution(/*timeIdx=*/0)[globalDofIdx];
336  if (priVars.adaptSwitchingVariable()) {
337  this->linearizer().markDofRed(globalDofIdx);
338  ++numSwitched_;
339  }
340  }
341 
342  // make sure that if there was a variable switch in an
343  // other partition we will also set the switch flag
344  // for our partition.
345  numSwitched_ = this->gridView_.comm().sum(numSwitched_);
346 
347  this->simulator_.model().newtonMethod().endIterMsg()
348  << ", num switched=" << numSwitched_;
349  }
350 
359  template <class DofEntity>
360  void serializeEntity(std::ostream &outstream, const DofEntity &dof)
361  {
362 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
363  int dofIdx = asImp_().dofMapper().index(dof);
364 #else
365  int dofIdx = asImp_().dofMapper().map(dof);
366 #endif
367 
368  // write phase state
369  if (!outstream.good()) {
370  OPM_THROW(std::runtime_error, "Could not serialize degree of freedom " << dofIdx);
371  }
372 
373  // write the primary variables
374  const auto& priVars = this->solution_[/*timeIdx=*/0][dofIdx];
375  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
376  outstream << priVars[eqIdx] << " ";
377 
378  // write the pseudo primary variables
379  outstream << priVars.switchingVarMeaning() << " ";
380  outstream << priVars.pvtRegionIndex() << " ";
381  }
382 
391  template <class DofEntity>
392  void deserializeEntity(std::istream &instream,
393  const DofEntity &dof)
394  {
395 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
396  int dofIdx = asImp_().dofMapper().index(dof);
397 #else
398  int dofIdx = asImp_().dofMapper().map(dof);
399 #endif
400 
401  // read in the "real" primary variables of the DOF
402  auto& priVars = this->solution_[/*timeIdx=*/0][dofIdx];
403  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
404  if (!instream.good())
405  OPM_THROW(std::runtime_error,
406  "Could not deserialize degree of freedom " << dofIdx);
407  instream >> priVars[eqIdx];
408  }
409 
410  // read the pseudo primary variables
411  int switchingVarMeaning;
412  instream >> switchingVarMeaning;
413 
414  int pvtRegionIdx;
415  instream >> pvtRegionIdx;
416 
417  if (!instream.good())
418  OPM_THROW(std::runtime_error,
419  "Could not deserialize degree of freedom " << dofIdx);
420 
421  typedef typename PrimaryVariables::SwitchingVarMeaning SVM;
422  priVars.setSwitchingVarMeaning(static_cast<SVM>(switchingVarMeaning));
423  priVars.setPvtRegionIndex(pvtRegionIdx);
424  }
425 
433  template <class Restarter>
434  void deserialize(Restarter &res)
435  {
436  ParentType::deserialize(res);
437 
438  // set the PVT indices of the primary variables. This is also done by writing
439  // them into the restart file and re-reading them, but it is better to calculate
440  // them from scratch because the input could have been changed in this regard...
441  ElementContext elemCtx(this->simulator_);
442  auto elemIt = this->gridView().template begin</*codim=*/0>();
443  auto elemEndIt = this->gridView().template end</*codim=*/0>();
444  for (; elemIt != elemEndIt; ++ elemIt) {
445  elemCtx.updateStencil(*elemIt);
446  for (int dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timIdx=*/0); ++dofIdx) {
447  int globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timIdx=*/0);
448  updatePvtRegionIndex_(this->solution_[/*timeIdx=*/0][globalDofIdx],
449  elemCtx, dofIdx, /*timeIdx=*/0);
450  }
451  }
452 
453  this->solution_[/*timeIdx=*/1] = this->solution_[/*timeIdx=*/0];
454  }
455 
456 
457 // HACK: this should be made private and the BaseModel should be
458 // declared to be a friend. Since C++-2003 (and more relevantly GCC
459 // 4.4) don't support friend typedefs, we need to make this method
460 // public until the oldest supported compiler supports friend
461 // typedefs...
462 
463 //protected:
464 // friend typename GET_PROP_TYPE(TypeTag, Discretization);
465 
466  template <class Context>
467  void supplementInitialSolution_(PrimaryVariables &priVars,
468  const Context &context, int dofIdx, int timeIdx)
469  { updatePvtRegionIndex_(priVars, context, dofIdx, timeIdx); }
470 
472  {
474 
475  // add the VTK output modules which make sense for the blackoil model
476  this->addOutputModule(new Ewoms::VtkBlackOilModule<TypeTag>(this->simulator_));
477  this->addOutputModule(new Ewoms::VtkCompositionModule<TypeTag>(this->simulator_));
478  }
479 
480 private:
481  Implementation &asImp_()
482  { return *static_cast<Implementation*>(this); }
483  const Implementation &asImp_() const
484  { return *static_cast<const Implementation*>(this); }
485 
486  template <class Context>
487  void updatePvtRegionIndex_(PrimaryVariables &priVars, const Context &context, int dofIdx, int timeIdx)
488  {
489  int regionIdx = context.problem().pvtRegionIndex(context, dofIdx, timeIdx);
490  priVars.setPvtRegionIndex(regionIdx);
491  }
492 
493  mutable Scalar referencePressure_;
494  int numSwitched_;
495 };
496 } // namespace Ewoms
497 
498 #endif
void registerOutputModules_()
Definition: multiphasebasemodel.hh:244
void finishInit()
Apply the initial conditions to the model.
Definition: blackoilmodel.hh:224
This template class contains the data which is required to calculate the fluxes of the fluid phases o...
Definition: blackoilextensivequantities.hh:47
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:38
static std::string name()
Definition: blackoilmodel.hh:234
Calculates the local residual of the black oil model.
Scalar eqWeight(int globalDofIdx, int eqIdx) const
Returns the relative weight of an equation.
Definition: blackoilmodel.hh:290
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmodule.hh:101
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: blackoilmodel.hh:307
void supplementInitialSolution_(PrimaryVariables &priVars, const Context &context, int dofIdx, int timeIdx)
Definition: blackoilmodel.hh:467
The multi-dimensional Newton method.
Definition: newtonmethod.hh:54
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
Implements a vector representing mass, molar or volumetric rates for the black oil model...
Implements a boundary vector for the fully implicit black-oil model.
BlackOilModel(Simulator &simulator)
Definition: blackoilmodel.hh:205
The primary variable and equation indices for the black-oil model.
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:45
Base class for all problems which use the black-oil model.
Definition: blackoilproblem.hh:40
std::string primaryVarName(int pvIdx) const
Given an primary variable index, return a human readable name.
Definition: blackoilmodel.hh:240
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
Implements a vector representing mass, molar or volumetric rates for the black oil model...
Definition: blackoilratevector.hh:48
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:45
void registerOutputModules_()
Definition: blackoilmodel.hh:471
Definition: blackoilmodel.hh:56
std::string eqName(int eqIdx) const
Given an equation index, return a human readable name.
Definition: blackoilmodel.hh:259
SET_PROP(NumericModel, ParameterTree)
Set the ParameterTree property.
Definition: basicproperties.hh:117
SET_INT_PROP(NumericModel, GridGlobalRefinements, 0)
void deserializeEntity(std::istream &instream, const DofEntity &dof)
Reads the current solution variables for a degree of freedom from a restart file. ...
Definition: blackoilmodel.hh:392
void deserialize(Restarter &res)
Deserializes the state of the model.
Definition: blackoilmodel.hh:434
Declares the properties required by the black oil model.
Implements a boundary vector for the fully implicit black-oil model.
Definition: blackoilboundaryratevector.hh:42
MultiPhaseBaseModel(Simulator &simulator)
Definition: multiphasebasemodel.hh:138
The primary variable and equation indices for the black-oil model.
Definition: blackoilindices.hh:37
SET_TYPE_PROP(NumericModel, Scalar, double)
Set the default type of scalar values to double.
Contains the quantities which are are constant within a finite volume in the black-oil model...
Definition: blackoilintensivequantities.hh:43
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: blackoilmodel.hh:212
Definition: baseauxiliarymodule.hh:35
void finishInit()
Apply the initial conditions to the model.
Definition: multiphasebasemodel.hh:157
VTK output module for the fluid composition.
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:145
Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: blackoilmodel.hh:274
VTK output module for the black oil model's parameters.
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:102
void switchPrimaryVars_()
Definition: blackoilmodel.hh:329
This template class contains the data which is required to calculate the fluxes of the fluid phases o...
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Represents the primary variables used by the black-oil model.
Contains the quantities which are are constant within a finite volume in the black-oil model...
A newton solver which is specific to the black oil model.
Definition: blackoilnewtonmethod.hh:39
NEW_TYPE_TAG(AuxModule)
Base class for all problems which use the black-oil model.
VTK output module for the black oil model's parameters.
Definition: vtkblackoilmodule.hh:70
#define INHERITS_FROM(...)
Syntactic sugar for NEW_TYPE_TAG.
Definition: propertysystem.hh:229
void serializeEntity(std::ostream &outstream, const DofEntity &dof)
Write the current solution for a degree of freedom to a restart file.
Definition: blackoilmodel.hh:360
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:74
A newton solver which is specific to the black oil model.
A fully-implicit black-oil flow model.
Definition: blackoilmodel.hh:53