pvsratevector.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_PVS_RATE_VECTOR_HH
29#define EWOMS_PVS_RATE_VECTOR_HH
30
31#include "pvsindices.hh"
32
34#include <opm/material/constraintsolvers/NcpFlash.hpp>
35#include <opm/material/common/Valgrind.hpp>
36
37#include <dune/common/fvector.hh>
38
39namespace Opm {
40
49template <class TypeTag>
51 : public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
52 getPropValue<TypeTag, Properties::NumEq>()>
53{
58
59 enum { conti0EqIdx = Indices::conti0EqIdx };
60 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
61 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
62 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
63
64 using ParentType = Dune::FieldVector<Evaluation, numEq>;
66
67public:
68 PvsRateVector() : ParentType()
69 { Opm::Valgrind::SetUndefined(*this); }
70
74 PvsRateVector(const Evaluation& value)
75 : ParentType(value)
76 {}
77
81 PvsRateVector(const PvsRateVector& value) : ParentType(value)
82 {}
83
87 void setMassRate(const ParentType& value)
88 {
89 // convert to molar rates
90 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
91 (*this)[conti0EqIdx + compIdx] = value[conti0EqIdx + compIdx];
92 (*this)[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
93 }
94 }
95
99 void setMolarRate(const ParentType& value)
100 { ParentType::operator=(value); }
101
105 template <class RhsEval>
106 void setEnthalpyRate(const RhsEval& rate)
107 { EnergyModule::setEnthalpyRate(*this, rate); }
108
112 template <class FluidState, class RhsEval>
113 void setVolumetricRate(const FluidState& fluidState, unsigned phaseIdx, const RhsEval& volume)
114 {
115 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
116 (*this)[conti0EqIdx + compIdx] =
117 fluidState.density(phaseIdx, compIdx)
118 * fluidState.moleFraction(phaseIdx, compIdx)
119 * volume;
120
121 EnergyModule::setEnthalpyRate(*this, fluidState, phaseIdx, volume);
122 }
123
127 template <class RhsEval>
128 PvsRateVector& operator=(const RhsEval& value)
129 {
130 for (unsigned i=0; i < this->size(); ++i)
131 (*this)[i] = value;
132 return *this;
133 }
134
139 {
140 for (unsigned i=0; i < this->size(); ++i)
141 (*this)[i] = other[i];
142 return *this;
143 }
144};
145
146} // namespace Opm
147
148#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:48
Implements a vector representing molar, mass or volumetric rates.
Definition: pvsratevector.hh:53
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: pvsratevector.hh:113
PvsRateVector(const Evaluation &value)
Definition: pvsratevector.hh:74
PvsRateVector(const PvsRateVector &value)
Copy constructor.
Definition: pvsratevector.hh:81
void setEnthalpyRate(const RhsEval &rate)
Set an enthalpy rate [J/As] where .
Definition: pvsratevector.hh:106
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition: pvsratevector.hh:99
PvsRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: pvsratevector.hh:128
PvsRateVector & operator=(const PvsRateVector &other)
Assignment operator from another rate vector.
Definition: pvsratevector.hh:138
PvsRateVector()
Definition: pvsratevector.hh:68
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition: pvsratevector.hh:87
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:242