flashboundaryratevector.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_BOUNDARY_RATE_VECTOR_HH
27 #define EWOMS_FLASH_BOUNDARY_RATE_VECTOR_HH
28 
29 #include "flashproperties.hh"
30 
32 #include <opm/material/common/Valgrind.hpp>
33 
34 namespace Ewoms {
35 
43 template <class TypeTag>
44 class FlashBoundaryRateVector : public GET_PROP_TYPE(TypeTag, RateVector)
45 {
46  typedef typename GET_PROP_TYPE(TypeTag, RateVector) ParentType;
47  typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
48  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
49  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
50  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
51  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
52 
53  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
54  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
55  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
56  enum { conti0EqIdx = Indices::conti0EqIdx };
57  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
58 
60  typedef Opm::MathToolbox<Evaluation> Toolbox;
61 
62 public:
63  FlashBoundaryRateVector() : ParentType()
64  {}
65 
70  FlashBoundaryRateVector(const Evaluation& value) : ParentType(value)
71  {}
72 
78  : ParentType(value)
79  {}
80 
84  template <class Context, class FluidState>
85  void setFreeFlow(const Context &context, int bfIdx, int timeIdx,
86  const FluidState &fluidState)
87  {
88  typename FluidSystem::ParameterCache paramCache;
89  paramCache.updateAll(fluidState);
90 
91  ExtensiveQuantities extQuants;
92  extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState, paramCache);
93  const auto &insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
94 
96  // advective fluxes of all components in all phases
98  (*this) = Toolbox::createConstant(0.0);
99  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
100  Evaluation meanMBoundary = 0;
101  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
102  meanMBoundary +=
103  fluidState.moleFraction(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
104 
105  Evaluation density;
106  if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
107  density = FluidSystem::density(fluidState, paramCache, phaseIdx);
108  else
109  density = insideIntQuants.fluidState().density(phaseIdx);
110 
111  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
112  Evaluation molarity;
113  if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
114  molarity = fluidState.moleFraction(phaseIdx, compIdx)*density/meanMBoundary;
115  else
116  molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
117 
118  // add advective flux of current component in current
119  // phase
120  (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
121  }
122 
123  if (enableEnergy) {
124  Evaluation specificEnthalpy;
125  if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
126  specificEnthalpy = FluidSystem::enthalpy(fluidState, paramCache, phaseIdx);
127  else
128  specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
129 
130  // currently we neglect heat conduction!
131  Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
132  EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
133  }
134  }
135 
136  EnergyModule::addToEnthalpyRate(*this, EnergyModule::heatConductionRate(extQuants));
137 
138 #ifndef NDEBUG
139  for (int i = 0; i < numEq; ++i) {
140  Valgrind::CheckDefined((*this)[i]);
141  }
142 #endif
143  }
144 
148  template <class Context, class FluidState>
149  void setInFlow(const Context &context, int bfIdx, int timeIdx,
150  const FluidState &fluidState)
151  {
152  this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
153 
154  // we only allow fluxes in the direction opposite to the outer unit normal
155  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
156  Evaluation &val = this->operator[](eqIdx);
157  val = Toolbox::min(0.0, val);
158  }
159  }
160 
164  template <class Context, class FluidState>
165  void setOutFlow(const Context &context, int bfIdx, int timeIdx,
166  const FluidState &fluidState)
167  {
168  this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
169 
170  // we only allow fluxes in the same direction as the outer unit normal
171  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
172  Evaluation &val = this->operator[](eqIdx);
173  val = Toolbox::max(0.0, val);
174  }
175  }
176 
180  void setNoFlow()
181  { (*this) = Toolbox::createConstant(0.0); }
182 };
183 
184 } // namespace Ewoms
185 
186 #endif
#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
void setOutFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: flashboundaryratevector.hh:165
void setInFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: flashboundaryratevector.hh:149
Implements a boundary vector for the fully implicit compositional multi-phase model which is based on...
Definition: flashboundaryratevector.hh:44
FlashBoundaryRateVector(const Evaluation &value)
Definition: flashboundaryratevector.hh:70
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Definition: baseauxiliarymodule.hh:35
FlashBoundaryRateVector()
Definition: flashboundaryratevector.hh:63
Declares the properties required by the compositional multi-phase model based on flash calculations...
void setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: flashboundaryratevector.hh:180
FlashBoundaryRateVector(const FlashBoundaryRateVector &value)
Definition: flashboundaryratevector.hh:77
void setFreeFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: flashboundaryratevector.hh:85
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...