flashlocalresidual.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 
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_FLASH_LOCAL_RESIDUAL_HH
27 #define EWOMS_FLASH_LOCAL_RESIDUAL_HH
28 
29 #include "flashproperties.hh"
30 
33 
34 namespace Ewoms {
41 template <class TypeTag>
42 class FlashLocalResidual: public GET_PROP_TYPE(TypeTag, DiscLocalResidual)
43 {
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  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
48  typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
49  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
50 
51  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
52  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
53  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
54  enum { conti0EqIdx = Indices::conti0EqIdx };
55 
56  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
58 
59  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
61 
62  typedef Opm::MathToolbox<Evaluation> Toolbox;
63 
64 public:
68  template <class LhsEval>
69  void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
70  const ElementContext &elemCtx,
71  int dofIdx,
72  int timeIdx,
73  int phaseIdx) const
74  {
75  const IntensiveQuantities &intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
76  const auto &fs = intQuants.fluidState();
77 
78  // compute storage term of all components within all phases
79  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
80  int eqIdx = conti0EqIdx + compIdx;
81  storage[eqIdx] +=
82  Toolbox::template toLhs<LhsEval>(fs.molarity(phaseIdx, compIdx))
83  * Toolbox::template toLhs<LhsEval>(fs.saturation(phaseIdx))
84  * Toolbox::template toLhs<LhsEval>(intQuants.porosity());
85  }
86 
87  EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
88  }
89 
93  template <class LhsEval>
94  void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
95  const ElementContext &elemCtx,
96  int dofIdx,
97  int timeIdx) const
98  {
99  storage = 0;
100  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
101  addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
102 
103  EnergyModule::addSolidHeatStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
104  }
105 
109  void computeFlux(RateVector &flux,
110  const ElementContext &elemCtx,
111  int scvfIdx,
112  int timeIdx) const
113  {
114  flux = 0.0;
115  addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
116  Valgrind::CheckDefined(flux);
117 
118  addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
119  Valgrind::CheckDefined(flux);
120  }
121 
125  void addAdvectiveFlux(RateVector &flux,
126  const ElementContext &elemCtx,
127  int scvfIdx,
128  int timeIdx) const
129  {
130  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
131 
132  int interiorIdx = extQuants.interiorIndex();
133  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
134  // data attached to upstream and the finite volume of the current phase
135  int upIdx = extQuants.upstreamIndex(phaseIdx);
136  const IntensiveQuantities &up = elemCtx.intensiveQuantities(upIdx, timeIdx);
137 
138  // this is a bit hacky because it is specific to the element-centered
139  // finite volume scheme. (N.B. that if finite differences are used to
140  // linearize the system of equations, it does not matter.)
141  if (upIdx == interiorIdx) {
142  Evaluation tmp =
143  up.fluidState().molarDensity(phaseIdx)
144  * extQuants.volumeFlux(phaseIdx);
145 
146  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
147  flux[conti0EqIdx + compIdx] +=
148  tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
149  }
150  }
151  else {
152  Evaluation tmp =
153  Toolbox::value(up.fluidState().molarDensity(phaseIdx))
154  * extQuants.volumeFlux(phaseIdx);
155 
156  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
157  flux[conti0EqIdx + compIdx] +=
158  tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
159  }
160  }
161  }
162 
163  EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
164  }
165 
169  void addDiffusiveFlux(RateVector &flux,
170  const ElementContext &elemCtx,
171  int scvfIdx,
172  int timeIdx) const
173  {
174  DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
175  EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
176  }
177 
181  void computeSource(RateVector &source,
182  const ElementContext &elemCtx,
183  int dofIdx,
184  int timeIdx) const
185  {
186  Valgrind::SetUndefined(source);
187  elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
188  Valgrind::CheckDefined(source);
189  }
190 };
191 
192 } // namespace Ewoms
193 
194 #endif
Classes required for molecular diffusion.
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:46
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: flashlocalresidual.hh:125
void computeSource(RateVector &source, const ElementContext &elemCtx, int dofIdx, int timeIdx) const
Calculate the source term of the equation.
Definition: flashlocalresidual.hh:181
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
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: flashlocalresidual.hh:94
Calculates the local residual of the compositional multi-phase model based on flash calculations...
Definition: flashlocalresidual.hh:42
#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: flashlocalresidual.hh:109
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, int scvfIdx, int timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: flashlocalresidual.hh:169
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Definition: baseauxiliarymodule.hh:35
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: flashlocalresidual.hh:69
Declares the properties required by the compositional multi-phase model based on flash calculations...
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...