pvsmodel.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) 2009-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_PVS_MODEL_HH
27 #define EWOMS_PVS_MODEL_HH
28 
29 #include <opm/material/localad/Math.hpp>
30 
31 #include "pvsproperties.hh"
32 #include "pvslocalresidual.hh"
33 #include "pvsnewtonmethod.hh"
34 #include "pvsprimaryvariables.hh"
35 #include "pvsratevector.hh"
36 #include "pvsboundaryratevector.hh"
39 #include "pvsindices.hh"
40 
47 #include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
48 #include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
49 
50 #include <iostream>
51 #include <sstream>
52 #include <string>
53 #include <vector>
54 
55 namespace Ewoms {
56 template <class TypeTag>
57 class PvsModel;
58 }
59 
60 namespace Ewoms {
61 namespace Properties {
64  VtkPhasePresence,
65  VtkComposition,
66  VtkEnergy,
67  VtkDiffusion));
68 
71  LocalResidual,
73 
76 
79 
82 
85 
88 
91 
94 
97 
98 // set the model to a medium verbosity
99 SET_INT_PROP(PvsModel, PvsVerbosity, 1);
100 
102 SET_BOOL_PROP(PvsModel, EnableEnergy, false);
103 
104 // disable molecular diffusion by default
105 SET_BOOL_PROP(PvsModel, EnableDiffusion, false);
106 
108 SET_SCALAR_PROP(PvsModel, PvsPressureBaseWeight, 1.0);
109 
111 SET_SCALAR_PROP(PvsModel, PvsSaturationsBaseWeight, 1.0);
112 
114 SET_SCALAR_PROP(PvsModel, PvsMoleFractionsBaseWeight, 1.0);
115 } // namespace Properties
116 
212 template <class TypeTag>
213 class PvsModel
214  : public MultiPhaseBaseModel<TypeTag>
215 {
216  typedef MultiPhaseBaseModel<TypeTag> ParentType;
217 
218  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
219  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
220  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
221  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
222 
223  typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
224  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
225  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
226  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
227 
228  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
229  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
230  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
231  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
232 
233  typedef typename GridView::template Codim<0>::Entity Element;
234  typedef typename GridView::template Codim<0>::Iterator ElementIterator;
235 
236  typedef Ewoms::EnergyModule<TypeTag, enableEnergy> EnergyModule;
237 
238 public:
239  PvsModel(Simulator &simulator)
240  : ParentType(simulator)
241  {
242  verbosity_ = EWOMS_GET_PARAM(TypeTag, int, PvsVerbosity);
243  numSwitched_ = 0;
244  }
245 
249  static void registerParameters()
250  {
252 
253  // register runtime parameters of the VTK output modules
256 
257  if (enableDiffusion)
259 
260  if (enableEnergy)
262 
263  EWOMS_REGISTER_PARAM(TypeTag, int, PvsVerbosity,
264  "The verbosity level of the primary variable "
265  "switching model");
266  }
267 
271  static std::string name()
272  { return "pvs"; }
273 
277  std::string primaryVarName(int pvIdx) const
278  {
279  std::string s;
280  if (!(s = EnergyModule::primaryVarName(pvIdx)).empty())
281  return s;
282 
283  std::ostringstream oss;
284  if (pvIdx == Indices::pressure0Idx)
285  oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
286  else if (Indices::switch0Idx <= pvIdx && pvIdx < Indices::switch0Idx
287  + numPhases - 1)
288  oss << "switch_" << pvIdx - Indices::switch0Idx;
289  else if (Indices::switch0Idx + numPhases - 1 <= pvIdx
290  && pvIdx < Indices::switch0Idx + numComponents - 1)
291  oss << "auxMoleFrac^" << FluidSystem::componentName(pvIdx);
292  else
293  assert(false);
294 
295  return oss.str();
296  }
297 
301  std::string eqName(int eqIdx) const
302  {
303  std::string s;
304  if (!(s = EnergyModule::eqName(eqIdx)).empty())
305  return s;
306 
307  std::ostringstream oss;
308  if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
309  + numComponents) {
310  int compIdx = eqIdx - Indices::conti0EqIdx;
311  oss << "continuity^" << FluidSystem::componentName(compIdx);
312  }
313  else
314  assert(false);
315 
316  return oss.str();
317  }
318 
323  {
324  ParentType::updateFailed();
325  numSwitched_ = 0;
326  }
327 
331  void updateBegin()
332  {
333  ParentType::updateBegin();
334 
335  // find the a reference pressure. The first degree of freedom
336  // might correspond to non-interior entities which would lead
337  // to an undefined value, so we have to iterate...
338  int nDof = this->numTotalDof();
339  for (int dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
340  if (this->dofTotalVolume(dofIdx) > 0) {
342  this->solution(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
343  break;
344  }
345  }
346  }
347 
351  Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
352  {
353  Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
354  if (tmp > 0)
355  // energy related quantity
356  return tmp;
357 
358  if (Indices::pressure0Idx == pvIdx) {
359  return 10 / referencePressure_;
360  }
361 
362  if (Indices::switch0Idx <= pvIdx && pvIdx < Indices::switch0Idx
363  + numPhases - 1) {
364  int phaseIdx = pvIdx - Indices::switch0Idx;
365 
366  if (!this->solution(/*timeIdx=*/0)[globalDofIdx].phaseIsPresent(phaseIdx))
367  // for saturations, the weight is always 1
368  return 1;
369 
370  // for saturations, the PvsMoleSaturationsBaseWeight
371  // property determines the weight
372  return GET_PROP_VALUE(TypeTag, PvsSaturationsBaseWeight);
373  }
374 
375  // for mole fractions, the PvsMoleFractionsBaseWeight
376  // property determines the weight
377  return GET_PROP_VALUE(TypeTag, PvsMoleFractionsBaseWeight);
378  }
379 
383  Scalar eqWeight(int globalDofIdx, int eqIdx) const
384  {
385  Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
386  if (tmp > 0)
387  // energy related equation
388  return tmp;
389 
390  int compIdx = eqIdx - Indices::conti0EqIdx;
391  assert(0 <= compIdx && compIdx <= numComponents);
392 
393  // make all kg equal
394  return FluidSystem::molarMass(compIdx);
395  }
396 
401  {
402  ParentType::advanceTimeLevel();
403  numSwitched_ = 0;
404  }
405 
410  bool switched() const
411  { return numSwitched_ > 0; }
412 
416  template <class DofEntity>
417  void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
418  {
419  // write primary variables
420  ParentType::serializeEntity(outstream, dofEntity);
421 
422 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
423  int dofIdx = this->dofMapper().index(dofEntity);
424 #else
425  int dofIdx = this->dofMapper().map(dofEntity);
426 #endif
427  if (!outstream.good())
428  OPM_THROW(std::runtime_error, "Could not serialize DOF " << dofIdx);
429 
430  outstream << this->solution(/*timeIdx=*/0)[dofIdx].phasePresence() << " ";
431  }
432 
436  template <class DofEntity>
437  void deserializeEntity(std::istream &instream, const DofEntity &dofEntity)
438  {
439  // read primary variables
440  ParentType::deserializeEntity(instream, dofEntity);
441 
442  // read phase presence
443 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
444  int dofIdx = this->dofMapper().index(dofEntity);
445 #else
446  int dofIdx = this->dofMapper().map(dofEntity);
447 #endif
448  if (!instream.good())
449  OPM_THROW(std::runtime_error,
450  "Could not deserialize DOF " << dofIdx);
451 
452  int tmp;
453  instream >> tmp;
454  this->solution(/*timeIdx=*/0)[dofIdx].setPhasePresence(tmp);
455  this->solution(/*timeIdx=*/1)[dofIdx].setPhasePresence(tmp);
456  }
457 
466  {
467  numSwitched_ = 0;
468 
469  std::vector<bool> visited(this->numGridDof(), false);
470  ElementContext elemCtx(this->simulator_);
471 
472  ElementIterator elemIt = this->gridView_.template begin<0>();
473  ElementIterator elemEndIt = this->gridView_.template end<0>();
474  for (; elemIt != elemEndIt; ++elemIt) {
475  const Element& elem = *elemIt;
476  if (elem.partitionType() != Dune::InteriorEntity)
477  continue;
478  elemCtx.updateStencil(elem);
479 
480  int numLocalDof = elemCtx.stencil(/*timeIdx=*/0).numPrimaryDof();
481  for (int dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
482  int globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
483 
484  if (visited[globalIdx])
485  continue;
486  visited[globalIdx] = true;
487 
488  // compute the intensive quantities of the current degree of freedom
489  auto &priVars = this->solution(/*timeIdx=*/0)[globalIdx];
490  elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
491  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
492 
493  // evaluate primary variable switch
494  short oldPhasePresence = priVars.phasePresence();
495 
496  // set the primary variables and the new phase state
497  // from the current fluid state
498  priVars.assignNaive(intQuants.fluidState());
499 
500  if (oldPhasePresence != priVars.phasePresence()) {
501  if (verbosity_ > 1)
502  printSwitchedPhases_(elemCtx,
503  dofIdx,
504  intQuants.fluidState(),
505  oldPhasePresence,
506  priVars);
507  this->linearizer().markDofRed(globalIdx);
508  ++numSwitched_;
509  }
510  }
511  }
512 
513  // make sure that if there was a variable switch in an
514  // other partition we will also set the switch flag
515  // for our partition.
516  numSwitched_ = this->gridView_.comm().sum(numSwitched_);
517 
518  if (verbosity_ > 0)
519  this->simulator_.model().newtonMethod().endIterMsg()
520  << ", num switched=" << numSwitched_;
521  }
522 
523  template <class FluidState>
524  void printSwitchedPhases_(const ElementContext &elemCtx,
525  int dofIdx,
526  const FluidState &fs,
527  int oldPhasePresence,
528  const PrimaryVariables &newPv) const
529  {
530  typedef Opm::MathToolbox<typename FluidState::Scalar> FsToolbox;
531 
532  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
533  bool oldPhasePresent = (oldPhasePresence & (1 << phaseIdx)) > 0;
534  bool newPhasePresent = newPv.phaseIsPresent(phaseIdx);
535  if (oldPhasePresent == newPhasePresent)
536  continue;
537 
538  const auto &pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
539  if (oldPhasePresent && !newPhasePresent) {
540  std::cout << "'" << FluidSystem::phaseName(phaseIdx)
541  << "' phase disappears at position " << pos
542  << ". saturation=" << fs.saturation(phaseIdx)
543  << std::flush;
544  }
545  else {
546  Scalar sumx = 0;
547  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
548  sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
549 
550  std::cout << "'" << FluidSystem::phaseName(phaseIdx)
551  << "' phase appears at position " << pos
552  << " sum x = " << sumx << std::flush;
553  }
554  }
555 
556  std::cout << ", new primary variables: ";
557  newPv.print();
558  std::cout << "\n" << std::flush;
559  }
560 
562  {
564 
565  // add the VTK output modules which are meaningful for the model
566  this->addOutputModule(new Ewoms::VtkPhasePresenceModule<TypeTag>(this->simulator_));
567  this->addOutputModule(new Ewoms::VtkCompositionModule<TypeTag>(this->simulator_));
568  if (enableDiffusion)
569  this->addOutputModule(new Ewoms::VtkDiffusionModule<TypeTag>(this->simulator_));
570  if (enableEnergy)
571  this->addOutputModule(new Ewoms::VtkEnergyModule<TypeTag>(this->simulator_));
572  }
573 
574  mutable Scalar referencePressure_;
575 
576  // number of switches of the phase state in the last Newton
577  // iteration
579 
580  // verbosity of the model
582 };
583 }
584 
585 #endif
void registerOutputModules_()
Definition: multiphasebasemodel.hh:244
VTK output module for quantities which make sense for models which assume thermal equilibrium...
std::string primaryVarName(int pvIdx) const
Given an primary variable index, return a human readable name.
Definition: pvsmodel.hh:277
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
bool switched() const
Return true if the primary variables were switched for at least one vertex after the last timestep...
Definition: pvsmodel.hh:410
Classes required for molecular diffusion.
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: pvsintensivequantities.hh:53
SET_SCALAR_PROP(NumericModel, EndTime,-1e100)
The default value for the simulation's end time.
int verbosity_
Definition: pvsmodel.hh:581
A newton solver which is specific to the compositional multi-phase PVS model.
A newton solver which is specific to the compositional multi-phase PVS model.
Definition: pvsnewtonmethod.hh:40
Implements a vector representing molar, mass or volumetric rates.
Definition: pvsratevector.hh:48
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: pvsmodel.hh:331
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: pvsmodel.hh:322
The multi-dimensional Newton method.
Definition: newtonmethod.hh:54
Represents the primary variables used in the primary variable switching compositional model...
Definition: pvsprimaryvariables.hh:54
VTK output module for quantities which make sense for models which assume thermal equilibrium...
Definition: vtkenergymodule.hh:70
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
The indices for the compositional multi-phase primary variable switching model.
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:45
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Definition: pvsboundaryratevector.hh:44
int numSwitched_
Definition: pvsmodel.hh:578
Contains all data which is required to calculate all fluxes at a flux integration point for the prima...
Definition: pvsextensivequantities.hh:48
void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
Write the current solution for a degree of freedom to a restart file.
Definition: pvsmodel.hh:417
static std::string name()
Definition: pvsmodel.hh:271
std::string eqName(int eqIdx) const
Given an equation index, return a human readable name.
Definition: pvsmodel.hh:301
SET_INT_PROP(NumericModel, GridGlobalRefinements, 0)
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:93
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkphasepresencemodule.hh:77
static void registerParameters()
Register all run-time parameters for the PVS compositional model.
Definition: pvsmodel.hh:249
Implements a vector representing molar, mass or volumetric rates.
SET_TYPE_PROP(NumericModel, Scalar, double)
Set the default type of scalar values to double.
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
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Contains all data which is required to calculate all fluxes at a flux integration point for the prima...
Definition: baseauxiliarymodule.hh:35
VTK output module for the fluid composition.
A generic compositional multi-phase model using primary-variable switching.
Definition: pvsmodel.hh:57
void registerOutputModules_()
Definition: pvsmodel.hh:561
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
Scalar eqWeight(int globalDofIdx, int eqIdx) const
Returns the relative weight of an equation.
Definition: pvsmodel.hh:383
Represents the primary variables used in the primary variable switching compositional model...
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
VTK output module for the fluid composition.
Definition: vtkphasepresencemodule.hh:54
Scalar referencePressure_
Definition: pvsmodel.hh:574
PvsModel(Simulator &simulator)
Definition: pvsmodel.hh:239
void printSwitchedPhases_(const ElementContext &elemCtx, int dofIdx, const FluidState &fs, int oldPhasePresence, const PrimaryVariables &newPv) const
Definition: pvsmodel.hh:524
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:44
NEW_TYPE_TAG(AuxModule)
The indices for the compositional multi-phase primary variable switching model.
Definition: pvsindices.hh:43
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Declares the properties required for the compositional multi-phase primary variable switching model...
SET_BOOL_PROP(FvBaseDiscretization, EnableVtkOutput, true)
Enable the VTK output by default.
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition: pvsmodel.hh:400
void switchPrimaryVars_()
Definition: pvsmodel.hh:465
void deserializeEntity(std::istream &instream, const DofEntity &dofEntity)
Reads the current solution variables for a degree of freedom from a restart file. ...
Definition: pvsmodel.hh:437
VTK output module for quantities which make sense for models which incorperate molecular diffusion...
Definition: vtkdiffusionmodule.hh:67
#define INHERITS_FROM(...)
Syntactic sugar for NEW_TYPE_TAG.
Definition: propertysystem.hh:229
Scalar primaryVarWeight(int globalDofIdx, int pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: pvsmodel.hh:351
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:98
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
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:95