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-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 "Gulder.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiEqModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& XiEqProperties,
48  const volScalarField& Su
49 )
50 :
51  XiEqModel(XiEqProperties, thermo, turbulence, Su),
52  XiEqCoef_(XiEqModelCoeffs_.lookup<scalar>("XiEqCoef")),
53  SuMin_(0.01*Su.average()),
54  uPrimeCoef_(XiEqModelCoeffs_.lookup<scalar>("uPrimeCoef")),
55  subGridSchelkin_
56  (
57  readBool(XiEqModelCoeffs_.lookup("subGridSchelkin"))
58  )
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
71 {
72  volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
73 
74  tmp<volScalarField> tepsilon = turbulence_.epsilon();
75  const volScalarField& epsilon = tepsilon();
76 
77  if (subGridSchelkin_)
78  {
79  up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
80  }
81 
82  volScalarField tauEta(sqrt(mag(thermo_.muu()/(thermo_.rhou()*epsilon))));
83 
84  volScalarField Reta
85  (
86  up
87  / (
88  sqrt(epsilon*tauEta)
89  + dimensionedScalar(up.dimensions(), 1e-8)
90  )
91  );
92 
93  return (1.0 + XiEqCoef_*sqrt(up/(Su_ + SuMin_))*Reta);
94 }
95 
96 
97 bool Foam::XiEqModels::Gulder::read(const dictionary& XiEqProperties)
98 {
99  XiEqModel::read(XiEqProperties);
100 
101  XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef_;
102  XiEqModelCoeffs_.lookup("uPrimeCoef") >> uPrimeCoef_;
103  XiEqModelCoeffs_.lookup("subGridSchelkin") >> subGridSchelkin_;
104 
105  return true;
106 }
107 
108 
109 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:58
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Definition: XiEqModel.C:68
Simple Gulder model for XiEq based on Gulders correlation with a linear correction function to give a...
Definition: Gulder.H:54
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Definition: Gulder.C:97
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Definition: Gulder.C:70
Gulder(const dictionary &XiEqProperties, const psiuMulticomponentThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su)
Construct from components.
Definition: Gulder.C:44
virtual ~Gulder()
Destructor.
Definition: Gulder.C:64
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
const scalar epsilon
addToRunTimeSelectionTable(XiEqModel, basicSubGrid, dictionary)
defineTypeNameAndDebug(basicSubGrid, 0)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
bool readBool(Istream &)
Definition: boolIO.C:60
const doubleScalar e
Definition: doubleScalar.H:105
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31