Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
chemistryModel< CompType, ThermoType > Class Template Reference

Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equation system and evaluation of chemical source terms. More...

Inheritance diagram for chemistryModel< CompType, ThermoType >:
Inheritance graph
[legend]
Collaboration diagram for chemistryModel< CompType, ThermoType >:
Collaboration graph
[legend]

Public Member Functions

 TypeName ("chemistryModel")
 Runtime type information. More...
 
 chemistryModel (const fvMesh &mesh, const word &phaseName)
 Construct from mesh. More...
 
virtual ~chemistryModel ()
 Destructor. More...
 
const PtrList< Reaction< ThermoType > > & reactions () const
 The reactions. More...
 
const PtrList< ThermoType > & specieThermo () const
 Thermodynamic data of the species. More...
 
label nSpecie () const
 The number of species. More...
 
label nReaction () const
 The number of reactions. More...
 
scalar Treact () const
 Temperature below which the reaction rates are assumed 0. More...
 
scalar & Treact ()
 Temperature below which the reaction rates are assumed 0. More...
 
virtual tmp< scalarFieldomega (const scalarField &c, const scalar T, const scalar p) const
 dc/dt = omega, rate of change in concentration, for each species More...
 
virtual scalar omega (const Reaction< ThermoType > &r, const scalarField &c, const scalar T, const scalar p, scalar &pf, scalar &cf, label &lRef, scalar &pr, scalar &cr, label &rRef) const
 Return the reaction rate for reaction r and the reference. More...
 
virtual scalar omegaI (label iReaction, const scalarField &c, const scalar T, const scalar p, scalar &pf, scalar &cf, label &lRef, scalar &pr, scalar &cr, label &rRef) const
 Return the reaction rate for iReaction and the reference. More...
 
virtual void calculate ()
 Calculates the reaction rates. More...
 
const DimensionedField< scalar, volMesh > & RR (const label i) const
 Return const access to the chemical source terms for specie, i. More...
 
virtual DimensionedField< scalar, volMesh > & RR (const label i)
 Return non const access to chemical source terms [kg/m3/s]. More...
 
virtual tmp< DimensionedField< scalar, volMesh > > calculateRR (const label reactionI, const label speciei) const
 Return reaction rate of the speciei in reactionI. More...
 
virtual scalar solve (const scalar deltaT)
 Solve the reaction system for the given time step. More...
 
virtual scalar solve (const scalarField &deltaT)
 Solve the reaction system for the given time step. More...
 
virtual tmp< volScalarFieldtc () const
 Return the chemical time scale. More...
 
virtual tmp< volScalarFieldSh () const
 Return source for enthalpy equation [kg/m/s3]. More...
 
virtual tmp< volScalarFielddQ () const
 Return the heat release, i.e. enthalpy/sec [kg/m2/s3]. More...
 
virtual label nEqns () const
 Number of ODE's to solve. More...
 
virtual void derivatives (const scalar t, const scalarField &c, scalarField &dcdt) const
 Calculate the derivatives in dydx. More...
 
