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 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_FLASH_RATE_VECTOR_HH
29#define EWOMS_FLASH_RATE_VECTOR_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/material/common/Valgrind.hpp>
34
36
37namespace Opm {
38
44template <class TypeTag>
46 : public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
47 getPropValue<TypeTag, Properties::NumEq>()>
48{
53
54 enum { conti0EqIdx = Indices::conti0EqIdx };
55 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
56 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
57
58 using ParentType = Dune::FieldVector<Evaluation, numEq>;
60
61public:
63 : ParentType()
64 { Valgrind::SetUndefined(*this); }
65
69 explicit FlashRateVector(const Evaluation& value)
70 : ParentType(value)
71 {}
72
78 : ParentType(value)
79 {}
80
84 void setMassRate(const ParentType& value)
85 {
86 // convert to molar rates
87 ParentType molarRate(value);
88 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
89 molarRate[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
90 }
91
92 setMolarRate(molarRate);
93 }
94
98 void setMolarRate(const ParentType& value)
99 { ParentType::operator=(value); }
100
104 void setEnthalpyRate(const Evaluation& rate)
105 { EnergyModule::setEnthalpyRate(*this, rate); }
106
110 template <class FluidState, class RhsEval>
111 void setVolumetricRate(const FluidState& fluidState, unsigned phaseIdx, const RhsEval& volume)
112 {
113 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
114 (*this)[conti0EqIdx + compIdx] =
115 fluidState.density(phaseIdx, compIdx) *
116 fluidState.moleFraction(phaseIdx, compIdx) *
117 volume;
118 }
119
120 EnergyModule::setEnthalpyRate(*this, fluidState, phaseIdx, volume);
121 }
122
126 template <class RhsEval>
127 FlashRateVector& operator=(const RhsEval& value)
128 {
129 for (unsigned i = 0; i < this->size(); ++i) {
130 (*this)[i] = value;
131 }
132 return *this;
133 }
134
139 {
140 for (unsigned i = 0; i < this->size(); ++i) {
141 (*this)[i] = other[i];
142 }
143 return *this;
144 }
145};
146
147} // namespace Opm
148
149#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Definition: flashratevector.hh:48
void setEnthalpyRate(const Evaluation &rate)
Set an enthalpy rate [J/As] where .
Definition: flashratevector.hh:104
FlashRateVector & operator=(const FlashRateVector &other)
Assignment operator from another rate vector.
Definition: flashratevector.hh:138
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition: flashratevector.hh:98
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: flashratevector.hh:111
FlashRateVector(const Evaluation &value)
Definition: flashratevector.hh:69
FlashRateVector(const FlashRateVector &value)
Definition: flashratevector.hh:77
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition: flashratevector.hh:84
FlashRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: flashratevector.hh:127
FlashRateVector()
Definition: flashratevector.hh:62
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:39
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:233