pvslocalresidual.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  Copyright (C) 2009-2012 by Klaus Mosthaf
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 2 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
27 #ifndef EWOMS_PVS_LOCAL_RESIDUAL_HH
28 #define EWOMS_PVS_LOCAL_RESIDUAL_HH
29 
30 #include "pvsproperties.hh"
31 
34 
35 namespace Ewoms {
36 
43 template <class TypeTag>
44 class PvsLocalResidual : public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
45 {
46  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
47  typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
48  typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
49  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
50  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
51  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
52 
53  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
54  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
55  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
56  enum { conti0EqIdx = Indices::conti0EqIdx };
57 
58  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
60 
61  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
63 
64  typedef Opm::MathToolbox<Evaluation> Toolbox;
65 
66 public:
70  template <class LhsEval>
71  void addPhaseStorage(Dune::FieldVector<LhsEval, numEq> &storage,
72  const ElementContext &elemCtx,
73  int dofIdx,
74  int timeIdx,
75  int phaseIdx) const
76  {
77  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
78  const auto &fs = intQuants.fluidState();
79 
80  // compute storage term of all components within all phases
81  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
82  int eqIdx = conti0EqIdx + compIdx;
83  storage[eqIdx] +=
84  Toolbox::template toLhs<LhsEval>(fs.molarity(phaseIdx, compIdx))
85  * Toolbox::template toLhs<LhsEval>(fs.saturation(phaseIdx))
86  * Toolbox::template toLhs<LhsEval>(intQuants.porosity());
87  }
88 
89  EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
90  }
91 
95  template <class LhsEval>
96  void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
97  const ElementContext& elemCtx,
98  int dofIdx,
99  int timeIdx) const
100  {
101  storage = 0.0;
102  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
103  addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
104 
105  EnergyModule::addSolidHeatStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
106  }
107 
111  void computeFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
112  {
113  flux = Toolbox::createConstant(0.0);
114  addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
115  Valgrind::CheckDefined(flux);
116 
117  addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
118  Valgrind::CheckDefined(flux);
119  }
120 
124  void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx,
125  int scvfIdx, int timeIdx) const
126  {
127  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
128 
129  int interiorIdx = extQuants.interiorIndex();
130  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
131  // data attached to upstream and the downstream DOFs
132  // of the current phase
133  int upIdx = extQuants.upstreamIndex(phaseIdx);
134  const IntensiveQuantities &up = elemCtx.intensiveQuantities(upIdx, timeIdx);
135 
136  // this is a bit hacky because it is specific to the element-centered
137  // finite volume scheme. (N.B. that if finite differences are used to
138  // linearize the system of equations, it does not matter.)
139  if (upIdx == interiorIdx) {
140  Evaluation tmp =
141  up.fluidState().molarDensity(phaseIdx)
142  * extQuants.volumeFlux(phaseIdx);
143 
144  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
145  flux[conti0EqIdx + compIdx] +=
146  tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
147  }
148  }
149  else {
150  Evaluation tmp =
151  Toolbox::value(up.fluidState().molarDensity(phaseIdx))
152  * extQuants.volumeFlux(phaseIdx);
153 
154  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
155  flux[conti0EqIdx + compIdx] +=
156  tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
157  }
158  }
159  }
160 
161  EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
162  }
163 
167  void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx,
168  int scvfIdx, int timeIdx) const
169  {
170  DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
171  EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
172  }
173 
177  void computeSource(RateVector &source,
178  const ElementContext &elemCtx,
179  int dofIdx,
180  int timeIdx) const
181  {
182  Valgrind::SetUndefined(source);
183  elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
184  Valgrind::CheckDefined(source);
185  }
186 };
187 
188 } // namespace Ewoms
189 
190 #endif
Classes required for molecular diffusion.
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: pvslocalresidual.hh:124
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:46
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume...
Definition: pvslocalresidual.hh:96
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: pvslocalresidual.hh:167
void computeSource(RateVector &source, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Calculate the source term of the equation.
Definition: pvslocalresidual.hh:177
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Definition: baseauxiliarymodule.hh:35
void computeFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume...
Definition: pvslocalresidual.hh:111
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:44
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, int dofIdx, int timeIdx, int phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase...
Definition: pvslocalresidual.hh:71
Declares the properties required for the compositional multi-phase primary variable switching model...
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...