ChemFoam

From OpenFOAMWiki

chemFoam

   Solver for chemistry problems, designed for use on single cell cases to
   provide comparison against other chemistry solvers, that uses a single cell
   mesh, and fields created from the initial conditions.

1 Solution Strategy

This solver provides an excellent starting point for those who want to have a first impression on the influence of the chemical reactions in the equations describing the evolution of the species concentration and also the evolution of temperature (i.e. the energy equation). Since the computational domain consists only of one cell, the only mechanism influencing the evolution of the spices concentration and of the temperature are the chemical reactions. The system of equations which is solved is the following:



    \frac{  \partial \rho Y_i}{\partial t}   =  \dot \omega_i (Y_i,T)
(1)

     h = u_0 + \frac{p}{\rho} + \int_0^t \frac{\dot q}{\rho} d\tau
(2)

     p = \frac{ \rho R T}{M_{avg}}
(3)

 Y_i, \omega_i, T, p, \rho, \dot q, R, M_{avg}, u_0, are the species mass mass fraction, the reaction rate, the temperature, the pressure, the density, the heat generated by the combustion, the gas constant, the average molar weight and the initial internal energy, respectively. Note that for a given enthalpy, pressure and density the temperature can be calculated.

Like most of the solvers present in OpenFOAM also this solver follows a segregated solution strategy. That means that for each quantity of interest one linear equation is solved and the coupling between the equations is achieved by explicit source terms. Briefly summarized the solution is achieved as follows:

1) Solve the chemistry: The purpose is to get the reaction rates for each species involved and the heat realized by the chemical reaction

2) Solve the species equation: The purpose it to get the species concentration at the new time step

3) Solve the energy equation: Here we get the enthalpy (over the thermodynamics also the temperature) at the new time step

4) Solve the pressure equation: Required to calculate the enthalpy


The source code can be found in chemFoam.C


 
 
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2010-2011 OpenCFD Ltd.
     \\/     M anipulation  |
-------------------------------------------------------------------------------
                            | Copyright (C) 2011-2017 OpenFOAM Foundation
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.
 
    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.
 
    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.
 
    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
 
Application
    chemFoam
 
Group
    grpCombustionSolvers
 
Description
    Solver for chemistry problems, designed for use on single cell cases to
    provide comparison against other chemistry solvers, that uses a single cell
    mesh, and fields created from the initial conditions.
 
\*---------------------------------------------------------------------------*/
 
#include "fvCFD.H"
#include "psiReactionThermo.H"
#include "BasicChemistryModel.H"
#include "reactingMixture.H"
#include "chemistrySolver.H"
#include "OFstream.H"
#include "thermoPhysicsTypes.H"
#include "basicSpecieMixture.H"
#include "hexCellFvMesh.H"
 
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
int main(int argc, char *argv[])
{
    argList::addNote
    (
        "Solver for chemistry problems, designed for use on single cell cases"
        " to provide comparison against other chemistry solvers"
    );
 
    argList::noParallel();
 
    #define CREATE_MESH createSingleCellMesh.H
    #define NO_CONTROL
    #include "postProcess.H"
 
    #include "setRootCaseLists.H"
    #include "createTime.H"
    #include "createSingleCellMesh.H"
    #include "createFields.H"
    #include "createFieldRefs.H"
    #include "readInitialConditions.H"
    #include "createControls.H"
 
    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
    Info<< "\nStarting time loop\n" << endl;
 
    while (runTime.run())
    {
        #include "readControls.H"
 
        #include "setDeltaT.H"
 
        ++runTime;
        Info<< "Time = " << runTime.timeName() << nl << endl;
 
        #include "solveChemistry.H"
        #include "YEqn.H"
        #include "hEqn.H"
        #include "pEqn.H"
 
        #include "output.H"
 
        runTime.printExecutionTime(Info);
    }
 
    Info << "Number of steps = " << runTime.timeIndex() << endl;
    Info << "End" << nl << endl;
 
    return 0;
}
 
 
// ************************************************************************* //
 
 

2 Chemical reactions

2.1 Chemical equations

