28#ifndef EWOMS_FLASH_RATE_VECTOR_HH
29#define EWOMS_FLASH_RATE_VECTOR_HH
31#include <dune/common/fvector.hh>
33#include <opm/material/common/Valgrind.hpp>
44template <
class TypeTag>
46 :
public Dune::FieldVector<GetPropType<TypeTag, Properties::Evaluation>,
47 getPropValue<TypeTag, Properties::NumEq>()>
54 enum { conti0EqIdx = Indices::conti0EqIdx };
58 using ParentType = Dune::FieldVector<Evaluation, numEq>;
64 { Valgrind::SetUndefined(*
this); }
87 ParentType molarRate(value);
88 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
89 molarRate[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
99 { ParentType::operator=(value); }
105 { EnergyModule::setEnthalpyRate(*
this, rate); }
110 template <
class Flu
idState,
class RhsEval>
113 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
114 (*this)[conti0EqIdx + compIdx] =
115 fluidState.density(phaseIdx, compIdx) *
116 fluidState.moleFraction(phaseIdx, compIdx) *
120 EnergyModule::setEnthalpyRate(*
this, fluidState, phaseIdx, volume);
126 template <
class RhsEval>
129 for (
unsigned i = 0; i < this->size(); ++i) {
138 FlashRateVector&
operator=(
const FlashRateVector& other)
140 for (
unsigned i = 0; i < this->size(); ++i) {
141 (*this)[i] = other[i];
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:54
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
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240