flashratevector.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_FLASH_RATE_VECTOR_HH
27 #define EWOMS_FLASH_RATE_VECTOR_HH
28 
29 #include <dune/common/fvector.hh>
30 
32 #include <opm/material/constraintsolvers/NcpFlash.hpp>
33 #include <opm/material/common/Valgrind.hpp>
34 
36 
37 namespace Ewoms {
38 
44 template <class TypeTag>
46  : public Dune::FieldVector<typename GET_PROP_TYPE(TypeTag, Evaluation),
47  GET_PROP_VALUE(TypeTag, NumEq)>
48 {
49  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
50  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
51  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
52  typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
53 
54  enum { conti0EqIdx = Indices::conti0EqIdx };
55  enum { numComponents = GET_PROP_VALUE(TypeTag, NumComponents) };
56  enum { numEq = GET_PROP_VALUE(TypeTag, NumEq) };
57 
58  typedef Dune::FieldVector<Evaluation, numEq> ParentType;
60 
61 public:
62  FlashRateVector() : ParentType()
63  { Valgrind::SetUndefined(*this); }
64 
68  FlashRateVector(const Evaluation& value) : ParentType(value)
69  {}
70 
75  FlashRateVector(const FlashRateVector& value) : ParentType(value)
76  {}
77 
81  void setMassRate(const ParentType& value)
82  {
83  // convert to molar rates
84  ParentType molarRate(value);
85  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
86  molarRate[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
87 
88  setMolarRate(molarRate);
89  }
90 
94  void setMolarRate(const ParentType &value)
95  { ParentType::operator=(value); }
96 
100  void setEnthalpyRate(const Evaluation& rate)
101  { EnergyModule::setEnthalpyRate(*this, rate); }
102 
106  template <class FluidState, class RhsEval>
107  void setVolumetricRate(const FluidState &fluidState, int phaseIdx, const RhsEval& volume)
108  {
109  for (int compIdx = 0; compIdx < numComponents; ++compIdx)
110  (*this)[conti0EqIdx + compIdx] =
111  fluidState.density(phaseIdx, compIdx)
112  * fluidState.moleFraction(phaseIdx, compIdx)
113  * volume;
114 
115  EnergyModule::setEnthalpyRate(*this, fluidState, phaseIdx, volume);
116  }
117 
121  template <class RhsEval>
122  FlashRateVector& operator=(const RhsEval& value)
123  {
124  for (unsigned i=0; i < this->size(); ++i)
125  (*this)[i] = value;
126  return *this;
127  }
128 
133  {
134  for (unsigned i=0; i < this->size(); ++i)
135  (*this)[i] = other[i];
136  return *this;
137  }
138 };
139 
140 } // namespace Ewoms
141 
142 #endif
Implements a vector representing rates of conserved quantities.
Definition: flashratevector.hh:45
FlashRateVector(const FlashRateVector &value)
Default constructor.
Definition: flashratevector.hh:75
FlashRateVector(const Evaluation &value)
Definition: flashratevector.hh:68
void setVolumetricRate(const FluidState &fluidState, int phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: flashratevector.hh:107
#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.
void setEnthalpyRate(const Evaluation &rate)
Set an enthalpy rate [J/As] where .
Definition: flashratevector.hh:100
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition: flashratevector.hh:81
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:54
Definition: baseauxiliarymodule.hh:35
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition: flashratevector.hh:94
FlashRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: flashratevector.hh:122
FlashRateVector & operator=(const FlashRateVector &other)
Assignment operator from another rate vector.
Definition: flashratevector.hh:132
FlashRateVector()
Definition: flashratevector.hh:62
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...