ncpintensivequantities.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_NCP_INTENSIVE_QUANTITIES_HH
27 #define EWOMS_NCP_INTENSIVE_QUANTITIES_HH
28 
29 #include "ncpproperties.hh"
30 
33 
34 #include <opm/material/constraintsolvers/NcpFlash.hpp>
35 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
36 #include <opm/material/constraintsolvers/CompositionFromFugacities.hpp>
37 
38 #include <dune/common/fvector.hh>
39 #include <dune/common/fmatrix.hh>
40 
41 namespace Ewoms {
49 template <class TypeTag>
51  : public GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities)
52  , public DiffusionIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableDiffusion) >
53  , public EnergyIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableEnergy) >
54  , public GET_PROP_TYPE(TypeTag, FluxModule)::FluxIntensiveQuantities
55 {
56  typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
57 
58  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
59  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
60  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
61  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
62  typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
63  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
64  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
65 
66  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
67  typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
68 
69  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
70  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
71  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
72  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
73  enum { fugacity0Idx = Indices::fugacity0Idx };
74  enum { saturation0Idx = Indices::saturation0Idx };
75  enum { pressure0Idx = Indices::pressure0Idx };
76  enum { dimWorld = GridView::dimensionworld };
77 
78  typedef Opm::CompositionFromFugacities<Scalar, FluidSystem, Evaluation>
79  CompositionFromFugacitiesSolver;
80  typedef Opm::CompositionalFluidState<Evaluation, FluidSystem,
81  /*storeEnthalpy=*/enableEnergy> FluidState;
82  typedef Dune::FieldVector<Evaluation, numComponents> ComponentVector;
83  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
86  typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
87 
88 public:
90  {}
91 
95  void update(const ElementContext &elemCtx,
96  int dofIdx,
97  int timeIdx)
98  {
99  ParentType::update(elemCtx,
100  dofIdx,
101  timeIdx);
102  ParentType::checkDefined();
103 
104  typename FluidSystem::ParameterCache paramCache;
105  const auto &priVars = elemCtx.primaryVars(dofIdx, timeIdx);
106 
107  // set the phase saturations
108  Evaluation sumSat = 0;
109  for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx) {
110  const Evaluation& val = priVars.makeEvaluation(saturation0Idx + phaseIdx, timeIdx);
111  fluidState_.setSaturation(phaseIdx, val);
112  sumSat += val;
113  }
114  fluidState_.setSaturation(numPhases - 1, 1.0 - sumSat);
115  Valgrind::CheckDefined(sumSat);
116 
117  // set the fluid phase temperature
118  EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
119 
120  // retrieve capillary pressure parameters
121  const auto &problem = elemCtx.problem();
122  const MaterialLawParams &materialParams =
123  problem.materialLawParams(elemCtx, dofIdx, timeIdx);
124  // calculate capillary pressures
125  Evaluation capPress[numPhases];
126  MaterialLaw::capillaryPressures(capPress, materialParams, fluidState_);
127  // add to the pressure of the first fluid phase
128  const Evaluation& pressure0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
129  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
130  fluidState_.setPressure(phaseIdx, pressure0 + (capPress[phaseIdx] - capPress[0]));
131 
132  ComponentVector fug;
133  // retrieve component fugacities
134  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
135  fug[compIdx] = priVars.makeEvaluation(fugacity0Idx + compIdx, timeIdx);
136 
137  // calculate phase compositions
138  const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
139  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
140  // initial guess
141  if (hint) {
142  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
143  // use the hint for the initial mole fraction!
144  const Evaluation& moleFracIJ = hint->fluidState().moleFraction(phaseIdx, compIdx);
145 
146  // set initial guess of the component's mole fraction
147  fluidState_.setMoleFraction(phaseIdx, compIdx, moleFracIJ);
148  }
149  }
150  else // !hint
151  CompositionFromFugacitiesSolver::guessInitial(fluidState_,
152  paramCache,
153  phaseIdx, fug);
154 
155  // calculate the phase composition from the component
156  // fugacities
157  CompositionFromFugacitiesSolver::solve(fluidState_, paramCache,
158  phaseIdx, fug);
159  }
160 
161  // porosity
162  porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
163  Valgrind::CheckDefined(porosity_);
164 
165  // relative permeabilities
166  MaterialLaw::relativePermeabilities(relativePermeability_,
167  materialParams, fluidState_);
168 
169  // dynamic viscosities
170  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
171  // viscosities
172  const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
173  fluidState_.setViscosity(phaseIdx, mu);
174 
175  mobility_[phaseIdx] = relativePermeability_[phaseIdx]/mu;
176  }
177 
178  // intrinsic permeability
179  intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
180 
181  // update the quantities specific for the velocity model
182  FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
183 
184  // energy related quantities
185  EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
186 
187  // update the diffusion specific quantities of the intensive quantities
188  DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
189 
190  checkDefined();
191  }
192 
196  const FluidState &fluidState() const
197  { return fluidState_; }
198 
202  const DimMatrix &intrinsicPermeability() const
203  { return intrinsicPerm_; }
204 
208  const Evaluation& relativePermeability(int phaseIdx) const
209  { return relativePermeability_[phaseIdx]; }
210 
214  const Evaluation& mobility(int phaseIdx) const
215  { return mobility_[phaseIdx]; }
216 
220  const Evaluation& porosity() const
221  { return porosity_; }
222 
226  void checkDefined() const
227  {
228 #if !defined NDEBUG && HAVE_VALGRIND
229  ParentType::checkDefined();
230 
231  Valgrind::CheckDefined(porosity_);
232  Valgrind::CheckDefined(relativePermeability_);
233 
234  fluidState_.checkDefined();
235 #endif
236  }
237 
238 private:
239  DimMatrix intrinsicPerm_;
240  FluidState fluidState_;
241  Evaluation porosity_;
242  Evaluation relativePermeability_[numPhases];
243  Evaluation mobility_[numPhases];
244 };
245 
246 } // namespace Ewoms
247 
248 #endif
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:530
Classes required for molecular diffusion.
const Evaluation & mobility(int phaseIdx) const
ImmiscibleIntensiveQuantities::mobility.
Definition: ncpintensivequantities.hh:214
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
NcpIntensiveQuantities()
Definition: ncpintensivequantities.hh:89
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
void checkDefined() const
IntensiveQuantities::checkDefined.
Definition: ncpintensivequantities.hh:226
const Evaluation & porosity() const
ImmiscibleIntensiveQuantities::porosity.
Definition: ncpintensivequantities.hh:220
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: ncpintensivequantities.hh:50
const FluidState & fluidState() const
ImmiscibleIntensiveQuantities::fluidState.
Definition: ncpintensivequantities.hh:196
Definition: baseauxiliarymodule.hh:35
void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
IntensiveQuantities::update.
Definition: ncpintensivequantities.hh:95
const DimMatrix & intrinsicPermeability() const
ImmiscibleIntensiveQuantities::intrinsicPermeability.
Definition: ncpintensivequantities.hh:202
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes...
Definition: diffusionmodule.hh:138
Declares the properties required for the NCP compositional multi-phase model.
const Evaluation & relativePermeability(int phaseIdx) const
ImmiscibleIntensiveQuantities::relativePermeability.
Definition: ncpintensivequantities.hh:208
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...