blackoillocalresidual.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_BLACK_OIL_LOCAL_RESIDUAL_HH
27 #define EWOMS_BLACK_OIL_LOCAL_RESIDUAL_HH
28 
29 #include "blackoilproperties.hh"
30 
31 namespace Ewoms {
37 template <class TypeTag>
38 class BlackOilLocalResidual : public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
39 {
40  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
41  typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
42  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
43  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
44  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
45  typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
46  typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
47 
48  enum { conti0EqIdx = Indices::conti0EqIdx };
49  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
50  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
51  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
52 
53  typedef Opm::MathToolbox<Evaluation> Toolbox;
54 
55 public:
59  template <class LhsEval>
60  void computeStorage(Dune::FieldVector<LhsEval, numEq> &storage,
61  const ElementContext &elemCtx,
62  int dofIdx,
63  int timeIdx) const
64  {
65 #ifndef NDEBUG
66  typedef Opm::MathToolbox<LhsEval> LhsToolbox;
67 #endif
68 
69  // retrieve the intensive quantities for the SCV at the specified point in time
70  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
71 
72  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
73  storage[eqIdx] = 0.0;
74 
75  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
76  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
77  storage[conti0EqIdx + compIdx] +=
78  Toolbox::template toLhs<LhsEval>(intQuants.porosity())
79  * Toolbox::template toLhs<LhsEval>(intQuants.fluidState().saturation(phaseIdx))
80  * Toolbox::template toLhs<LhsEval>(intQuants.fluidState().molarity(phaseIdx, compIdx));
81  assert(std::isfinite(LhsToolbox::value(storage[conti0EqIdx + compIdx])));
82  }
83  }
84  }
85 
89  void computeFlux(RateVector &flux,
90  const ElementContext &elemCtx,
91  int scvfIdx,
92  int timeIdx) const
93  {
94  assert(timeIdx == 0);
95 
96  const ExtensiveQuantities &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
97 
98  for (int eqIdx=0; eqIdx < numEq; eqIdx++)
99  flux[eqIdx] = 0;
100 
101  int interiorIdx = extQuants.interiorIndex();
102  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
103  int upIdx = extQuants.upstreamIndex(phaseIdx);
104  const IntensiveQuantities &up = elemCtx.intensiveQuantities(upIdx, /*timeIdx=*/0);
105 
106  // this is a bit hacky because it is specific to the element-centered
107  // finite volume scheme. (N.B. that if finite differences are used to
108  // linearize the system of equations, it does not matter.)
109  if (upIdx == interiorIdx) {
110  Evaluation tmp =
111  up.fluidState().molarDensity(phaseIdx)
112  * extQuants.volumeFlux(phaseIdx);
113 
114  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
115  flux[conti0EqIdx + compIdx] +=
116  tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
117  }
118  }
119  else {
120  Evaluation tmp =
121  Toolbox::value(up.fluidState().molarDensity(phaseIdx))
122  * extQuants.volumeFlux(phaseIdx);
123 
124  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
125  flux[conti0EqIdx + compIdx] +=
126  tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
127  }
128  }
129  }
130  }
131 
135  void computeSource(RateVector &source,
136  const ElementContext &elemCtx,
137  int dofIdx,
138  int timeIdx) const
139  {
140  // retrieve the source term intrinsic to the problem
141  elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
142  }
143 };
144 
145 } // namespace Ewoms
146 
147 #endif
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:38
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: blackoillocalresidual.hh:89
#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
Declares the properties required by the black oil model.
Definition: baseauxiliarymodule.hh:35
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: blackoillocalresidual.hh:60
void computeSource(RateVector &source, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Calculate the source term of the equation.
Definition: blackoillocalresidual.hh:135