multiphasebasemodel.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) 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_MULTI_PHASE_BASE_MODEL_HH
27 #define EWOMS_MULTI_PHASE_BASE_MODEL_HH
28 
29 #include <opm/material/localad/Math.hpp>
31 #include "multiphasebaseproblem.hh"
33 
36 
37 #include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
38 #include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
39 #include <opm/material/heatconduction/DummyHeatConductionLaw.hpp>
40 
41 #include <dune/common/unused.hh>
42 
43 namespace Ewoms {
44 template <class TypeTag>
46 }
47 
48 namespace Ewoms {
49 namespace Properties {
51 NEW_TYPE_TAG(MultiPhaseBaseModel, INHERITS_FROM(VtkMultiPhase, VtkTemperature));
52 
54 SET_SPLICES(MultiPhaseBaseModel, SpatialDiscretizationSplice);
55 
59 SET_TAG_PROP(MultiPhaseBaseModel, SpatialDiscretizationSplice, VcfvDiscretization);
60 
62 SET_INT_PROP(MultiPhaseBaseModel, NumEq, GET_PROP_TYPE(TypeTag, Indices)::numEq);
64 SET_INT_PROP(MultiPhaseBaseModel, NumPhases, GET_PROP_TYPE(TypeTag, FluidSystem)::numPhases);
66 SET_INT_PROP(MultiPhaseBaseModel, NumComponents, GET_PROP_TYPE(TypeTag, FluidSystem)::numComponents);
67 
70 
73 
78 {
79 private:
80  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
81  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
82  typedef Opm::NullMaterialTraits<Scalar, FluidSystem::numPhases> Traits;
83 
84 public:
85  typedef Opm::NullMaterial<Traits> type;
86 };
87 
93  MaterialLawParams,
94  typename GET_PROP_TYPE(TypeTag, MaterialLaw)::Params);
95 
98  HeatConductionLaw,
99  Opm::DummyHeatConductionLaw<typename GET_PROP_TYPE(TypeTag, Scalar)>);
100 
104  HeatConductionLawParams,
105  typename GET_PROP_TYPE(TypeTag, HeatConductionLaw)::Params);
106 
108 SET_BOOL_PROP(MultiPhaseBaseModel, EnableGravity, false);
109 
110 } // namespace Properties
111 
117 template <class TypeTag>
118 class MultiPhaseBaseModel : public GET_PROP_TYPE(TypeTag, Discretization)
119 {
120  typedef typename GET_PROP_TYPE(TypeTag, Discretization) ParentType;
121  typedef typename GET_PROP_TYPE(TypeTag, Model) Implementation;
122  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
123  typedef typename GET_PROP_TYPE(TypeTag, ThreadManager) ThreadManager;
124  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
125  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
126  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
127  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
128  typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
129  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
130 
131  typedef typename GridView::template Codim<0>::Iterator ElementIterator;
132  typedef typename GridView::template Codim<0>::Entity Element;
133 
134  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
135  enum { numComponents = FluidSystem::numComponents };
136 
137 public:
139  : ParentType(simulator)
140  { }
141 
145  static void registerParameters()
146  {
147  ParentType::registerParameters();
148 
149  // register runtime parameters of the VTK output modules
152  }
153 
157  void finishInit()
158  {
159  ParentType::finishInit();
160 
161  // Reduce the Newton tolerance by the volume of smallest degree of freedom. (For
162  // large grids the tolerance needs to be reduced because the total mass lost
163  // would be too large.)
164  Scalar minDofVolume = 1e100;
165  for (size_t globalDofIdx = 0; globalDofIdx < this->numGridDof(); ++ globalDofIdx)
166  if (this->isLocalDof(globalDofIdx))
167  minDofVolume = std::min(minDofVolume, this->dofTotalVolume(globalDofIdx));
168  minDofVolume = this->gridView().comm().min(minDofVolume);
169 
170  Scalar newtonTolerance = EWOMS_GET_PARAM(TypeTag, Scalar, NewtonRawTolerance);
171  newtonTolerance /= std::sqrt(minDofVolume);
172  this->newtonMethod().setTolerance(newtonTolerance);
173  }
174 
180  bool phaseIsConsidered(int phaseIdx) const
181  { return true; }
182 
190  void globalPhaseStorage(EqVector &storage, int phaseIdx)
191  {
192  assert(0 <= phaseIdx && phaseIdx < numPhases);
193 
194  storage = 0;
195 
196  ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView());
197  OmpMutex addMutex;
198 #ifdef _OPENMP
199 #pragma omp parallel
200 #endif
201  {
202  // Attention: the variables below are thread specific and thus cannot be
203  // moved in front of the #pragma!
204  int threadId = ThreadManager::threadId();
205  ElementContext elemCtx(this->simulator_);
206  ElementIterator elemIt = this->gridView().template begin</*codim=*/0>();
207  EqVector tmp;
208 
209  for (threadedElemIt.beginParallel(elemIt);
210  !threadedElemIt.isFinished(elemIt);
211  threadedElemIt.increment(elemIt))
212  {
213  const Element& elem = *elemIt;
214  if (elem.partitionType() != Dune::InteriorEntity)
215  continue; // ignore ghost and overlap elements
216 
217  elemCtx.updateStencil(elem);
218  elemCtx.updateIntensiveQuantities(/*timeIdx=*/0);
219 
220  const auto &stencil = elemCtx.stencil(/*timeIdx=*/0);
221 
222  for (int dofIdx = 0; dofIdx < elemCtx.numDof(/*timeIdx=*/0); ++dofIdx) {
223  const auto &scv = stencil.subControlVolume(dofIdx);
224  const auto &intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
225 
226  tmp = 0;
227  this->localResidual(threadId).addPhaseStorage(tmp,
228  elemCtx,
229  dofIdx,
230  /*timeIdx=*/0,
231  phaseIdx);
232  tmp *= scv.volume()*intQuants.extrusionFactor();
233 
234  OmpMutex addLock(addMutex);
235  storage += tmp;
236  addLock.unlock();
237  }
238  }
239  }
240 
241  storage = this->gridView_.comm().sum(storage);
242  }
243 
245  {
246  ParentType::registerOutputModules_();
247 
248  // add the VTK output modules which make sense for all multi-phase models
249  this->addOutputModule(new Ewoms::VtkMultiPhaseModule<TypeTag>(this->simulator_));
250  this->addOutputModule(new Ewoms::VtkTemperatureModule<TypeTag>(this->simulator_));
251  }
252 
253 private:
254  const Implementation &asImp_() const
255  { return *static_cast<const Implementation *>(this); }
256 };
257 } // namespace Ewoms
258 
259 #endif
void registerOutputModules_()
Definition: multiphasebasemodel.hh:244
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkmultiphasemodule.hh:131
void beginParallel(EntityIterator &threadPrivateIt)
Definition: threadedentityiterator.hh:53
void globalPhaseStorage(EqVector &storage, int phaseIdx)
Compute the total storage inside one phase of all conservation quantities.
Definition: multiphasebasemodel.hh:190
SET_TAG_PROP(FvBaseDiscretization, LinearSolverSplice, ParallelIterativeLinearSolver)
use a parallel iterative linear solver by default
bool phaseIsConsidered(int phaseIdx) const
Returns true iff a fluid phase is used by the model.
Definition: multiphasebasemodel.hh:180
Specifies a flux module which uses the Darcy relation.
Definition: darcyfluxmodule.hh:58
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:45
void unlock()
Definition: locks.hh:41
This file contains the necessary classes to calculate the velocity out of a pressure potential gradie...
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
This class calculates the pressure potential gradients and the filter velocities for multi-phase flow...
SET_PROP(NumericModel, ParameterTree)
Set the ParameterTree property.
Definition: basicproperties.hh:117
SET_INT_PROP(NumericModel, GridGlobalRefinements, 0)
VTK output module for quantities which make sense for all models which deal with multiple fluid phase...
Definition: vtkmultiphasemodule.hh:93
MultiPhaseBaseModel(Simulator &simulator)
Definition: multiphasebasemodel.hh:138
SET_TYPE_PROP(NumericModel, Scalar, double)
Set the default type of scalar values to double.
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
Definition: baseauxiliarymodule.hh:35
void finishInit()
Apply the initial conditions to the model.
Definition: multiphasebasemodel.hh:157
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:145
The base class for the vertex centered finite volume discretization scheme.
Definition: vcfvdiscretization.hh:42
VTK output module for the temperature in which assume thermal equilibrium.
Definition: vtktemperaturemodule.hh:59
Provides an STL-iterator like interface to iterate over the enties of a GridView in OpenMP threaded a...
Definition: threadedentityiterator.hh:39
static int threadId()
Return the index of the current OpenMP thread.
Definition: threadmanager.hh:123
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:55
Implements a shallow wrapper around the "raw" locks provided by OpenMP.
Definition: locks.hh:35
The base class for the vertex centered finite volume discretization scheme.
NEW_TYPE_TAG(AuxModule)
void increment(EntityIterator &threadPrivateIt)
Definition: threadedentityiterator.hh:68
Defines the common properties required by the porous medium multi-phase models.
bool isFinished(const EntityIterator &threadPrivateIt) const
Definition: threadedentityiterator.hh:63
SET_BOOL_PROP(FvBaseDiscretization, EnableVtkOutput, true)
Enable the VTK output by default.
#define INHERITS_FROM(...)
Syntactic sugar for NEW_TYPE_TAG.
Definition: propertysystem.hh:229
SET_SPLICES(FvBaseDiscretization, LinearSolverSplice, LocalLinearizerSplice)
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtktemperaturemodule.hh:82
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:95