flashintensivequantities.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_INTENSIVE_QUANTITIES_HH
27 #define EWOMS_FLASH_INTENSIVE_QUANTITIES_HH
28 
29 #include "flashproperties.hh"
30 #include "flashindices.hh"
31 
34 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
35 
36 #include <dune/common/fvector.hh>
37 #include <dune/common/fmatrix.hh>
38 
39 namespace Ewoms {
40 
47 template <class TypeTag>
49  : public GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities)
50  , public DiffusionIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableDiffusion) >
51  , public EnergyIntensiveQuantities<TypeTag, GET_PROP_VALUE(TypeTag, EnableEnergy) >
52  , public GET_PROP_TYPE(TypeTag, FluxModule)::FluxIntensiveQuantities
53 {
54  typedef typename GET_PROP_TYPE(TypeTag, DiscIntensiveQuantities) ParentType;
55 
56  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
57  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
58  typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
59  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
60  typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
61  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
62  typedef typename GET_PROP_TYPE(TypeTag, ThreadManager) ThreadManager;
63 
64  // primary variable indices
65  enum { cTot0Idx = Indices::cTot0Idx };
66  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
67  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
68  enum { enableDiffusion = GET_PROP_VALUE(TypeTag, EnableDiffusion) };
69  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
70  enum { dimWorld = GridView::dimensionworld };
71 
72  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
73  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
74  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
75  typedef typename GET_PROP_TYPE(TypeTag, FlashSolver) FlashSolver;
76 
77  typedef Dune::FieldVector<Evaluation, numComponents> ComponentVector;
78  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
79 
80  typedef typename FluxModule::FluxIntensiveQuantities FluxIntensiveQuantities;
83 
84 public:
86  typedef Opm::CompositionalFluidState<Evaluation, FluidSystem, enableEnergy> FluidState;
87 
91  void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
92  {
93  ParentType::update(elemCtx,
94  dofIdx,
95  timeIdx);
96  EnergyIntensiveQuantities::updateTemperatures_(fluidState_, elemCtx, dofIdx, timeIdx);
97 
98  const auto &priVars = elemCtx.primaryVars(dofIdx, timeIdx);
99  const auto &problem = elemCtx.problem();
100  Scalar flashTolerance = EWOMS_GET_PARAM(TypeTag, Scalar, FlashTolerance);
101  if (flashTolerance <= 0)
102  flashTolerance = std::numeric_limits<Scalar>::epsilon()*1e7;
103 
104  // extract the total molar densities of the components
105  ComponentVector cTotal;
106  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
107  cTotal[compIdx] = priVars.makeEvaluation(cTot0Idx + compIdx, timeIdx);
108 
109  typename FluidSystem::ParameterCache paramCache;
110  const auto *hint = elemCtx.thermodynamicHint(dofIdx, timeIdx);
111  if (hint) {
112  // use the same fluid state as the one of the hint, but
113  // make sure that we don't overwrite the temperature
114  // specified by the primary variables
115  const Evaluation& T = fluidState_.temperature(/*phaseIdx=*/0);
116  fluidState_.assign(hint->fluidState());
117  fluidState_.setTemperature(T);
118  }
119  else
120  FlashSolver::guessInitial(fluidState_, paramCache, cTotal);
121 
122  // compute the phase compositions, densities and pressures
123  const MaterialLawParams &materialParams =
124  problem.materialLawParams(elemCtx, dofIdx, timeIdx);
125  FlashSolver::template solve<MaterialLaw>(fluidState_,
126  paramCache,
127  materialParams,
128  cTotal,
129  flashTolerance);
130 
131  // calculate relative permeabilities
132  MaterialLaw::relativePermeabilities(relativePermeability_,
133  materialParams, fluidState_);
134  Valgrind::CheckDefined(relativePermeability_);
135 
136  // set the phase viscosities
137  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
138  const Evaluation& mu = FluidSystem::viscosity(fluidState_, paramCache, phaseIdx);
139  fluidState_.setViscosity(phaseIdx, mu);
140 
141  mobility_[phaseIdx] = relativePermeability_[phaseIdx] / mu;
142  Valgrind::CheckDefined(mobility_[phaseIdx]);
143  }
144 
146  // calculate the remaining quantities
148 
149  // porosity
150  porosity_ = problem.porosity(elemCtx, dofIdx, timeIdx);
151  Valgrind::CheckDefined(porosity_);
152 
153  // intrinsic permeability
154  intrinsicPerm_ = problem.intrinsicPermeability(elemCtx, dofIdx, timeIdx);
155 
156  // update the quantities specific for the velocity model
157  FluxIntensiveQuantities::update_(elemCtx, dofIdx, timeIdx);
158 
159  // energy related quantities
160  EnergyIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
161 
162  // update the diffusion specific quantities of the intensive quantities
163  DiffusionIntensiveQuantities::update_(fluidState_, paramCache, elemCtx, dofIdx, timeIdx);
164  }
165 
169  const FluidState &fluidState() const
170  { return fluidState_; }
171 
175  const DimMatrix &intrinsicPermeability() const
176  { return intrinsicPerm_; }
177 
181  const Evaluation& relativePermeability(int phaseIdx) const
182  { return relativePermeability_[phaseIdx]; }
183 
187  const Evaluation& mobility(int phaseIdx) const
188  {
189  return mobility_[phaseIdx];
190  }
191 
195  const Evaluation& porosity() const
196  { return porosity_; }
197 
198 private:
199  DimMatrix intrinsicPerm_;
200  FluidState fluidState_;
201  Evaluation porosity_;
202  Evaluation relativePermeability_[numPhases];
203  Evaluation mobility_[numPhases];
204 };
205 
206 } // namespace Ewoms
207 
208 #endif
Provides the volumetric quantities required for the energy equation.
Definition: energymodule.hh:530
Classes required for molecular diffusion.
const Evaluation & porosity() const
Returns the average porosity within the control volume.
Definition: flashintensivequantities.hh:195
Opm::CompositionalFluidState< Evaluation, FluidSystem, enableEnergy > FluidState
The type of the object returned by the fluidState() method.
Definition: flashintensivequantities.hh:86
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
Contains the intensive quantities of the flash-based compositional multi-phase model.
Definition: flashintensivequantities.hh:48
const DimMatrix & intrinsicPermeability() const
Returns the intrinsic permeability tensor a degree of freedom.
Definition: flashintensivequantities.hh:175
#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: flashintensivequantities.hh:169
void update(const ElementContext &elemCtx, int dofIdx, int timeIdx)
Definition: flashintensivequantities.hh:91
Simplifies multi-threaded capabilities.
Definition: threadmanager.hh:48
Definition: baseauxiliarymodule.hh:35
const Evaluation & mobility(int phaseIdx) const
Returns the effective mobility of a given phase within the control volume.
Definition: flashintensivequantities.hh:187
Provides the volumetric quantities required for the calculation of molecular diffusive fluxes...
Definition: diffusionmodule.hh:138
Defines the primary variable and equation indices for the compositional multi-phase model based on fl...
Declares the properties required by the compositional multi-phase model based on flash calculations...
const Evaluation & relativePermeability(int phaseIdx) const
Returns the relative permeability of a given phase within the control volume.
Definition: flashintensivequantities.hh:181
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:95