1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 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 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 "psiuReactionThermo.H"
#include "turbulentFluidThermoModel.H"
#include "multivariateSurfaceInterpolationScheme.H"
#include "fvcDiv.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class XiModel Declaration
\*---------------------------------------------------------------------------*/
class XiModel
{
protected:
// Protected data
dictionary XiModelCoeffs_;
const psiuReactionThermo& 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 psiuReactionThermo& 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 psiuReactionThermo& 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 psiuReactionThermo& 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;
//- Write fields related to Xi model
virtual void writeFields() = 0;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
|