richardsintensivequantities.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_RICHARDS_INTENSIVE_QUANTITIES_HH
27 #define EWOMS_RICHARDS_INTENSIVE_QUANTITIES_HH
28 
29 #include "richardsproperties.hh"
30 
31 #include <opm/material/fluidstates/ImmiscibleFluidState.hpp>
32 
33 #include <dune/common/fvector.hh>
34 #include <dune/common/fmatrix.hh>
35 
36 namespace Ewoms {
37 
44 template <class TypeTag>
46  : public GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities)
47  , public GET_PROP_TYPE(TypeTag, FluxModule)::FluxIntensiveQuantities
48 {
49  typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
50  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
51  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
52  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
53  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
54  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
55  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
56  typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
57 
58  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
59  enum { pressureWIdx = Indices::pressureWIdx };
60  enum { numPhases = FluidSystem::numPhases };
61  enum { liquidPhaseIdx = GET_PROP_VALUE(TypeTag, LiquidPhaseIndex) };
62  enum { gasPhaseIdx = GET_PROP_VALUE(TypeTag, GasPhaseIndex) };
63  enum { dimWorld = GridView::dimensionworld };
64 
65  typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
66  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
67  typedef Dune::FieldVector<Scalar, numPhases> ScalarPhaseVector;
68  typedef Dune::FieldVector<Evaluation, numPhases> PhaseVector;
69  typedef Opm::MathToolbox<Evaluation> Toolbox;
70 
71 public:
73  typedef Opm::ImmiscibleFluidState<Evaluation, FluidSystem> FluidState;
74 
78  void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
79  {
80  ParentType::update(elemCtx, dofIdx, timeIdx);
81 
82  const auto& T = elemCtx.problem().temperature(elemCtx, dofIdx, timeIdx);
83  fluidState_.setTemperature(T);
84 
85  // material law parameters
86  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
87  const auto &problem = elemCtx.problem();
88  const typename MaterialLaw::Params &materialParams =
89  problem.materialLawParams(elemCtx, dofIdx, timeIdx);
90  const auto &priVars = elemCtx.primaryVars(dofIdx, timeIdx);
91 
93  // calculate the pressures
95 
96  // first, we have to find the minimum capillary pressure (i.e. Sw = 0)
97  fluidState_.setSaturation(liquidPhaseIdx, 1.0);
98  fluidState_.setSaturation(gasPhaseIdx, 0.0);
99  ScalarPhaseVector pC;
100  MaterialLaw::capillaryPressures(pC, materialParams, fluidState_);
101 
102  // non-wetting pressure can be larger than the
103  // reference pressure if the medium is fully
104  // saturated by the wetting phase
105  const Evaluation& pW = priVars.makeEvaluation(pressureWIdx, timeIdx);
106  Evaluation pN =
107  Toolbox::max(elemCtx.problem().referencePressure(elemCtx, dofIdx, /*timeIdx=*/0),
108  pW + (pC[gasPhaseIdx] - pC[liquidPhaseIdx]));
109 
111  // calculate the saturations
113  fluidState_.setPressure(liquidPhaseIdx, pW);
114  fluidState_.setPressure(gasPhaseIdx, pN);
115 
116  PhaseVector sat;
117  MaterialLaw::saturations(sat, materialParams, fluidState_);
118  fluidState_.setSaturation(liquidPhaseIdx, sat[liquidPhaseIdx]);
119  fluidState_.setSaturation(gasPhaseIdx, sat[gasPhaseIdx]);
120 
121  typename FluidSystem::ParameterCache paramCache;
122  paramCache.updateAll(fluidState_);
123 
124  // compute and set the wetting phase viscosity
125  const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, liquidPhaseIdx);
126  fluidState_.setViscosity(liquidPhaseIdx, mu);
127  fluidState_.setViscosity(gasPhaseIdx, 1e-20);
128 
129  // compute and set the wetting phase density
130  const Evaluation& rho = FluidSystem::density(fluidState_, paramCache, liquidPhaseIdx);
131  fluidState_.setDensity(liquidPhaseIdx, rho);
132  fluidState_.setDensity(gasPhaseIdx, 1e-20);
133 
134  // relperms
135  MaterialLaw::relativePermeabilities(relativePermeability_, materialParams, fluidState_);
136 
137  // mobilities
138  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
139  mobility_[phaseIdx] = relativePermeability_[phaseIdx]/fluidState_.viscosity(phaseIdx);
140 
141  // porosity
142  porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
143 
144  // intrinsic permeability
145  intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
146 
147  // update the quantities specific for the velocity model
148  FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
149  }
150 
154  const FluidState &fluidState() const
155  { return fluidState_; }
156 
160  const Evaluation& porosity() const
161  { return porosity_; }
162 
166  const DimMatrix &intrinsicPermeability() const
167  { return intrinsicPerm_; }
168 
172  const Evaluation& relativePermeability(int phaseIdx) const
173  { return relativePermeability_[phaseIdx]; }
174 
178  const Evaluation& mobility(int phaseIdx) const
179  { return mobility_[phaseIdx]; }
180 
181 private:
182  FluidState fluidState_;
183  DimMatrix intrinsicPerm_;
184  Evaluation relativePermeability_[numPhases];
185  Evaluation mobility_[numPhases];
186  Evaluation porosity_;
187 };
188 
189 } // namespace Ewoms
190 
191 #endif
void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
Definition: richardsintensivequantities.hh:78
Opm::ImmiscibleFluidState< Evaluation, FluidSystem > FluidState
The type returned by the fluidState() method.
Definition: richardsintensivequantities.hh:73
#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
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: richardsintensivequantities.hh:160
const Evaluation & mobility(int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: richardsintensivequantities.hh:178
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: richardsintensivequantities.hh:166
Intensive quantities required by the Richards model.
Definition: richardsintensivequantities.hh:45
Definition: baseauxiliarymodule.hh:35
const Evaluation & relativePermeability(int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: richardsintensivequantities.hh:172
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: richardsintensivequantities.hh:154
Contains the property declarations for the Richards model.