This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/combustion/PDRFoam/XiModels/XiModel/XiModel.H

254 lines
6.9 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright held by original author
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
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 inestability 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 inestability.
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 inestability. It is a constant (default 2.5).
SourceFiles
XiModel.C
\*---------------------------------------------------------------------------*/
#ifndef XiModel_H
#define XiModel_H
#include "IOdictionary.H"
#include "hhuCombustionThermo.H"
#include "RASModel.H"
#include "multivariateSurfaceInterpolationScheme.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class XiModel Declaration
\*---------------------------------------------------------------------------*/
class XiModel
{
protected:
// Protected data
dictionary XiModelCoeffs_;
const hhuCombustionThermo& thermo_;
const compressible::RASModel& turbulence_;
const volScalarField& Su_;
const volScalarField& rho_;
const volScalarField& b_;
const surfaceScalarField& phi_;
//- Flame wrinking field
volScalarField Xi_;
private:
// Private Member Functions
//- Disallow copy construct
XiModel(const XiModel&);
//- Disallow default bitwise assignment
void operator=(const XiModel&);
public:
//- Runtime type information
TypeName("XiModel");
// Declare run-time constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
XiModel,
dictionary,
(
const dictionary& XiProperties,
const hhuCombustionThermo& thermo,
const compressible::RASModel& turbulence,
const volScalarField& Su,
const volScalarField& rho,
const volScalarField& b,
const surfaceScalarField& phi
),
(
XiProperties,
thermo,
turbulence,
Su,
rho,
b,
phi
)
);
// Selectors
//- Return a reference to the selected Xi model
static autoPtr<XiModel> New
(
const dictionary& XiProperties,
const hhuCombustionThermo& thermo,
const compressible::RASModel& turbulence,
const volScalarField& Su,
const volScalarField& rho,
const volScalarField& b,
const surfaceScalarField& phi
);
// Constructors
//- Construct from components
XiModel
(
const dictionary& XiProperties,
const hhuCombustionThermo& 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-wrinking Xi
virtual const volScalarField& Xi() const
{
return Xi_;
}
//- Return the flame diffusivity
virtual tmp<volScalarField> Db() const
{
return turbulence_.muEff();
}
//- Add Xi to the multivariateSurfaceInterpolationScheme table
// if required
virtual void addXi
(
multivariateSurfaceInterpolationScheme<scalar>::fieldTable&
)
{}
//- Correct the flame-wrinking Xi
virtual void correct() = 0;
//- Correct the flame-wrinking 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;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //