pvsprimaryvariables.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) 2011-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_PRIMARY_VARIABLES_HH
27 #define EWOMS_PVS_PRIMARY_VARIABLES_HH
28 
31 
32 #include <opm/material/constraintsolvers/NcpFlash.hpp>
33 #include <opm/material/fluidstates/CompositionalFluidState.hpp>
34 
35 #include "pvsindices.hh"
36 #include "pvsproperties.hh"
37 
38 #include <dune/common/fvector.hh>
39 
40 #include <iostream>
41 
42 namespace Ewoms {
43 
53 template <class TypeTag>
55 {
56  typedef FvBasePrimaryVariables<TypeTag> ParentType;
58  typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) Implementation;
59 
60  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
61  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
62  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
63  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
64  typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
65  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
66 
67  // primary variable indices
68  enum { pressure0Idx = Indices::pressure0Idx };
69  enum { switch0Idx = Indices::switch0Idx };
70 
71  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
72  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
73  enum { enableEnergy = GET_PROP_VALUE(TypeTag, EnableEnergy) };
74 
75  typedef typename Opm::MathToolbox<Evaluation> Toolbox;
76  typedef Dune::FieldVector<Scalar, numComponents> ComponentVector;
78  typedef Opm::NcpFlash<Scalar, FluidSystem> NcpFlash;
79 
80 public:
81  PvsPrimaryVariables() : ParentType()
82  { Valgrind::SetDefined(*this); }
83 
87  explicit PvsPrimaryVariables(Scalar value) : ParentType(value)
88  {
89  Valgrind::CheckDefined(value);
90  Valgrind::SetDefined(*this);
91 
92  phasePresence_ = 0;
93  }
94 
99  PvsPrimaryVariables(const PvsPrimaryVariables &value) : ParentType(value)
100  {
101  Valgrind::SetDefined(*this);
102 
103  phasePresence_ = value.phasePresence_;
104  }
105 
109  template <class FluidState>
110  void assignMassConservative(const FluidState &fluidState,
111  const MaterialLawParams &matParams,
112  bool isInEquilibrium = false)
113  {
114 #ifndef NDEBUG
115  // make sure the temperature is the same in all fluid phases
116  for (int phaseIdx = 1; phaseIdx < numPhases; ++phaseIdx) {
117  assert(std::abs(fluidState.temperature(0) - fluidState.temperature(phaseIdx)) < 1e-30);
118  }
119 #endif // NDEBUG
120 
121  // for the equilibrium case, we don't need complicated
122  // computations.
123  if (isInEquilibrium) {
124  assignNaive(fluidState);
125  return;
126  }
127 
128  // use a flash calculation to calculate a fluid state in
129  // thermodynamic equilibrium
130  typename FluidSystem::ParameterCache paramCache;
131  Opm::CompositionalFluidState<Scalar, FluidSystem> fsFlash;
132 
133  // use the externally given fluid state as initial value for
134  // the flash calculation
135  fsFlash.assign(fluidState);
136 
137  // calculate the phase densities
138  paramCache.updateAll(fsFlash);
139  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
140  fsFlash.setDensity(phaseIdx, FluidSystem::density(fsFlash, paramCache,
141  phaseIdx));
142  }
143 
144  // calculate the "global molarities"
145  ComponentVector globalMolarities(0.0);
146  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
147  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
148  globalMolarities[compIdx] +=
149  fsFlash.saturation(phaseIdx) * fsFlash.molarity(phaseIdx, compIdx);
150  }
151  }
152 
153  // run the flash calculation
154  // NcpFlash::guessInitial(fsFlash, paramCache, globalMolarities);
155  NcpFlash::template solve<MaterialLaw>(fsFlash, paramCache, matParams,
156  globalMolarities);
157 
158  // use the result to assign the primary variables
159  assignNaive(fsFlash);
160  }
161 
166  short phasePresence() const
167  { return phasePresence_; }
168 
175  void setPhasePresence(short value)
176  { phasePresence_ = value; }
177 
185  void setPhasePresent(int phaseIdx, bool yesno = true)
186  {
187  if (yesno)
188  setPhasePresence(phasePresence_ | (1 << phaseIdx));
189  else
190  setPhasePresence(phasePresence_ & ~(1 << phaseIdx));
191  }
192 
197  unsigned char implicitSaturationIdx() const
198  { return lowestPresentPhaseIdx(); }
199 
208  static bool phaseIsPresent(int phaseIdx, short phasePresence)
209  { return phasePresence & (1 << phaseIdx); }
210 
217  bool phaseIsPresent(int phaseIdx) const
218  { return phasePresence_ & (1 << phaseIdx); }
219 
224  { return ffs(phasePresence_) - 1; }
225 
229  ThisType &operator=(const Implementation &value)
230  {
231  ParentType::operator=(value);
232  phasePresence_ = value.phasePresence_;
233 
234  return *this;
235  }
236 
240  ThisType &operator=(Scalar value)
241  {
242  ParentType::operator=(value);
243 
244  phasePresence_ = 0;
245  return *this;
246  }
247 
255  Evaluation explicitSaturationValue(int phaseIdx, int timeIdx) const
256  {
257  if (!phaseIsPresent(phaseIdx) || phaseIdx == lowestPresentPhaseIdx())
258  // non-present phases have saturation 0
259  return Toolbox::createConstant(0.0);
260 
261  int varIdx = switch0Idx + phaseIdx - 1;
262  if (timeIdx != 0)
263  Toolbox::createConstant((*this)[varIdx]);
264  return Toolbox::createVariable((*this)[varIdx], varIdx);
265  }
266 
270  template <class FluidState>
271  void assignNaive(const FluidState &fluidState)
272  {
273  typedef Opm::MathToolbox<typename FluidState::Scalar> FsToolbox;
274 
275  // assign the phase temperatures. this is out-sourced to
276  // the energy module
277  EnergyModule::setPriVarTemperatures(*this, fluidState);
278 
279  // set the pressure of the first phase
280  (*this)[pressure0Idx] = FsToolbox::value(fluidState.pressure(/*phaseIdx=*/0));
281  Valgrind::CheckDefined((*this)[pressure0Idx]);
282 
283  // determine the phase presence.
284  phasePresence_ = 0;
285  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
286  // use a NCP condition to determine if the phase is
287  // present or not
288  Scalar a = 1;
289  for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
290  a -= FsToolbox::value(fluidState.moleFraction(phaseIdx, compIdx));
291  }
292  Scalar b = FsToolbox::value(fluidState.saturation(phaseIdx));
293 
294  if (b > a)
295  phasePresence_ |= (1 << phaseIdx);
296  }
297 
298  // some phase must be present
299  if (phasePresence_ == 0)
300  OPM_THROW(Opm::NumericalProblem,
301  "Phase state was 0, i.e., no fluid is present");
302 
303  // set the primary variables which correspond to mole
304  // fractions of the present phase which has the lowest index.
305  int lowestPhaseIdx = lowestPresentPhaseIdx();
306  for (int switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
307  int phaseIdx = switchIdx;
308  int compIdx = switchIdx + 1;
309  if (switchIdx >= lowestPhaseIdx)
310  ++phaseIdx;
311 
312  if (phaseIsPresent(phaseIdx)) {
313  (*this)[switch0Idx + switchIdx] = FsToolbox::value(fluidState.saturation(phaseIdx));
314  Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
315  }
316  else {
317  (*this)[switch0Idx + switchIdx] =
318  FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx));
319  Valgrind::CheckDefined((*this)[switch0Idx + switchIdx]);
320  }
321  }
322 
323  // set the mole fractions in of the remaining components in
324  // the phase with the lowest index
325  for (int compIdx = numPhases - 1; compIdx < numComponents - 1;
326  ++compIdx) {
327  (*this)[switch0Idx + compIdx] =
328  FsToolbox::value(fluidState.moleFraction(lowestPhaseIdx, compIdx + 1));
329  Valgrind::CheckDefined((*this)[switch0Idx + compIdx]);
330  }
331  }
332 
336  void print(std::ostream &os = std::cout) const
337  {
338  os << "(p_" << FluidSystem::phaseName(0) << " = "
339  << this->operator[](pressure0Idx);
340  int lowestPhaseIdx = lowestPresentPhaseIdx();
341  for (int switchIdx = 0; switchIdx < numPhases - 1; ++switchIdx) {
342  int phaseIdx = switchIdx;
343  int compIdx = switchIdx + 1;
344  if (phaseIdx >= lowestPhaseIdx)
345  ++phaseIdx; // skip the saturation of the present
346  // phase with the lowest index
347 
348  if (phaseIsPresent(phaseIdx)) {
349  os << ", S_" << FluidSystem::phaseName(phaseIdx) << " = "
350  << (*this)[switch0Idx + switchIdx];
351  }
352  else {
353  os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
354  << FluidSystem::componentName(compIdx) << " = "
355  << (*this)[switch0Idx + switchIdx];
356  }
357  }
358  for (int compIdx = numPhases - 1; compIdx < numComponents - 1;
359  ++compIdx) {
360  os << ", x_" << FluidSystem::phaseName(lowestPhaseIdx) << "^"
361  << FluidSystem::componentName(compIdx + 1) << " = "
362  << (*this)[switch0Idx + compIdx];
363  }
364  os << ")";
365  os << ", phase presence: " << static_cast<int>(phasePresence_) << std::flush;
366  }
367 
368 private:
369  unsigned char phasePresence_;
370 };
371 
372 } // namespace Ewoms
373 
374 #endif
Evaluation explicitSaturationValue(int phaseIdx, int timeIdx) const
Returns an explcitly stored saturation for a given phase.
Definition: pvsprimaryvariables.hh:255
bool phaseIsPresent(int phaseIdx) const
Returns true iff a phase is present for the current phase presence.
Definition: pvsprimaryvariables.hh:217
static bool phaseIsPresent(int phaseIdx, short phasePresence)
Returns true iff a phase is present for a given phase presence.
Definition: pvsprimaryvariables.hh:208
Represents the primary variables used by the a model.
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &matParams, bool isInEquilibrium=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: pvsprimaryvariables.hh:110
Represents the primary variables used in the primary variable switching compositional model...
Definition: pvsprimaryvariables.hh:54
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
ThisType & operator=(Scalar value)
Assignment operator from a scalar value.
Definition: pvsprimaryvariables.hh:240
short phasePresence() const
Return the fluid phases which are present in a given control volume.
Definition: pvsprimaryvariables.hh:166
The indices for the compositional multi-phase primary variable switching model.
ThisType & operator=(const Implementation &value)
Assignment operator from an other primary variables object.
Definition: pvsprimaryvariables.hh:229
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: pvsprimaryvariables.hh:271
PvsPrimaryVariables(const PvsPrimaryVariables &value)
Default constructor.
Definition: pvsprimaryvariables.hh:99
Represents the primary variables used by the a model.
Definition: fvbaseprimaryvariables.hh:41
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
void setPhasePresent(int phaseIdx, bool yesno=true)
Set whether a given indivividual phase should be present or not.
Definition: pvsprimaryvariables.hh:185
Definition: baseauxiliarymodule.hh:35
int lowestPresentPhaseIdx() const
Returns the phase with the lowest index that is present.
Definition: pvsprimaryvariables.hh:223
void setPhasePresence(short value)
Set which fluid phases are present in a given control volume.
Definition: pvsprimaryvariables.hh:175
unsigned char implicitSaturationIdx() const
Returns the index of the phase with's its saturation is determined by the closure condition of satura...
Definition: pvsprimaryvariables.hh:197
PvsPrimaryVariables()
Definition: pvsprimaryvariables.hh:81
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: pvsprimaryvariables.hh:336
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...
PvsPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: pvsprimaryvariables.hh:87