virtual void jacobian (const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
 Calculate the Jacobian of the system. More...
 
virtual void solve (scalarField &c, scalar &T, scalar &p, scalar &deltaT, scalar &subDeltaT) const
 
template<class DeltaTType >
Foam::scalar solve (const DeltaTType &deltaT)
 
- Public Member Functions inherited from ODESystem
 ODESystem ()
 Construct null. More...
 
virtual ~ODESystem ()
 Destructor. More...
 

Protected Types

typedef ThermoType thermoType
 

Protected Member Functions

PtrList< DimensionedField< scalar, volMesh > > & RR ()
 Write access to chemical source terms. More...
 

Protected Attributes

PtrList< volScalarField > & Y_
 Reference to the field of specie mass fractions. More...
 
const PtrList< Reaction< ThermoType > > & reactions_
 Reactions. More...
 
const PtrList< ThermoType > & specieThermo_
 Thermodynamic data of the species. More...
 
label nSpecie_
 Number of species. More...
 
label nReaction_
 Number of reactions. More...
 
scalar Treact_
 Temperature below which the reaction rates are assumed 0. More...
 
PtrList< DimensionedField< scalar, volMesh > > RR_
 List of reaction rate per specie [kg/m3/s]. More...
 

Detailed Description

template<class CompType, class ThermoType>
class Foam::chemistryModel< CompType, ThermoType >

Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equation system and evaluation of chemical source terms.

Source files

Definition at line 59 of file chemistryModel.H.

Member Typedef Documentation

typedef ThermoType thermoType
protected

Definition at line 80 of file chemistryModel.H.

Constructor & Destructor Documentation

chemistryModel ( const fvMesh mesh,
const word phaseName 
)

Construct from mesh.

Definition at line 35 of file chemistryModel.C.

References Foam::dimMass, Foam::dimTime, Foam::dimVolume, Foam::endl(), forAll, Foam::Info, mesh, Foam::name(), fvMesh::time(), and Time::timeName().

Here is the call graph for this function:

~chemistryModel ( )
virtual

Destructor.

Definition at line 88 of file chemistryModel.C.

References chemistryModel< CompType, ThermoType >::omega().

Here is the call graph for this function:

Member Function Documentation

Foam::PtrList< Foam::DimensionedField< Foam::scalar, Foam::volMesh > > & RR ( )
inlineprotected

Write access to chemical source terms.

(e.g. for multi-chemistry model)

Definition at line 33 of file chemistryModelI.H.

Referenced by chemistryModel< CompType, ThermoType >::RR(), and chemistryModel< CompType, ThermoType >::Treact().

Here is the caller graph for this function:

TypeName ( "chemistryModel< CompType, ThermoType >"  )

Runtime type information.

const Foam::PtrList< Foam::Reaction< ThermoType > > & reactions ( ) const
inline

The reactions.

Definition at line 41 of file chemistryModelI.H.

const Foam::PtrList< ThermoType > & specieThermo ( ) const
inline

Thermodynamic data of the species.

Definition at line 49 of file chemistryModelI.H.

Foam::label nSpecie ( ) const
inline

The number of species.

Definition at line 57 of file chemistryModelI.H.

Foam::label nReaction ( ) const
inline

The number of reactions.

Definition at line 65 of file chemistryModelI.H.

Foam::scalar Treact ( ) const
inline

Temperature below which the reaction rates are assumed 0.

Definition at line 73 of file chemistryModelI.H.

Foam::scalar & Treact ( )
inline

Temperature below which the reaction rates are assumed 0.

Definition at line 81 of file chemistryModelI.H.

References chemistryModel< CompType, ThermoType >::RR().

Here is the call graph for this function:

Foam::tmp< Foam::scalarField > omega ( const scalarField c,
const scalar  T,
const scalar  p 
) const
virtual

dc/dt = omega, rate of change in concentration, for each species

Definition at line 97 of file chemistryModel.C.

References forAll, Reaction< ReactionThermo >::lhs(), chemistryModel< CompType, ThermoType >::omegaI(), R, tmp< T >::ref(), Reaction< ReactionThermo >::rhs(), s(), and scalarField().

Referenced by chemistryModel< CompType, ThermoType >::omegaI(), and chemistryModel< CompType, ThermoType >::~chemistryModel().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::scalar omega ( const Reaction< ThermoType > &  r,
const scalarField c,
const scalar  T,
const scalar  p,
scalar &  pf,
scalar &  cf,
label lRef,
scalar &  pr,
scalar &  cr,
label rRef 
) const
virtual

Return the reaction rate for reaction r and the reference.

species and charateristic times

Definition at line 161 of file chemistryModel.C.

References Foam::constant::physicoChemical::c2, chemistryModel< CompType, ThermoType >::derivatives(), Foam::exp(), Reaction< ReactionThermo >::kf(), Reaction< ReactionThermo >::kr(), Reaction< ReactionThermo >::lhs(), Foam::max(), Foam::pow(), Reaction< ReactionThermo >::rhs(), and s().

Here is the call graph for this function:

Foam::scalar omegaI ( label  iReaction,
const scalarField c,
const scalar  T,
const scalar  p,
scalar &  pf,
scalar &  cf,
label lRef,
scalar &  pr,
scalar &  cr,
label rRef 
) const
virtual

Return the reaction rate for iReaction and the reference.

species and charateristic times

Definition at line 139 of file chemistryModel.C.

References chemistryModel< CompType, ThermoType >::omega(), and R.

Referenced by chemistryModel< CompType, ThermoType >::omega().

Here is the call graph for this function:

Here is the caller graph for this function:

void calculate ( )
virtual

Calculates the reaction rates.

Definition at line 712 of file chemistryModel.C.

References Foam::constant::universal::c, forAll, mesh, Foam::constant::mathematical::pi(), rho, thermo, and timeName.

Here is the call graph for this function:

const Foam::DimensionedField< Foam::scalar, Foam::volMesh > & RR ( const label  i) const
inline

Return const access to the chemical source terms for specie, i.

Definition at line 90 of file chemistryModelI.H.

References chemistryModel< CompType, ThermoType >::RR().

Here is the call graph for this function:

Foam::DimensionedField< Foam::scalar, Foam::volMesh > & RR ( const label  i)
virtual

Return non const access to chemical source terms [kg/m3/s].

Definition at line 100 of file chemistryModelI.H.

Foam::tmp< Foam::DimensionedField< Foam::scalar, Foam::volMesh > > calculateRR ( const label  reactionI,
const label  speciei 
) const
virtual

Return reaction rate of the speciei in reactionI.

Definition at line 632 of file chemistryModel.C.

References Foam::constant::universal::c, Foam::dimMass, Foam::dimTime, Foam::dimVolume, forAll, mesh, Foam::constant::mathematical::pi(), tmp< T >::ref(), rho, Foam::constant::thermodynamic::RR, thermo, and timeName.

Referenced by chemistryModel< CompType, ThermoType >::nEqns().

Here is the call graph for this function:

Here is the caller graph for this function:

Foam::scalar solve ( const scalar  deltaT)
virtual

Solve the reaction system for the given time step.

and return the characteristic time

Definition at line 844 of file chemistryModel.C.

References Foam::min(), and solve().

Here is the call graph for this function:

Foam::scalar solve ( const scalarField deltaT)
virtual

Solve the reaction system for the given time step.

and return the characteristic time

Definition at line 859 of file chemistryModel.C.

Foam::tmp< Foam::volScalarField > tc ( ) const
virtual

Return the chemical time scale.

Definition at line 465 of file chemistryModel.C.

References Foam::constant::universal::c, Foam::dimTime, forAll, mesh, Foam::constant::mathematical::pi(), tmp< T >::ref(), rho, Reaction< ReactionThermo >::rhs(), s(), thermo, and timeName.

Here is the call graph for this function:

Foam::tmp< Foam::volScalarField > Sh ( ) const
virtual

Return source for enthalpy equation [kg/m/s3].

Definition at line 551 of file chemistryModel.C.

References Foam::dimEnergy, Foam::dimTime, Foam::dimVolume, forAll, tmp< T >::ref(), and Sh.

Here is the call graph for this function:

Foam::tmp< Foam::volScalarField > dQ ( ) const
virtual

Return the heat release, i.e. enthalpy/sec [kg/m2/s3].

Definition at line 591 of file chemistryModel.C.

References Foam::dimEnergy, Foam::dimTime, dQ, tmp< T >::ref(), GeometricField< Type, PatchField, GeoMesh >::ref(), and Sh.

Here is the call graph for this function:

Foam::label nEqns ( ) const
virtual

Number of ODE's to solve.

Implements ODESystem.

Definition at line 622 of file chemistryModel.C.

References chemistryModel< CompType, ThermoType >::calculateRR().

Here is the call graph for this function:

void derivatives ( const scalar  x,
const scalarField y,
scalarField dydx 
) const
virtual

Calculate the derivatives in dydx.

Implements ODESystem.

Definition at line 279 of file chemistryModel.C.

References cp, chemistryModel< CompType, ThermoType >::jacobian(), and rho.

Referenced by chemistryModel< CompType, ThermoType >::omega().

Here is the call graph for this function:

Here is the caller graph for this function:

void jacobian ( const scalar  x,
const scalarField y,
scalarField dfdx,
scalarSquareMatrix dfdy 
) const
virtual

Calculate the Jacobian of the system.

Need by the stiff-system solvers

Implements ODESystem.

Definition at line 324 of file chemistryModel.C.

References Foam::constant::physicoChemical::c2, delta, forAll, Reaction< ReactionThermo >::kf(), Reaction< ReactionThermo >::kr(), Reaction< ReactionThermo >::lhs(), Foam::max(), nSpecie, Foam::pow(), and Reaction< ReactionThermo >::rhs().

Referenced by chemistryModel< CompType, ThermoType >::derivatives().

Here is the call graph for this function:

Here is the caller graph for this function:

void solve ( scalarField c,
scalar &  T,
scalar &  p,
scalar &  deltaT,
scalar &  subDeltaT 
) const
virtual

Definition at line 869 of file chemistryModel.C.

References NotImplemented.

Foam::scalar solve ( const DeltaTType &  deltaT)

Definition at line 762 of file chemistryModel.C.

References Foam::constant::universal::c, correct, forAll, mesh, Foam::min(), Foam::constant::mathematical::pi(), rho, solve(), thermo, and timeName.

Here is the call graph for this function:

Member Data Documentation

PtrList<volScalarField>& Y_
protected

Reference to the field of specie mass fractions.

Definition at line 85 of file chemistryModel.H.

const PtrList<Reaction<ThermoType> >& reactions_
protected

Reactions.

Definition at line 88 of file chemistryModel.H.

const PtrList<ThermoType>& specieThermo_
protected

Thermodynamic data of the species.

Definition at line 91 of file chemistryModel.H.

label nSpecie_
protected

Number of species.

Definition at line 94 of file chemistryModel.H.

label nReaction_
protected

Number of reactions.

Definition at line 97 of file chemistryModel.H.

scalar Treact_
protected

Temperature below which the reaction rates are assumed 0.

Definition at line 100 of file chemistryModel.H.

PtrList<DimensionedField<scalar, volMesh> > RR_
protected

List of reaction rate per specie [kg/m3/s].

Definition at line 103 of file chemistryModel.H.


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