28 #ifndef EWOMS_FLASH_RATE_VECTOR_HH
29 #define EWOMS_FLASH_RATE_VECTOR_HH
31 #include <dune/common/fvector.hh>
34 #include <opm/material/constraintsolvers/NcpFlash.hpp>
35 #include <opm/common/Valgrind.hpp>
46 template <
class TypeTag>
48 :
public Dune::FieldVector<typename GET_PROP_TYPE(TypeTag, Evaluation),
49 GET_PROP_VALUE(TypeTag, NumEq)>
51 typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
52 typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
53 typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
54 typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
56 enum { conti0EqIdx = Indices::conti0EqIdx };
60 typedef Dune::FieldVector<Evaluation, numEq> ParentType;
65 { Opm::Valgrind::SetUndefined(*
this); }
86 ParentType molarRate(value);
87 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
88 molarRate[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
97 { ParentType::operator=(value); }
103 { EnergyModule::setEnthalpyRate(*
this, rate); }
108 template <
class Flu
idState,
class RhsEval>
111 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
112 (*
this)[conti0EqIdx + compIdx] =
113 fluidState.density(phaseIdx, compIdx)
114 * fluidState.moleFraction(phaseIdx, compIdx)
117 EnergyModule::setEnthalpyRate(*
this, fluidState, phaseIdx, volume);
123 template <
class RhsEval>
126 for (
unsigned i=0; i < this->size(); ++i)
136 for (
unsigned i=0; i < this->size(); ++i)
137 (*
this)[i] = other[i];
void setEnthalpyRate(const Evaluation &rate)
Set an enthalpy rate [J/As] where .
Definition: flashratevector.hh:102
FlashRateVector(const FlashRateVector &value)
Default constructor.
Definition: flashratevector.hh:77
Contains the intensive quantities of the flash-based compositional multi-phase model.
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition: flashratevector.hh:109
FlashRateVector & operator=(const FlashRateVector &other)
Assignment operator from another rate vector.
Definition: flashratevector.hh:134
Provides the auxiliary methods required for consideration of the energy equation. ...
Definition: energymodule.hh:59
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:469
Contains the classes required to consider energy as a conservation quantity in a multi-phase module...
Implements a vector representing rates of conserved quantities.
Definition: flashratevector.hh:47
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition: flashratevector.hh:83
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition: flashratevector.hh:96
FlashRateVector(const Evaluation &value)
Definition: flashratevector.hh:70
FlashRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition: flashratevector.hh:124