Files
OpenFOAM-12/applications/solvers/combustion/PDRFoam/XiModels/XiModel/XiModel.H
Henry Weller b7ea5fcc29 solvers::XiFluid: New solver module for compressible premixed/partially-premixed combustion
executed with foamRun for single region simulations of foamMultiRun for
multi-region simulations.  Replaces XiFoam and all the corresponding
tutorials have been updated and moved to tutorials/modules/XiFluid.

Class
    Foam::solvers::XiFluid

Description
    Solver module for compressible premixed/partially-premixed combustion with
    turbulence modelling.

    Combusting RANS code using the b-Xi two-equation model.
    Xi may be obtained by either the solution of the Xi transport
    equation or from an algebraic expression.  Both approaches are
    based on Gulder's flame speed correlation which has been shown
    to be appropriate by comparison with the results from the
    spectral model.

    Strain effects are encorporated directly into the Xi equation
    but not in the algebraic approximation.  Further work need to be
    done on this issue, particularly regarding the enhanced removal rate
    caused by flame compression.  Analysis using results of the spectral
    model will be required.

    For cases involving very lean Propane flames or other flames which are
    very strain-sensitive, a transport equation for the laminar flame
    speed is present.  This equation is derived using heuristic arguments
    involving the strain time scale and the strain-rate at extinction.
    the transport velocity is the same as that for the Xi equation.

    Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
    pseudo-transient and steady simulations.

    Optional fvModels and fvConstraints are provided to enhance the simulation
    in many ways including adding various sources, chemical reactions,
    combustion, Lagrangian particles, radiation, surface film etc. and
    constraining or limiting the solution.

    Reference:
    \verbatim
        Greenshields, C. J., & Weller, H. G. (2022).
        Notes on Computational Fluid Dynamics: General Principles.
        CFD Direct Ltd.: Reading, UK.
    \endverbatim

SourceFiles
    XiFluid.C

See also
    Foam::solvers::fluidSolver
    Foam::solvers::isothermalFluid
2022-12-29 23:53:33 +00:00

258 lines
7.2 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
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/>.
Class
Foam::XiModel
Description
Base-class for all Xi models used by the b-Xi combustion model.
See Technical Report SH/RE/01R for details on the PDR modelling.
Xi is given through an algebraic expression (\link algebraic.H \endlink),
by solving a transport equation (\link transport.H \endlink) or a
fixed value (\link fixed.H \endlink).
See report TR/HGW/10 for details on the Weller two equations model.
In the algebraic and transport methods \f$\Xi_{eq}\f$ is calculated in
similar way. In the algebraic approach, \f$\Xi_{eq}\f$ is the value used in
the \f$ b \f$ transport equation.
\f$\Xi_{eq}\f$ is calculated as follows:
\f$\Xi_{eq} = 1 + (1 + 2\Xi_{coeff}(0.5 - \dwea{b}))(\Xi^* - 1)\f$
where:
\f$ \dwea{b} \f$ is the regress variable.
\f$ \Xi_{coeff} \f$ is a model constant.
\f$ \Xi^* \f$ is the total equilibrium wrinkling combining the effects
of the flame instability and turbulence interaction and is given by
\f[
\Xi^* = \frac {R}{R - G_\eta - G_{in}}
\f]
where:
\f$ G_\eta \f$ is the generation rate of wrinkling due to turbulence
interaction.
\f$ G_{in} = \kappa \rho_{u}/\rho_{b} \f$ is the generation
rate due to the flame instability.
By adding the removal rates of the two effects:
\f[
R = G_\eta \frac{\Xi_{\eta_{eq}}}{\Xi_{\eta_{eq}} - 1}
+ G_{in} \frac{\Xi_{{in}_{eq}}}{\Xi_{{in}_{eq}} - 1}
\f]
where:
\f$ R \f$ is the total removal.
\f$ G_\eta \f$ is a model constant.
\f$ \Xi_{\eta_{eq}} \f$ is the flame wrinkling due to turbulence.
\f$ \Xi_{{in}_{eq}} \f$ is the equilibrium level of the flame wrinkling
generated by instability. It is a constant (default 2.5).
SourceFiles
XiModel.C
\*---------------------------------------------------------------------------*/
#ifndef XiModel_H
#define XiModel_H
#include "IOdictionary.H"
#include "psiuMulticomponentThermo.H"
#include "compressibleMomentumTransportModels.H"
#include "multivariateSurfaceInterpolationScheme.H"
#include "fvcDiv.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class XiModel Declaration
\*---------------------------------------------------------------------------*/
class XiModel
{
protected:
// Protected data
dictionary XiModelCoeffs_;
const psiuMulticomponentThermo& thermo_;
const compressible::RASModel& turbulence_;
const volScalarField& Su_;
const volScalarField& rho_;
const volScalarField& b_;
const surfaceScalarField& phi_;
//- Flame wrinkling field
volScalarField Xi_;
public:
//- Runtime type information
TypeName("XiModel");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
XiModel,
dictionary,
(
const dictionary& XiProperties,
const psiuMulticomponentThermo& thermo,
const compressible::RASModel& turbulence,
const volScalarField& Su,
const volScalarField& rho,
const volScalarField& b,
const surfaceScalarField& phi
),
(
XiProperties,
thermo,
turbulence,
Su,
rho,
b,
phi
)
);
// Constructors
//- Construct from components
XiModel
(
const dictionary& XiProperties,
const psiuMulticomponentThermo& thermo,
const compressible::RASModel& turbulence,
const volScalarField& Su,
const volScalarField& rho,
const volScalarField& b,
const surfaceScalarField& phi
);
//- Disallow default bitwise copy construction
XiModel(const XiModel&);
// Selectors
//- Return a reference to the selected Xi model
static autoPtr<XiModel> New
(
const dictionary& XiProperties,
const psiuMulticomponentThermo& thermo,
const compressible::RASModel& turbulence,
const volScalarField& Su,
const volScalarField& rho,
const volScalarField& b,
const surfaceScalarField& phi
);
//- Destructor
virtual ~XiModel();
// Member Functions
//- Return the flame-wrinkling Xi
virtual const volScalarField& Xi() const
{
return Xi_;
}
//- Return the flame diffusivity
virtual tmp<volScalarField> Db() const
{
return volScalarField::New
(
"Db",
turbulence_.rho()*turbulence_.nuEff()
);
}
//- Add Xi to the multivariateSurfaceInterpolationScheme table
// if required
virtual void addXi
(
multivariateSurfaceInterpolationScheme<scalar>::fieldTable&
)
{}
//- Correct the flame-wrinkling Xi
virtual void correct() = 0;
//- Correct the flame-wrinkling Xi using the given convection scheme
virtual void correct(const fv::convectionScheme<scalar>&)
{
correct();
}
//- Update properties from given dictionary
virtual bool read(const dictionary& XiProperties) = 0;
//- Write fields related to Xi model
virtual void writeFields() = 0;
// Member Operators
//- Disallow default bitwise assignment
void operator=(const XiModel&) = delete;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //