Gulder.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-2024 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 "Gulder.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiEqModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
46 
47  XiEqCoeff_.readIfPresent(dict);
48  SuMin_.readIfPresent(dict);
49 
50  return true;
51 }
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const dictionary& dict,
60  const volScalarField& Su
61 )
62 :
64  XiEqCoeff_("XiEqCoeff", dimless, 0.62),
65  SuMin_("SuMin", 0.01*Su.average())
66 {
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
80 {
81  const volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
82 
83  tmp<volScalarField> tepsilon = turbulence_.epsilon();
84  const volScalarField& epsilon = tepsilon();
85 
86  const volScalarField tauEta
87  (
88  sqrt(mag(thermo_.muu()/(thermo_.rhou()*epsilon)))
89  );
90 
91  const volScalarField Reta
92  (
93  up
94  /(
95  sqrt(epsilon*tauEta)
96  + dimensionedScalar(up.dimensions(), 1e-8)
97  )
98  );
99 
100  return (1 + XiEqCoeff_*sqrt(up/max(Su_, SuMin_))*Reta);
101 }
102 
103 
104 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Templated abstract base class for thermophysical transport models.
Base class for equilibrium flame wrinkling XiEq models.
Definition: XiEqModel.H:52
virtual bool readCoeffs(const dictionary &dict)=0
Update coefficients from given dictionary.
Definition: XiEqModel.C:39
Simple model for the equilibrium flame wrinkling XiEq based on Gulder's turbulent flame speed correla...
Definition: Gulder.H:74
Gulder(const dictionary &dict, const psiuMulticomponentThermo &thermo, const fluidThermoThermophysicalTransportModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: Gulder.C:56
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
Definition: Gulder.C:43
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Definition: Gulder.C:79
virtual ~Gulder()
Destructor.
Definition: Gulder.C:73
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
bool readIfPresent(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type> if found in the dictionary.
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
const scalar epsilon
defineTypeNameAndDebug(constant, 0)
addToRunTimeSelectionTable(XiEqModel, constant, dictionary)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
const dimensionSet dimless
dimensioned< Type > average(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31