ABLATE Source Documentation  0.12.33
ablate::eos::EOS Class Referenceabstract

#include <eos.hpp>

+ Inheritance diagram for ablate::eos::EOS:

Public Member Functions

 EOS (std::string typeIn)
 
virtual void View (std::ostream &stream) const =0
 
virtual ThermodynamicFunction GetThermodynamicFunction (ThermodynamicProperty property, const std::vector< domain::Field > &fields) const =0
 
virtual ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction (ThermodynamicProperty property, const std::vector< domain::Field > &fields) const =0
 
virtual EOSFunction GetFieldFunctionFunction (const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const =0
 
virtual const std::vector< std::string > & GetSpeciesVariables () const =0
 
virtual const std::vector< std::string > & GetProgressVariables () const =0
 
virtual std::vector< std::string > GetFieldTags () const
 
virtual std::vector< std::shared_ptr< domain::FieldDescriptor > > GetAdditionalFields () const
 
virtual const std::vector< std::string > & GetFieldFunctionProperties () const
 

Static Public Attributes

static const std::string YI = "Yi"
 
static const std::string PROGRESS = "Progress"
 

Protected Attributes

const std::string type
 

Friends

std::ostream & operator<< (std::ostream &out, const EOS &eos)
 

Detailed Description

The equation of state is designed to to compute thermodynamic properties based upon the conserved field variables being solved. This can range from euler & densityYi and euler & progresses variables with ChemTab.

Member Function Documentation

◆ GetAdditionalFields()

virtual std::vector<std::shared_ptr<domain::FieldDescriptor> > ablate::eos::EOS::GetAdditionalFields ( ) const
inlinevirtual

Helper function that would allow setting up of any additional required fields needed by the eos

Returns

◆ GetFieldFunctionFunction()

virtual EOSFunction ablate::eos::EOS::GetFieldFunctionFunction ( const std::string &  field,
ThermodynamicProperty  property1,
ThermodynamicProperty  property2,
std::vector< std::string >  otherProperties 
) const
pure virtual

Single function to produce fieldFunction function for any two properties, velocity, and species mass fractions. These calls can be slower and should be used for init/output only

Parameters
field
property1
property2

Implemented in ablate::eos::zerorkEOS, ablate::eos::TwoPhase, ablate::eos::TChemSoot, ablate::eos::TChem, ablate::eos::StiffenedGas, ablate::eos::PerfectGas, and ablate::eos::ChemTab.

◆ GetFieldFunctionProperties()

virtual const std::vector<std::string>& ablate::eos::EOS::GetFieldFunctionProperties ( ) const
inlinevirtual

Properties known by this equation of state used for the FieldFunction calculations. This can be the same as the GetSpeciesVariables. species model functions

Returns

Reimplemented in ablate::eos::TwoPhase, and ablate::eos::ChemTab.

◆ GetFieldTags()

virtual std::vector<std::string> ablate::eos::EOS::GetFieldTags ( ) const
inlinevirtual

Returns any additional field tags that can be used when setting up the field

Returns

Reimplemented in ablate::eos::ChemTab.

◆ GetProgressVariables()

virtual const std::vector<std::string>& ablate::eos::EOS::GetProgressVariables ( ) const
pure virtual

Returns a vector of all extra variables required to utilize the equation of state

Returns

Implemented in ablate::eos::zerorkEOS, ablate::eos::TwoPhase, ablate::eos::TChemSoot, ablate::eos::TChemBase, ablate::eos::StiffenedGas, ablate::eos::PerfectGas, and ablate::eos::ChemTab.

◆ GetSpeciesVariables()

virtual const std::vector<std::string>& ablate::eos::EOS::GetSpeciesVariables ( ) const
pure virtual

Required species to utilize the equation of state species model functions

Returns

Implemented in ablate::eos::zerorkEOS, ablate::eos::TwoPhase, ablate::eos::TChemBase, ablate::eos::StiffenedGas, ablate::eos::PerfectGas, and ablate::eos::ChemTab.

◆ GetThermodynamicFunction()

virtual ThermodynamicFunction ablate::eos::EOS::GetThermodynamicFunction ( ThermodynamicProperty  property,
const std::vector< domain::Field > &  fields 
) const
pure virtual

Single function to produce thermodynamic function for any property based upon the available fields

Parameters
property
fields
Returns

Implemented in ablate::eos::zerorkEOS, ablate::eos::TwoPhase, ablate::eos::TChemSoot, ablate::eos::TChem, ablate::eos::StiffenedGas, ablate::eos::PerfectGas, and ablate::eos::ChemTab.

◆ GetThermodynamicTemperatureFunction()

virtual ThermodynamicTemperatureFunction ablate::eos::EOS::GetThermodynamicTemperatureFunction ( ThermodynamicProperty  property,
const std::vector< domain::Field > &  fields 
) const
pure virtual

Single function to produce thermodynamic function for any property based upon the available fields and temperature

Parameters
property
fields
Returns

Implemented in ablate::eos::zerorkEOS, ablate::eos::TwoPhase, ablate::eos::TChemSoot, ablate::eos::TChem, ablate::eos::StiffenedGas, ablate::eos::PerfectGas, and ablate::eos::ChemTab.

◆ View()

virtual void ablate::eos::EOS::View ( std::ostream &  stream) const
pure virtual

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const EOS eos 
)
friend

Support function for printing any eos

Parameters
out
eos
Returns

The documentation for this class was generated from the following file: