28 #ifndef EWOMS_BLACK_OIL_RATE_VECTOR_HH
29 #define EWOMS_BLACK_OIL_RATE_VECTOR_HH
31 #include <dune/common/fvector.hh>
33 #include <opm/common/Valgrind.hpp>
34 #include <opm/material/constraintsolvers/NcpFlash.hpp>
49 template <
class TypeTag>
51 :
public Dune::FieldVector<typename GET_PROP_TYPE(TypeTag, Evaluation),
52 GET_PROP_VALUE(TypeTag, NumEq)>
54 typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
55 typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
56 typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
57 typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
61 enum { conti0EqIdx = Indices::conti0EqIdx };
63 typedef Opm::MathToolbox<Evaluation> Toolbox;
64 typedef Dune::FieldVector<Evaluation, numEq> ParentType;
68 { Opm::Valgrind::SetUndefined(*
this); }
76 template <
class Eval = Evaluation>
77 BlackOilRateVector(
const typename std::enable_if<std::is_same<Eval, Evaluation>::value, Evaluation>::type& value) : ParentType(value)
91 { ParentType::operator=(value); }
99 ParentType::operator=(value);
102 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
103 (*
this)[conti0EqIdx + compIdx] *= FluidSystem::molarMass(compIdx);
109 template <
class Flu
idState,
class RhsEval>
112 const RhsEval& volume)
114 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
115 (*
this)[conti0EqIdx + compIdx] =
116 fluidState.density(phaseIdx)
117 * fluidState.massFraction(phaseIdx, compIdx)
124 template <
class RhsEval>
127 for (
unsigned i=0; i < this->size(); ++i)
137 for (
unsigned i=0; i < this->size(); ++i)
138 (*
this)[i] = other[i];
BlackOilRateVector(const BlackOilRateVector &value)
Default constructor.
Definition: blackoilratevector.hh:84
Contains the quantities which are are constant within a finite volume in the black-oil model...
Implements a vector representing mass, molar or volumetric rates for the black oil model...
Definition: blackoilratevector.hh:50
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:469
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition: blackoilratevector.hh:90
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: blackoilratevector.hh:110
BlackOilRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: blackoilratevector.hh:125
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition: blackoilratevector.hh:96
BlackOilRateVector & operator=(const BlackOilRateVector &other)
Assignment operator from another rate vector.
Definition: blackoilratevector.hh:135
BlackOilRateVector(Scalar value)
Definition: blackoilratevector.hh:73