blackoilboundaryratevector.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_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
27 #define EWOMS_BLACK_OIL_BOUNDARY_RATE_VECTOR_HH
28 
29 #include <opm/material/common/Valgrind.hpp>
30 #include <opm/material/constraintsolvers/NcpFlash.hpp>
31 
33 
34 namespace Ewoms {
35 
41 template <class TypeTag>
42 class BlackOilBoundaryRateVector : public GET_PROP_TYPE(TypeTag, RateVector)
43 {
44  typedef typename GET_PROP_TYPE(TypeTag, RateVector) ParentType;
45  typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) ExtensiveQuantities;
46  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
47  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
48  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
49 
50  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
51  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
52  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
53  enum { conti0EqIdx = Indices::conti0EqIdx };
54 
55 public:
59  BlackOilBoundaryRateVector() : ParentType()
60  {}
61 
65  BlackOilBoundaryRateVector(Scalar value) : ParentType(value)
66  {}
67 
72  : ParentType(value)
73  {}
74 
75 
79  template <class Context, class FluidState>
80  void setFreeFlow(const Context &context, int bfIdx, int timeIdx,
81  const FluidState &fluidState)
82  {
83  typename FluidSystem::ParameterCache paramCache;
84  paramCache.updateAll(fluidState);
85 
86  ExtensiveQuantities extQuants;
87  extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState, paramCache);
88  const auto &insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
89 
91  // advective fluxes of all components in all phases
93  (*this) = 0.0;
94  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
95  Scalar meanMBoundary = 0;
96  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
97  meanMBoundary += fluidState.moleFraction(phaseIdx, compIdx)
98  * FluidSystem::molarMass(compIdx);
99 
100  Scalar density;
101  if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
102  density = FluidSystem::density(fluidState, paramCache, phaseIdx);
103  else
104  density = insideIntQuants.fluidState().density(phaseIdx);
105 
106  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
107  Scalar molarity;
108  if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx))
109  molarity =
110  fluidState.moleFraction(phaseIdx, compIdx)
111  * density
112  / meanMBoundary;
113  else
114  molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
115 
116  // add advective flux of current component in current
117  // phase
118  (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx) * molarity;
119  }
120  }
121 
122 #ifndef NDEBUG
123  for (int i = 0; i < numEq; ++i) {
124  Valgrind::CheckDefined((*this)[i]);
125  }
126  Valgrind::CheckDefined(*this);
127 #endif
128  }
129 
133  template <class Context, class FluidState>
134  void setInFlow(const Context &context, int bfIdx, int timeIdx,
135  const FluidState &fluidState)
136  {
137  this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
138 
139  // we only allow fluxes in the direction opposite to the outer
140  // unit normal
141  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
142  Scalar &val = this->operator[](eqIdx);
143  val = std::min<Scalar>(0.0, val);
144  }
145  }
146 
150  template <class Context, class FluidState>
151  void setOutFlow(const Context &context, int bfIdx, int timeIdx,
152  const FluidState &fluidState)
153  {
154  this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
155 
156  // we only allow fluxes in the same direction as the outer
157  // unit normal
158  for (int eqIdx = 0; eqIdx < numEq; ++eqIdx) {
159  Scalar &val = this->operator[](eqIdx);
160  val = std::max<Scalar>(0.0, val);
161  }
162  }
163 
167  void setNoFlow()
168  { (*this) = 0.0; }
169 };
170 
171 } // namespace Ewoms
172 
173 #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 setNoFlow()
Specify a no-flow boundary for all conserved quantities.
Definition: blackoilboundaryratevector.hh:167
BlackOilBoundaryRateVector(const BlackOilBoundaryRateVector &value)
Copy constructor.
Definition: blackoilboundaryratevector.hh:71
Implements a boundary vector for the fully implicit black-oil model.
Definition: blackoilboundaryratevector.hh:42
void setFreeFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify a free-flow boundary.
Definition: blackoilboundaryratevector.hh:80
BlackOilBoundaryRateVector(Scalar value)
Definition: blackoilboundaryratevector.hh:65
Definition: baseauxiliarymodule.hh:35
void setOutFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an outflow boundary.
Definition: blackoilboundaryratevector.hh:151
BlackOilBoundaryRateVector()
Default constructor.
Definition: blackoilboundaryratevector.hh:59
Contains the quantities which are are constant within a finite volume in the black-oil model...
void setInFlow(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState)
Specify an inflow boundary.
Definition: blackoilboundaryratevector.hh:134