28#ifndef EWOMS_FLASH_RATE_VECTOR_HH
29#define EWOMS_FLASH_RATE_VECTOR_HH
31#include <dune/common/fvector.hh>
34#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>;
63 { Opm::Valgrind::SetUndefined(*
this); }
84 ParentType molarRate(value);
85 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
86 molarRate[conti0EqIdx + compIdx] /= FluidSystem::molarMass(compIdx);
95 { ParentType::operator=(value); }
101 { EnergyModule::setEnthalpyRate(*
this, rate); }
106 template <
class Flu
idState,
class RhsEval>
109 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
110 (*
this)[conti0EqIdx + compIdx] =
111 fluidState.density(phaseIdx, compIdx)
112 * fluidState.moleFraction(phaseIdx, compIdx)
115 EnergyModule::setEnthalpyRate(*
this, fluidState, phaseIdx, volume);
121 template <
class RhsEval>
124 for (
unsigned i=0; i < this->size(); ++i)
134 for (
unsigned i=0; i < this->size(); ++i)
135 (*
this)[i] = other[i];
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:50
Implements a vector representing rates of conserved quantities.
Definition flashratevector.hh:48
void setEnthalpyRate(const Evaluation &rate)
Set an enthalpy rate [J/As] where .
Definition flashratevector.hh:100
FlashRateVector & operator=(const FlashRateVector &other)
Assignment operator from another rate vector.
Definition flashratevector.hh:132
void setMolarRate(const ParentType &value)
Set a molar rate of the conservation quantities.
Definition flashratevector.hh:94
void setVolumetricRate(const FluidState &fluidState, unsigned phaseIdx, const RhsEval &volume)
Set a volumetric rate of a phase.
Definition flashratevector.hh:107
FlashRateVector(const Evaluation &value)
Definition flashratevector.hh:68
FlashRateVector(const FlashRateVector &value)
Definition flashratevector.hh:75
void setMassRate(const ParentType &value)
Set a mass rate of the conservation quantities.
Definition flashratevector.hh:81
FlashRateVector & operator=(const RhsEval &value)
Assignment operator from a scalar or a function evaluation.
Definition flashratevector.hh:122
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 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:235
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242