ncplocalresidual.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) 2011 by Philipp Nuske
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_NCP_LOCAL_RESIDUAL_HH
28 #define EWOMS_NCP_LOCAL_RESIDUAL_HH
29 
30 #include "ncpproperties.hh"
31 
34 
35 namespace Ewoms {
42 template <class TypeTag>
43 class NcpLocalResidual : public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
44 {
45  typedef typename GET_PROP_TYPE(TypeTag, DiscLocalResidual) ParentType;
46  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
47  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
48  typedef typename GET_PROP_TYPE(TypeTag, EqVector) EqVector;
49  typedef typename GET_PROP_TYPE(TypeTag, RateVector) RateVector;
50  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
51  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
52  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
53 
54  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
55  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
56  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
57  enum { ncp0EqIdx = Indices::ncp0EqIdx };
58  enum { conti0EqIdx = Indices::conti0EqIdx };
59 
60  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
62 
63  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
65 
66  typedef Dune::FieldVector<Evaluation, numEq> EvalEqVector;
67  typedef Dune::BlockVector<EvalEqVector> ElemEvalEqVector;
68  typedef Opm::MathToolbox<Evaluation> Toolbox;
69 
70 public:
74  template <class LhsEval>
75  void addPhaseStorage(Dune::FieldVector<LhsEval, numEq> &storage,
76  const ElementContext &elemCtx,
77  int dofIdx,
78  int timeIdx,
79  int phaseIdx) const
80  {
81  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
82  const auto &fluidState = intQuants.fluidState();
83 
84  // compute storage term of all components within all phases
85  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
86  int eqIdx = conti0EqIdx + compIdx;
87  storage[eqIdx] +=
88  Toolbox::template toLhs<LhsEval>(fluidState.molarity(phaseIdx, compIdx))
89  * Toolbox::template toLhs<LhsEval>(fluidState.saturation(phaseIdx))
90  * Toolbox::template toLhs<LhsEval>(intQuants.porosity());
91  }
92 
93  EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
94  }
95 
99  template <class LhsEval>
100  void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
101  const ElementContext &elemCtx,
102  int dofIdx,
103  int timeIdx) const
104  {
105  storage = 0;
106  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
107  addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
108 
109  EnergyModule::addSolidHeatStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
110  }
111 
115  void computeFlux(RateVector &flux, const ElementContext &elemCtx,
116  int scvfIdx, int timeIdx) const
117  {
118  flux = Toolbox::createConstant(0.0);
119  addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
120  Valgrind::CheckDefined(flux);
121 
122  addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
123  Valgrind::CheckDefined(flux);
124  }
125 
129  void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx,
130  int scvfIdx, int timeIdx) const
131  {
132  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
133 
134  int interiorIdx = extQuants.interiorIndex();
135  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136  // data attached to upstream and the downstream DOFs
137  // of the current phase
138  int upIdx = extQuants.upstreamIndex(phaseIdx);
139  const IntensiveQuantities &up = elemCtx.intensiveQuantities(upIdx, timeIdx);
140 
141  // this is a bit hacky because it is specific to the element-centered
142  // finite volume scheme. (N.B. that if finite differences are used to
143  // linearize the system of equations, it does not matter.)
144  if (upIdx == interiorIdx) {
145  Evaluation tmp =
146  up.fluidState().molarDensity(phaseIdx)
147  * extQuants.volumeFlux(phaseIdx);
148 
149  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
150  flux[conti0EqIdx + compIdx] +=
151  tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
152  }
153  }
154  else {
155  Evaluation tmp =
156  Toolbox::value(up.fluidState().molarDensity(phaseIdx))
157  * extQuants.volumeFlux(phaseIdx);
158 
159  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
160  flux[conti0EqIdx + compIdx] +=
161  tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
162  }
163  }
164  }
165 
166  EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
167  }
168 
172  void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx,
173  int scvfIdx, int timeIdx) const
174  {
175  DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
176  EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
177  }
178 
185  void computeSource(RateVector &source,
186  const ElementContext &elemCtx,
187  int dofIdx,
188  int timeIdx) const
189  {
190  Valgrind::SetUndefined(source);
191  elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
192  Valgrind::CheckDefined(source);
193  }
194 
195 
199  void evalConstraints_(ElemEvalEqVector& residual,
200  ElemEvalEqVector& storageTerm,
201  const ElementContext& elemCtx, int timeIdx) const
202  {
203  // set the auxiliary functions
204  for (int dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
205  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
206  residual[dofIdx][ncp0EqIdx + phaseIdx] =
207  phaseNcp_(elemCtx, dofIdx, timeIdx, phaseIdx);
208  }
209  }
210 
211  // overwrite the constraint equations
212  ParentType::evalConstraints_(residual, storageTerm, elemCtx, timeIdx);
213  }
214 
218  Evaluation phaseNcp_(const ElementContext &elemCtx,
219  int dofIdx,
220  int timeIdx,
221  int phaseIdx) const
222  {
223  const auto &fluidState = elemCtx.intensiveQuantities(dofIdx, timeIdx).fluidState();
224 
225  const Evaluation& a = phaseNotPresentIneq_(fluidState, phaseIdx);
226  const Evaluation& b = phasePresentIneq_(fluidState, phaseIdx);
227  return Toolbox::min(a, b);
228  }
229 
234  template <class FluidState>
235  Evaluation phasePresentIneq_(const FluidState &fluidState, int phaseIdx) const
236  { return fluidState.saturation(phaseIdx); }
237 
242  template <class FluidState>
243  Evaluation phaseNotPresentIneq_(const FluidState &fluidState, int phaseIdx) const
244  {
245  // difference of sum of mole fractions in the phase from 100%
246  Evaluation a = 1;
247  for (int i = 0; i < numComponents; ++i)
248  a -= fluidState.moleFraction(phaseIdx, i);
249  return a;
250  }
251 };
252 
253 } // namespace Ewoms
254 
255 #endif
Evaluation phasePresentIneq_(const FluidState &fluidState, int phaseIdx) const
Returns the value of the inequality where a phase is present.
Definition: ncplocalresidual.hh:235
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: ncplocalresidual.hh:129
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: ncplocalresidual.hh:115
Classes required for molecular diffusion.
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:46
Evaluation phaseNotPresentIneq_(const FluidState &fluidState, int phaseIdx) const
Returns the value of the inequality where a phase is not present.
Definition: ncplocalresidual.hh:243
void evalConstraints_(ElemEvalEqVector &residual, ElemEvalEqVector &storageTerm, const ElementContext &elemCtx, int timeIdx) const
Set the values of the constraint volumes of the current element.
Definition: ncplocalresidual.hh:199
#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
Evaluation phaseNcp_(const ElementContext &elemCtx, int dofIdx, int timeIdx, int phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: ncplocalresidual.hh:218
void computeSource(RateVector &source, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Calculate the source term of the equation.
Definition: ncplocalresidual.hh:185
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: ncplocalresidual.hh:100
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Definition: baseauxiliarymodule.hh:35
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: ncplocalresidual.hh:172
Details needed to calculate the local residual in the compositional multi-phase NCP-model ...
Definition: ncplocalresidual.hh:43
Declares the properties required for the NCP compositional multi-phase model.
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: ncplocalresidual.hh:75
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...