pvsintensivequantities.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_PVS_INTENSIVE_QUANTITIES_HH
27 #define EWOMS_PVS_INTENSIVE_QUANTITIES_HH
28 
29 #include "pvsproperties.hh"
30 
33 
34 #include <opm/material/constraintsolvers/ComputeFromReferencePhase.hpp>
35 #include <opm/material/constraintsolvers/MiscibleMultiPhaseComposition.hpp>
36 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
37 
38 #include <dune/common/fvector.hh>
39 #include <dune/common/fmatrix.hh>
40 
41 #include <iostream>
42 
43 namespace Ewoms {
52 template <class TypeTag>
54  : public GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities)
55  , public DiffusionIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableDiffusion) >
56  , public EnergyIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableEnergy) >
57  , public GET_PROP_TYPE(TypeTag, FluxModule)::FluxIntensiveQuantities
58 {
59  typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
60 
61  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
62  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
63  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
64  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
65  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
66  typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
67  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
68  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
69  typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
70 
71  enum { switch0Idx = Indices::switch0Idx };
72  enum { pressure0Idx = Indices::pressure0Idx };
73  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
74  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
75  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
76  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
77  enum { dimWorld = GridView::dimensionworld };
78 
79  typedef Opm::MathToolbox<Evaluation> Toolbox;
80  typedef Opm::MiscibleMultiPhaseComposition<Scalar, FluidSystem, Evaluation> MiscibleMultiPhaseComposition;
81  typedef Opm::ComputeFromReferencePhase<Scalar, FluidSystem, Evaluation> ComputeFromReferencePhase;
82 
83  typedef Dune::FieldVector<Scalar, numPhases> PhaseVector;
84  typedef Dune::FieldVector<Evaluation, numPhases> EvalPhaseVector;
85  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
86 
87  typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
90 
91 public:
93  typedef Opm::CompositionalFluidState<Evaluation, FluidSystem> FluidState;
94 
98  void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
99  {
100  ParentType::update(elemCtx, dofIdx, timeIdx);
101  EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
102 
103  const auto &priVars = elemCtx.primaryVars(dofIdx, timeIdx);
104  const auto &problem = elemCtx.problem();
105 
107  // set the saturations
109  Evaluation sumSat = Toolbox::createConstant(0.0);
110  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
111  fluidState_.setSaturation(phaseIdx, priVars.explicitSaturationValue(phaseIdx, timeIdx));
112  Valgrind::CheckDefined(fluidState_.saturation(phaseIdx));
113  sumSat += fluidState_.saturation(phaseIdx);
114  }
115  Valgrind::CheckDefined(priVars.implicitSaturationIdx());
116  Valgrind::CheckDefined(sumSat);
117  fluidState_.setSaturation(priVars.implicitSaturationIdx(), 1.0 - sumSat);
118 
120  // set the pressures of the fluid phases
122 
123  // calculate capillary pressure
124  const MaterialLawParams &materialParams =
125  problem.materialLawParams(elemCtx, dofIdx, timeIdx);
126  EvalPhaseVector pC;
127  MaterialLaw::capillaryPressures(pC, materialParams, fluidState_);
128 
129  // set the absolute phase pressures in the fluid state
130  const Evaluation& p0 = priVars.makeEvaluation(pressure0Idx, timeIdx);
131  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
132  fluidState_.setPressure(phaseIdx, p0 + (pC[phaseIdx] - pC[0]));
133 
135  // calculate the phase compositions
137 
138  typename FluidSystem::ParameterCache paramCache;
139  int lowestPresentPhaseIdx = priVars.lowestPresentPhaseIdx();
140  int numNonPresentPhases = 0;
141  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
142  if (!priVars.phaseIsPresent(phaseIdx))
143  ++numNonPresentPhases;
144  }
145 
146  // now comes the tricky part: calculate phase compositions
147  if (numNonPresentPhases == numPhases - 1) {
148  // only one phase is present, i.e. the primary variables
149  // contain the complete composition of the phase
150  Evaluation sumx = Toolbox::createConstant(0.0);
151  for (int compIdx = 1; compIdx < numComponents; ++compIdx) {
152  const Evaluation& x = priVars.makeEvaluation(switch0Idx + compIdx - 1, timeIdx);
153  fluidState_.setMoleFraction(lowestPresentPhaseIdx, compIdx, x);
154  sumx += x;
155  }
156 
157  // set the mole fraction of the first component
158  fluidState_.setMoleFraction(lowestPresentPhaseIdx, 0, 1 - sumx);
159 
160  // calculate the composition of the remaining phases (as
161  // well as the densities of all phases). this is the job
162  // of the "ComputeFromReferencePhase" constraint solver
163  ComputeFromReferencePhase::solve(fluidState_, paramCache,
164  lowestPresentPhaseIdx,
165  /*setViscosity=*/true,
166  /*setEnthalpy=*/false);
167  }
168  else {
169  // create the auxiliary constraints
170  int numAuxConstraints = numComponents + numNonPresentPhases - numPhases;
171  Opm::MMPCAuxConstraint<Evaluation> auxConstraints[numComponents];
172 
173  int auxIdx = 0;
174  int switchIdx = 0;
175  for (; switchIdx < numPhases - 1; ++switchIdx) {
176  int compIdx = switchIdx + 1;
177  int switchPhaseIdx = switchIdx;
178  if (switchIdx >= lowestPresentPhaseIdx)
179  switchPhaseIdx += 1;
180 
181  if (!priVars.phaseIsPresent(switchPhaseIdx)) {
182  auxConstraints[auxIdx].set(lowestPresentPhaseIdx, compIdx,
183  priVars.makeEvaluation(switch0Idx + switchIdx, timeIdx));
184  ++auxIdx;
185  }
186  }
187 
188  for (; auxIdx < numAuxConstraints; ++auxIdx, ++switchIdx) {
189  int compIdx = numPhases - numNonPresentPhases + auxIdx;
190  auxConstraints[auxIdx].set(lowestPresentPhaseIdx, compIdx,
191  priVars.makeEvaluation(switch0Idx + switchIdx, timeIdx));
192  }
193 
194  // both phases are present, i.e. phase compositions are a result of the the
195  // gas <-> liquid equilibrium. This is the job of the
196  // "MiscibleMultiPhaseComposition" constraint solver
197  MiscibleMultiPhaseComposition::solve(fluidState_, paramCache,
198  priVars.phasePresence(),
199  auxConstraints,
200  numAuxConstraints,
201  /*setViscosity=*/true,
202  /*setEnthalpy=*/false);
203  }
204 
205 #ifndef NDEBUG
206  // make valgrind happy and set the enthalpies to NaN
207  if (!enableEnergy) {
208  Scalar myNan = std::numeric_limits<Scalar>::quiet_NaN();
209  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
210  fluidState_.setEnthalpy(phaseIdx, myNan);
211  }
212 #endif
213 
215  // calculate the remaining quantities
217 
218  // calculate relative permeabilities
219  MaterialLaw::relativePermeabilities(relativePermeability_,
220  materialParams, fluidState_);
221  Valgrind::CheckDefined(relativePermeability_);
222 
223  // mobilities
224  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
225  mobility_[phaseIdx] =
226  relativePermeability_[phaseIdx] / fluidState().viscosity(phaseIdx);
227 
228  // porosity
229  porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
230  Valgrind::CheckDefined(porosity_);
231 
232  // intrinsic permeability
233  intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
234 
235  // update the quantities specific for the velocity model
236  FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
237 
238  // energy related quantities
239  EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
240 
241  // update the diffusion specific quantities of the intensive quantities
242  DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
243 
244  fluidState_.checkDefined();
245  }
246 
250  const FluidState &fluidState() const
251  { return fluidState_; }
252 
256  const DimMatrix &intrinsicPermeability() const
257  { return intrinsicPerm_; }
258 
262  const Evaluation& relativePermeability(int phaseIdx) const
263  { return relativePermeability_[phaseIdx]; }
264 
268  const Evaluation& mobility(int phaseIdx) const
269  { return mobility_[phaseIdx]; }
270 
274  const Evaluation& porosity() const
275  { return porosity_; }
276 
277 private:
278  FluidState fluidState_;
279  Evaluation porosity_;
280  DimMatrix intrinsicPerm_;
281  Evaluation relativePermeability_[numPhases];
282  Evaluation mobility_[numPhases];
283 };
284 
285 } // namespace Ewoms
286 
287 #endif
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:530
Classes required for molecular diffusion.
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: pvsintensivequantities.hh:53
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
Definition: pvsintensivequantities.hh:98
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
const FluidState & fluidState() const
Returns the phase state for the control-volume.
Definition: pvsintensivequantities.hh:250
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: pvsintensivequantities.hh:274
Definition: baseauxiliarymodule.hh:35
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: pvsintensivequantities.hh:256
const Evaluation & mobility(int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: pvsintensivequantities.hh:268
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes...
Definition: diffusionmodule.hh:138
const Evaluation & relativePermeability(int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: pvsintensivequantities.hh:262
Opm::CompositionalFluidState< Evaluation, FluidSystem > FluidState
The type of the object returned by the fluidState() method.
Definition: pvsintensivequantities.hh:93
Declares the properties required for the compositional multi-phase primary variable switching model...
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...