richardslocalresidual.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_RICHARDS_LOCAL_RESIDUAL_HH
27 #define EWOMS_RICHARDS_LOCAL_RESIDUAL_HH
28 
30 
32 
33 namespace Ewoms {
34 
39 template <class TypeTag>
40 class RichardsLocalResidual : public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
41 {
42  typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
43  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
44  typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
45  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
46  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
47  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
48 
49  enum { contiEqIdx = Indices::contiEqIdx };
50  enum { liquidPhaseIdx = GET_PROP_VALUE(TypeTag, LiquidPhaseIndex) };
51  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
52  typedef Opm::MathToolbox<Evaluation> Toolbox;
53 
54 public:
58  template <class LhsEval>
59  void computeStorage(Dune::FieldVector<LhsEval, numEq> &storage,
60  const ElementContext &elemCtx,
61  int dofIdx,
62  int timeIdx) const
63  {
64  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
65 
66  // partial time derivative of the wetting phase mass
67  storage[contiEqIdx] =
68  Toolbox::template toLhs<LhsEval>(intQuants.fluidState().density(liquidPhaseIdx))
69  *Toolbox::template toLhs<LhsEval>(intQuants.fluidState().saturation(liquidPhaseIdx))
70  *Toolbox::template toLhs<LhsEval>(intQuants.porosity());
71  }
72 
76  void computeFlux(RateVector &flux, const ElementContext &elemCtx,
77  int scvfIdx, int timeIdx) const
78  {
79  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
80 
81  int interiorIdx = extQuants.interiorIndex();
82  int upIdx = extQuants.upstreamIndex(liquidPhaseIdx);
83 
84  const IntensiveQuantities &up = elemCtx.intensiveQuantities(upIdx, timeIdx);
85 
86  // compute advective mass flux of the liquid phase. This is slightly hacky
87  // because it is specific to the element-centered finite volume method.
88  const Evaluation& rho = up.fluidState().density(liquidPhaseIdx);
89  if (interiorIdx == upIdx)
90  flux[contiEqIdx] = extQuants.volumeFlux(liquidPhaseIdx)*rho;
91  else
92  flux[contiEqIdx] = extQuants.volumeFlux(liquidPhaseIdx)*Toolbox::value(rho);
93  }
94 
98  void computeSource(RateVector &source,
99  const ElementContext &elemCtx,
100  int dofIdx,
101  int timeIdx) const
102  { elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx); }
103 };
104 
105 } // namespace Ewoms
106 
107 #endif
Element-wise calculation of the residual for the Richards model.
Definition: richardslocalresidual.hh:40
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
void computeSource(RateVector &source, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Calculate the source term of the equation.
Definition: richardslocalresidual.hh:98
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
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: richardslocalresidual.hh:76
Intensive quantities required by the Richards model.
Calculates and stores the data which is required to calculate the flux of fluid over a face of a fini...
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: richardslocalresidual.hh:59