algebraic.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "algebraic.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& XiProperties,
48  const volScalarField& Su,
49  const volScalarField& rho,
50  const volScalarField& b,
51  const surfaceScalarField& phi
52 )
53 :
54  XiModel(XiProperties, thermo, turbulence, Su, rho, b, phi),
55  XiShapeCoef(XiModelCoeffs_.lookup<scalar>("XiShapeCoef")),
56  XiEqModel_(XiEqModel::New(XiProperties, thermo, turbulence, Su)),
57  XiGModel_(XiGModel::New(XiProperties, thermo, turbulence, Su))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
70 {
71  return XiGModel_->Db();
72 }
73 
74 
76 {
77  volScalarField XiEqEta(XiEqModel_->XiEq());
78  volScalarField GEta(XiGModel_->G());
79 
80  volScalarField R(GEta*XiEqEta/(XiEqEta - 0.999));
81 
82  volScalarField XiEqStar(R/(R - GEta));
83 
84  Xi_ == 1.0 + (1.0 + (2*XiShapeCoef)*(0.5 - b_))*(XiEqStar - 1.0);
85 }
86 
87 
89 {
90  XiModel::read(XiProperties);
91 
92  XiModelCoeffs_.lookup("XiShapeCoef") >> XiShapeCoef;
93 
94  return true;
95 }
96 
97 
98 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:58
Base-class for all Xi generation models used by the b-Xi combustion model. See Technical Report SH/RE...
Definition: XiGModel.H:55
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:109
virtual bool read(const dictionary &XiProperties)=0
Update properties from given dictionary.
Definition: XiModel.C:86
Simple algebraic model for Xi based on Gulders correlation with a linear correction function to give ...
Definition: algebraic.H:58
virtual ~algebraic()
Destructor.
Definition: algebraic.C:63
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: algebraic.C:75
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
Definition: algebraic.C:88
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: algebraic.C:69
algebraic(const dictionary &XiProperties, const psiuMulticomponentThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Construct from components.
Definition: algebraic.C:44
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
volScalarField & b
Definition: createFields.H:27
addToRunTimeSelectionTable(XiModel, algebraic, dictionary)
defineTypeNameAndDebug(algebraic, 0)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31