Elementary chemical reaction involving  K species in  I reactions can be expressed in the following form (see also https://www3.nd.edu/~powers/ame.60636/chemkin2000.pdf or the book of [1]):


  \sum_{k=1}^K  \nu_{ki}' \chi_k \rightleftharpoons \sum_{k=1}^K  \nu_{ki}'' \chi_k  \ \ \ i=1...I
(4)

  \nu_{ki}',  \nu_{ki}'' and   \chi_k are the forward, the backward stoichiometric coefficient and the chemical symbol of the specie   k , respectively.

In order to make the notation more clear to those who start learning about the usage of CFD to simulate combustion process, Zel'dovich mechanism for the oxidation of N0 (see https://en.wikipedia.org/wiki/Zeldovich_mechanism or also the book of [2]) is used as an example:


  N_2 + O \rightleftharpoons  NO + N

N + O_2 \rightleftharpoons NO + O

Based on the above example the forward   \nu_{ki}' , the backward stoichiometric coefficient    \nu_{ki}'' and the chemical symbol of the specie k   \chi_k can be written as follows:


 \nu_{ki}' = \begin{pmatrix}1 & 1 & 0 & 0 & 0 \\0 & 0 & 0 & 1 & 1  \\\end{pmatrix},   \nu_{ki}'' = \begin{pmatrix}0 & 0 & 1 & 1 & 0 \\0 & 1 & 1 & 0 & 0  \\\end{pmatrix} ,   \chi_k = \begin{pmatrix}N_2 \\ O \\ NO \\ N \\ O_2  \\\end{pmatrix}

The reactions are specified in the file constant/thermophysicalProperties where the chemistry reader, the file containing the reactions and the thermodynamic properties have to be specified:


 
/*--------------------------------*- C++ -*----------------------------------*\
| =========                 |                                                 |
| \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox           |
|  \\    /   O peration     | Version:  v1912                                 |
|   \\  /    A nd           | Website:  www.openfoam.com                      |
|    \\/     M anipulation  |                                                 |
\*---------------------------------------------------------------------------*/
FoamFile
{
    version     2.0;
    format      ascii;
    class       dictionary;
    location    "constant";
    object      thermophysicalProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
thermoType
{
    type            hePsiThermo;
    mixture         reactingMixture;
    transport       sutherland;
    thermo          janaf;
    energy          sensibleEnthalpy;
    equationOfState perfectGas;
    specie          specie;
}
 
 
chemistryReader foamChemistryReader;
foamChemistryFile "<constant>/reactions";
foamChemistryThermoFile "<constant>/thermo";
 
//CHEMKINFile         "<case>/chemkin/chem.inp";
//CHEMKINThermoFile   "<case>/chemkin/therm.dat";
//CHEMKINTransportFile "<case>/chemkin/transportProperties";
 
// ************************************************************************* //

For the above example the reactions look as follows:


 
elements        2(O N);
 
species         5(O O2 N2 N NO);
 
reactions
{
    un-named-reaction-0
    {
        type            reversibleArrheniusReaction;
        reaction        "NO^1 + N^1 = N2^1 + O^1";
        A               2.1e+10;
        beta            0;
        Ta              0;
    }
    un-named-reaction-1
    {
        type            reversibleArrheniusReaction;
        reaction        "N^1 + O2^1 = NO^1 + O^1";
        A               5.83e+06;
        beta            1.01;
        Ta              3120;
    }
}

In our example both reactions are reversible Arrhenius reactions. A is the per-exponent factor, Ta the activation temperature and beta is the exponent of the temperature in the calculation of the reaction rate (see below for further explanations). The values are taken from the book of [3]. Note that the values of the constants are calculated in the OpenFOAM native unit system of kmol, m3, s, K.

2.2 Chemical reaction rate

Having established the chemical reactions which are desired to be solved, the next step is to describe the velocity at which the chemical reaction occur, i.e. the change of the concentration of the single species with time. The reaction rate  \omega_i of the reaction i can be written as (see also the book [4]):


\omega_i = k_{fi} \prod_k [X]_k^{\nu_{ki}'} - k_{ri} \prod_k [X]_k^{\nu_{ki}''}

 k_{fi} and  k_{ri} are the forward and reverse rate constant reaction i. Note that in OpenFoam the default value of the exponents in the above expression are the stoichiometric coefficients   \nu_{ki}',  \nu_{ki}'' are used. However the user can explicitly set it by adding the hat sign ^ after the chemical symbol (see also https://openfoamwiki.net/images/8/8a/Christ_OFW6_slides.pdf).   [X]_k is the molar concentration of the species k. To continue with our example the reaction rates of the Zeldovich mechanism can be written:


 \omega_1 = k_{f1}[NO][N] - k_{r1}[N2][O]

\omega_2 = k_{f2}[N][O2] - k_{r2}[NO][O]

The calculation of the reaction rates of each reaction is done in the function omega of the file $FOAM_SRC/thermophysicalModels/chemistryModel/chemistryModel/StandardChemistryModel/StandardChemistryModel.C:


 
template<class ReactionThermo, class ThermoType>
Foam::scalar Foam::StandardChemistryModel<ReactionThermo, ThermoType>::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
{
    const scalar kf = R.kf(p, T, c);
    const scalar kr = R.kr(kf, p, T, c);
 
    pf = 1.0;
    pr = 1.0;
 
    const label Nl = R.lhs().size();
    const label Nr = R.rhs().size();
 
    label slRef = 0;
    lRef = R.lhs()[slRef].index;
 
    pf = kf;
    for (label s = 1; s < Nl; s++)
    {
        const label si = R.lhs()[s].index;
 
        if (c[si] < c[lRef])
        {
            const scalar exp = R.lhs()[slRef].exponent;
            pf *= pow(max(c[lRef], 0.0), exp);
            lRef = si;
            slRef = s;
        }
        else
        {
            const scalar exp = R.lhs()[s].exponent;
            pf *= pow(max(c[si], 0.0), exp);
        }
    }
    cf = max(c[lRef], 0.0);
 
    {
        const scalar exp = R.lhs()[slRef].exponent;
        if (exp < 1.0)
        {
            if (cf > SMALL)
            {
                pf *= pow(cf, exp - 1.0);
            }
            else
            {
                pf = 0.0;
            }
        }
        else
        {
            pf *= pow(cf, exp - 1.0);
        }
    }
 
    label srRef = 0;
    rRef = R.rhs()[srRef].index;
 
    // Find the matrix element and element position for the rhs
    pr = kr;
    for (label s = 1; s < Nr; s++)
    {
        const label si = R.rhs()[s].index;
        if (c[si] < c[rRef])
        {
            const scalar exp = R.rhs()[srRef].exponent;
            pr *= pow(max(c[rRef], 0.0), exp);
            rRef = si;
            srRef = s;
        }
        else
        {
            const scalar exp = R.rhs()[s].exponent;
            pr *= pow(max(c[si], 0.0), exp);
        }
    }
    cr = max(c[rRef], 0.0);
 
    {
        const scalar exp = R.rhs()[srRef].exponent;
        if (exp < 1.0)
        {
            if (cr>SMALL)
            {
                pr *= pow(cr, exp - 1.0);
            }
            else
            {
                pr = 0.0;
            }
        }
        else
        {
            pr *= pow(cr, exp - 1.0);
        }
    }
 
    return pf*cf - pr*cr;
}


In order to get the reaction rate of the spices k, the reaction rates of each reaction containing the species k has to be multiplied with the stoichiometric coefficient and summed with a positive sign if the species is a product in with a negative sign if the species is a reactant (see also e.g. https://en.wikipedia.org/wiki/Reaction_rate). It can be written as (see also https://www3.nd.edu/~powers/ame.60636/chemkin2000.pdf and the books cited here):


\omega_k = \sum_i \omega_i ({\nu_{ki}'} - {\nu_{ki}''})

For our example reaction we get:


 \omega_k = \frac{d [X_k]}{dt} =  \begin{pmatrix}\frac{d[N_2]}{dt} \\ \frac{d[O]}{dt} \\ \frac{d[NO]}{dt} \\ \frac{d[N]}{dt} \\ \frac{d[O_2]}{dt}    \\\end{pmatrix} = \begin{pmatrix} -\omega_1 \\ -\omega_1  + \omega_2   \\   \omega_1  + \omega_2     \\  \omega_1 - \omega_2       \\  \omega_2    \\\end{pmatrix} =  \begin{pmatrix} -(k_{f1}[NO][N] - k_{r1}[N2][O]) \\ -(k_{f1}[NO][N] - k_{r1}[N2][O])  + k_{f2}[N][O2] - k_{r2}[NO][O]  \\   k_{f1}[NO][N] - k_{r1}[N2][O]  + k_{f2}[N][O2] - k_{r2}[NO][O]     \\  k_{f1}[NO][N] - k_{r1}[N2][O] - (k_{f2}[N][O2] - k_{r2}[NO][O])       \\  k_{f2}[N][O2] - k_{r2}[NO][O]    \\\end{pmatrix}

Looking at the above equation we see that the equation describing the evolution of the species concentration is a coupled non--linear differential equation. The calculation of the reaction rates of each species is done in the function of the file $FOAM_SRC/thermophysicalModels/chemistryModel/chemistryModel/StandardChemistryModel/StandardChemistryModel.C:


 
template<class ReactionThermo, class ThermoType>
void Foam::StandardChemistryModel<ReactionThermo, ThermoType>::omega
(
    const scalarField& c,
    const scalar T,
    const scalar p,
    scalarField& dcdt
) const
{
    scalar pf, cf, pr, cr;
    label lRef, rRef;
 
    dcdt = Zero;
 
    forAll(reactions_, i)
    {
        const Reaction<ThermoType>& R = reactions_[i];
 
        scalar omegai = omega
        (
            R, c, T, p, pf, cf, lRef, pr, cr, rRef
        );
 
        forAll(R.lhs(), s)
        {
            const label si = R.lhs()[s].index;
            const scalar sl = R.lhs()[s].stoichCoeff;
            dcdt[si] -= sl*omegai;
        }
 
        forAll(R.rhs(), s)
        {
            const label si = R.rhs()[s].index;
            const scalar sr = R.rhs()[s].stoichCoeff;
            dcdt[si] += sr*omegai;
        }
    }
}

2.3 Chemical reaction rate constants

The remaining two quantities to define are the forward and backward reaction rate constants  k_f and  k_r .

2.3.1 Forward rate constants

The forward rate constant is calculated for Arrhenius type reactions as follows (see also [5]):


k_f = A T^\beta exp(-T_a/T)

A is the per-exponent factor, Ta the activation temperature and beta is the exponent of the temperature in the calculation of the reaction rate. The soure code can be find in $FOAM_SRC/thermophysicalModels/specie/reaction/reactionRate/ArrheniusReactionRate/ArrheniusReactionRateI.H:


 
inline Foam::scalar Foam::ArrheniusReactionRate::operator()
(
    const scalar p,
    const scalar T,
    const scalarField&
) const
{
    scalar ak = A_;
 
    if (mag(beta_) > VSMALL)
    {
        ak *= pow(T, beta_);
    }
 
    if (mag(Ta_) > VSMALL)
    {
        ak *= exp(-Ta_/T);
    }
 
    return ak;
}

2.3.2 Reverse rate constants

The reverse rate of the ith reaction is calculated as follows (see also https://www3.nd.edu/~powers/ame.60636/chemkin2000.pdf) :


k_{ri} = \frac{k_{fi}}{K_{ci}}

The equilibrium constant is calculated as follows:


{K_{ci}} = K_{pi} (\frac{p_{atm}}{R T})^{\sum_k \nu_{ki}}

\sum \nu_{ki} = \sum \nu_{ki}' -  \sum \nu_{ki}'' is the sum of the stoichiometric coefficient of the ith reaction. Furthermore



 K_{pi} = exp \left( \sum_k \nu_{ki}( \frac{S^0_k}{R} -  \frac{H^0_k}{RT}) \right)

 H^0_k and  S^0_k are the molar enthalpy and entropy of the species k .

The reverse rate constant is calculated in the file $FOAM_SRC/thermophysicalModels/specie/reaction/Reactions/ReversibleReaction/ReversibleReaction.C:


 
template
<
    template<class> class ReactionType,
    class ReactionThermo,
    class ReactionRate
>
Foam::scalar Foam::ReversibleReaction
<
    ReactionType,
    ReactionThermo,
    ReactionRate
>::kr
(
    const scalar kfwd,
    const scalar p,
    const scalar T,
    const scalarField& c
) const
{
    return kfwd/max(this->Kc(p, T), 1e-6);
}

The calculation of Kc can be found in the file $FOAM_SRC/thermophysicalModels/specie/thermo/thermo/thermoI.H.


 
template<class Thermo, template<class> class Type>
inline Foam::scalar
Foam::species::thermo<Thermo, Type>::K(const scalar p, const scalar T) const
{
    scalar arg = -this->Y()*this->G(Pstd, T)/(RR*T);
 
    if (arg < 600)
    {
        return exp(arg);
    }
    else
    {
        return VGREAT;
    }
}
 
 
template<class Thermo, template<class> class Type>
inline Foam::scalar
Foam::species::thermo<Thermo, Type>::Kp(const scalar p, const scalar T) const
{
    return K(p, T);
}
 
 
template<class Thermo, template<class> class Type>
inline Foam::scalar
Foam::species::thermo<Thermo, Type>::Kc(const scalar p, const scalar T) const
{
    const scalar nm = this->Y()/this->W();
 
    if (equal(nm, SMALL))
    {
        return Kp(p, T);
    }
    else
    {
        return Kp(p, T)*pow(Pstd/(RR*T), nm);
    }
}

It is not so easily seen how the sum of the stoichiometric coefficients enters the calculation of the equilibrium constant. For this reason it is explained in a bit more detail (see also the thread https://www.cfd-online.com/Forums/openfoam/224970-understanding-reverse-rate-calculation-chemical-reactions.html). In order to calculate the equilibrium constant for each reaction a pseudo mixture is set where the thermodynamic properties of the mixture are weighted with the stoichiometric coefficient. This is done in the constructor of the Reaction object by calling the function setThermo:


 
template<class ReactionThermo>
Foam::Reaction<ReactionThermo>::Reaction
(
    const speciesTable& species,
    const List<specieCoeffs>& lhs,
    const List<specieCoeffs>& rhs,
    const HashPtrTable<ReactionThermo>& thermoDatabase,
    bool initReactionThermo
)
:
    ReactionThermo::thermoType(*thermoDatabase[species[0]]),
    name_("un-named-reaction-" + Foam::name(getNewReactionID())),
    species_(species),
    lhs_(lhs),
    rhs_(rhs)
{
    if (initReactionThermo)
    {
        setThermo(thermoDatabase);
    }
}
 
template<class ReactionThermo>
void Foam::Reaction<ReactionThermo>::setThermo
(
    const HashPtrTable<ReactionThermo>& thermoDatabase
)
{
 
    typename ReactionThermo::thermoType rhsThermo
    (
        rhs_[0].stoichCoeff
        *(*thermoDatabase[species_[rhs_[0].index]]).W()
        *(*thermoDatabase[species_[rhs_[0].index]])
    );
 
    for (label i=1; i<rhs_.size(); ++i)
    {
        rhsThermo +=
            rhs_[i].stoichCoeff
        *(*thermoDatabase[species_[rhs_[i].index]]).W()
        *(*thermoDatabase[species_[rhs_[i].index]]);
    }
 
    typename ReactionThermo::thermoType lhsThermo
    (
        lhs_[0].stoichCoeff
       *(*thermoDatabase[species_[lhs_[0].index]]).W()
       *(*thermoDatabase[species_[lhs_[0].index]])
    );
 
    for (label i=1; i<lhs_.size(); ++i)
    {
        lhsThermo +=
            lhs_[i].stoichCoeff
           *(*thermoDatabase[species_[lhs_[i].index]]).W()
           *(*thermoDatabase[species_[lhs_[i].index]]);
    }
 
    ReactionThermo::thermoType::operator=(lhsThermo == rhsThermo);
}

To undestand how the equilibrium constant Kc is calcualted we have to have a look at the operators +=, * and == defined in the thermodynamic packages. In the above listing a spieces is constructed where the thermodynimc properties are first multiplied with the stoichiometric coefficient and the molare weight, then the sum of the spieces at the right and left hand side of the chemical reaction is performed and as last step the == operator is applied. As an illustrative example to explain this operations we take the constant properties termodynamic package. The conclusions drawn here can be applied to all other packages as well.


 
template<class EquationOfState>
inline Foam::hConstThermo<EquationOfState> Foam::operator*
(
    const scalar s,
    const hConstThermo<EquationOfState>& ct
)
{
    return hConstThermo<EquationOfState>
    (
        s*static_cast<const EquationOfState&>(ct),
        ct.Cp_,
        ct.Hf_
    );
}
 
template<class EquationOfState>
inline void Foam::hConstThermo<EquationOfState>::operator+=
(
    const hConstThermo<EquationOfState>& ct
)
{
    scalar Y1 = this->Y();
 
    EquationOfState::operator+=(ct);
 
    if (mag(this->Y()) > SMALL)
    {
        Y1 /= this->Y();
        scalar Y2 = ct.Y()/this->Y();
 
        Cp_ = Y1*Cp_ + Y2*ct.Cp_;
        Hf_ = Y1*Hf_ + Y2*ct.Hf_;
    }
}
 
template<class EquationOfState>
inline Foam::hConstThermo<EquationOfState> Foam::operator==
(
    const hConstThermo<EquationOfState>& ct1,
    const hConstThermo<EquationOfState>& ct2
)
{
    EquationOfState eofs
    (
        static_cast<const EquationOfState&>(ct1)
     == static_cast<const EquationOfState&>(ct2)
    );
 
    return hConstThermo<EquationOfState>
    (
        eofs,
        ct2.Y()/eofs.Y()*ct2.Cp_
      - ct1.Y()/eofs.Y()*ct1.Cp_,
        ct2.Y()/eofs.Y()*ct2.Hf_
      - ct1.Y()/eofs.Y()*ct1.Hf_
    );
}

The first operator analized the * operator. When a hConstThermo object is muliplied with a scalar s the operator returns a new object very only the equationOfState attribute is muliplied by this scalar. For this reason we have to look at the operator * in the equationOfState object. As an axemple we look at the perfectGas equation of state:


 
template<class Specie>
inline Foam::perfectGas<Specie> Foam::operator*
(
    const scalar s,
    const perfectGas<Specie>& pg
)
{
    return perfectGas<Specie>(s*static_cast<const Specie&>(pg));
}

It returns a perfectGas object were the Specie attribut is muliplied with this scalar. For the Specie attribut the * operator returns a specie with the Y_ attribute multiplied by the scalar s:


 
inline specie operator*(const scalar s, const specie& st)
{
    return specie
    (
        s*st.Y_,
        st.molWeight_
    );
}

So the application of the operator leads only to a change of the Y_ attribute in both the left and right hand side of the chemical reaction:


Y' = \nu' W'

Y'' = \nu'' W''

The next operator we're going to examine is the += operator. We see that first the Y_ attribute is stored to the variable Y1. Than the += operator of the equation of state is executed. After that the weighted sum of the thermodynamic properties of the two involved species is taken.


 
template<class Specie>
inline void Foam::perfectGas<Specie>::operator+=(const perfectGas<Specie>& pg)
{
    Specie::operator+=(pg);
}



 
inline void specie::operator+=(const specie& st)
{
    const scalar sumY = Y_ + st.Y_;
    if (mag(sumY) > SMALL)
    {
        molWeight_ = sumY/(Y_/molWeight_ + st.Y_/st.molWeight_);
    }
 
    Y_ = sumY;
}
 
 

It's easyly seen, that for the attributes Y and molWeight we obtain for the left hand side and right hand side of the chemical equation:


Y' = \sum_k \nu_k' W_k'

Y'' = \sum_k \nu_k'' W_k''

molWeight' = \frac{\sum_k \nu_k' W_k'}{\sum_k \nu_k'}

molWeight'' =  \frac{\sum_k \nu_k'' W_k''}{\sum_k \nu_k''}

For the thermodynamic properties Hf (chemical enthalpy) and the specific heat capacity we get:


Hf' =  \frac{\sum_k \nu_k' W_k'Hf'_k}{\sum_k \nu_k'W_k'}

Hf'' =  \frac{\sum_k \nu_k'' W_k''Hf''_k}{\sum_k \nu_k''W_k''}

cp' =  \frac{\sum_k \nu_k' W_k'cp'_k}{\sum_k \nu_k'W_k'}

cp'' =  \frac{\sum_k \nu_k'' W_k''cp''_k}{\sum_k \nu_k''W_k''}

The final operator we have a look is the == operator. Also here first the equation of state == operator is executed. Also as for the other operators the equation of state == operator executes the == of the species:


 
template<class Specie>
inline Foam::perfectGas<Specie> Foam::operator==
(
    const perfectGas<Specie>& pg1,
    const perfectGas<Specie>& pg2
)
{
    return perfectGas<Specie>
    (
        static_cast<const Specie&>(pg1) == static_cast<const Specie&>(pg2)
    );
}


 
inline specie operator==(const specie& st1, const specie& st2)
{
    scalar diffY = st2.Y_ - st1.Y_;
    if (mag(diffY) < SMALL)
    {
        diffY = SMALL;
    }
 
    const scalar diffRW = st2.Y_/st2.molWeight_ - st1.Y_/st1.molWeight_;
 
    #ifdef __clang__
    // Using intermediate volatile bool to prevent compiler optimising out the
    // if block (above) - CLANG 3.7.1
    volatile const bool valid = (mag(diffRW) > SMALL);
    const scalar molWeight = valid ? diffY/diffRW : GREAT;
    #else
    scalar molWeight = GREAT;
    if (mag(diffRW) > SMALL)
    {
         molWeight = diffY/diffRW;
    }
    #endif
 
    return specie(diffY, molWeight);
}


So after the application of the == we get for the Y and molWeight attribut of the specie:


Y = \sum_k \nu_k' W_k' - \sum_k \nu_k'' W_k''

molWeight' = \frac{\sum_k \nu_k' W_k' - \sum_k \nu_k'' W_k''}{\sum_k \nu_k' - \sum_k \nu_k''}

For the thermodynamic properties we get:


cp = \frac{\sum_k \nu_k' W_k'cp_k' - \sum_k \nu_k'' W_k''cp_k''}{\sum_k \nu_k' W_k' - \sum_k \nu_k'' W_k'' }

Hf = \frac{\sum_k \nu_k' W_k'Hf_k' - \sum_k \nu_k'' W_k''Hf_k''}{\sum_k \nu_k' W_k' - \sum_k \nu_k'' W_k'' }

Finnaly if we look how the enthalpy and the entropy are computed for the combination hconst and perfectGas we can veryfy that the same definition for the equilibrium constant  K_{ci} as for chemkin is obtained:


 
template<class Thermo, template<class> class Type>
inline Foam::scalar
Foam::species::thermo<Thermo, Type>::G(const scalar p, const scalar T) const
{
    return this->Ha(p, T) - T*this->S(p, T);
}


 
template<class EquationOfState>
inline Foam::scalar Foam::hConstThermo<EquationOfState>::Ha
(
    const scalar p, const scalar T
) const
{
    return Hs(p, T) + Hc();
}
 
 
template<class EquationOfState>
inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hs
(
    const scalar p, const scalar T
) const
{
    return Cp_*T + EquationOfState::H(p, T);
}
 
 
template<class EquationOfState>
inline Foam::scalar Foam::hConstThermo<EquationOfState>::Hc() const
{
    return Hf_;
}


 
template<class Specie>
inline Foam::scalar Foam::perfectGas<Specie>::S(scalar p, scalar T) const
{
    return -this->R()*log(p/Pstd);
}


 
template<class Specie>
inline Foam::scalar Foam::perfectGas<Specie>::H(scalar p, scalar T) const
{
    return 0;
}

2.4 Heat released by the chemical reactions

The heat released by the chemical reactions is returned by the function Qdot() in the file $FOAM_SRC/thermophysicalModels/chemistryModel/chemistryModel/StandardChemistryModel/StandardChemistryModel.C (if the standard chemistry model is chosen):


 
template<class ReactionThermo, class ThermoType>
Foam::tmp<Foam::volScalarField>
Foam::StandardChemistryModel<ReactionThermo, ThermoType>::Qdot() const
{
    tmp<volScalarField> tQdot
    (
        new volScalarField
        (
            IOobject
            (
                "Qdot",
                this->mesh_.time().timeName(),
                this->mesh_,
                IOobject::NO_READ,
                IOobject::NO_WRITE,
                false
            ),
            this->mesh_,
            dimensionedScalar(dimEnergy/dimVolume/dimTime, Zero)
        )
    );
 
    if (this->chemistry_)
    {
        scalarField& Qdot = tQdot.ref();
 
        forAll(Y_, i)
        {
            forAll(Qdot, celli)
            {
                const scalar hi = specieThermo_[i].Hc();
                Qdot[celli] -= hi*RR_[i][celli];
            }
        }
    }
 
    return tQdot;
}

Hc is the chemical enthalpy and depends on the thermodynamic model chosen.

3 Evolution of species


 
{
    forAll(Y, specieI)
    {
        volScalarField& Yi = Y[specieI];
 
        solve
        (
            fvm::ddt(rho, Yi) - chemistry.RR(specieI),
            mesh.solver("Yi")
        );
    }
}

4 Evolution of energy

The energy equation is derived from the first law of thermodynamics (see also https://en.wikipedia.org/wiki/Enthalpy): For closed systems the change of internal energy is equal to the heat supplied to the system. For open systems the change of enthalpy is equal to the heat supplied to the system. For the first case we get:

_
u = u_0 + q = u_0 + \int \dot q dt
(x)

If we add on both sides the pressure divided by the density we get:


_
u + \frac{p}{\rho}  = u_0 + \frac{p}{\rho} + \int \dot q dt \ \ , \ \  h = u + \frac{p}{\rho}
(x)

 u, u_0 are the internal energy and the internal energy at the time 0.


 
{
    volScalarField& h = thermo.he();
 
    if (constProp == "volume")
    {
        h[0] = u0 + p[0]/rho[0] + integratedHeat;
    }
    else
    {
        h[0] = h0 + integratedHeat;
    }
 
    thermo.correct();
}

With the above equation the evolution of the enthalpy with time can be calculated. Knowing the functional dependency of the heat capacity at constant pressure  c_p a relation between the enthalpy and the temperature can be established. For an example of how it is down in OpenFOAM with a polynomial dependency of the heat capacity from the temperature can be found here https://caefn.com/openfoam/temperature-calculation.

5 Evolution of the pressure

In order to calculate the pressure Daltons model is used: It says that the thermodynamic pressure can be calculated as the sum of the partial pressures of the individual species (see e.g. https://en.wikipedia.org/wiki/Dalton%27s_law):


_
p = \sum p_i = \sum n_i \frac{R T}{V} = \sum \frac{m_i}{M_i} \frac{R T}{V}  = \sum \frac{y_i}{M_i} \underbrace{ \frac{m_{tot}}{V}}_{\rho_0} R T
(x)

 p_i, n_i, m_i, M_i, y_i, m_{tot}, \rho_0 are the partial pressure, the number of molecules, the mass, the molar mass, the mass faction of the species i, the total mass and the initial density, respectively.


 
{
    rho = thermo.rho();
    if (constProp == "volume")
    {
        scalar invW = 0.0;
        forAll(Y, i)
        {
            invW += Y[i][0]/specieData[i].W();
        }
 
        Rspecific[0] = 1000.0*constant::physicoChemical::R.value()*invW;
 
        p[0] = rho0*Rspecific[0]*thermo.T()[0];
        rho[0] = rho0;
    }
}

Note that the constant 1000 in the above code comes from the fact that in OpenFOAM the molar weight is specified in [kg/kmol]. See the file $FOAM_SRC/thermophysicalModels/specie/specie/specie.H.

6 References

  1. Powers, Joseph. Combustion thermodynamics and dynamics. Cambridge University Press, 2016.
  2. Powers, Joseph. Combustion thermodynamics and dynamics. Cambridge University Press, 2016.
  3. Powers, Joseph. Combustion thermodynamics and dynamics. Cambridge University Press, 2016.
  4. Poinsot, Thierry, and Denis Veynante. Theoretical and numerical combustion. RT Edwards, Inc., 2005.
  5. Poinsot, Thierry, and Denis Veynante. Theoretical and numerical combustion. RT Edwards, Inc., 2005.