Gulder.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 {
35  defineTypeNameAndDebug(Gulder, 0);
36  addToRunTimeSelectionTable(XiEqModel, Gulder, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiEqModels::Gulder::Gulder
44 (
45  const dictionary& XiEqProperties,
46  const psiuReactionThermo& thermo,
48  const volScalarField& Su
49 )
50 :
51  XiEqModel(XiEqProperties, thermo, turbulence, Su),
52  XiEqCoef_(readScalar(XiEqModelCoeffs_.lookup("XiEqCoef"))),
53  SuMin_(0.01*Su.average()),
54  uPrimeCoef_(readScalar(XiEqModelCoeffs_.lookup("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  const volScalarField& epsilon = turbulence_.epsilon();
74 
75  if (subGridSchelkin_)
76  {
77  up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
78  }
79 
81 
82  volScalarField Reta
83  (
84  up
85  / (
86  sqrt(epsilon*tauEta)
87  + dimensionedScalar("1e-8", up.dimensions(), 1e-8)
88  )
89  );
90 
91  return (1.0 + XiEqCoef_*sqrt(up/(Su_ + SuMin_))*Reta);
92 }
93 
94 
95 bool Foam::XiEqModels::Gulder::read(const dictionary& XiEqProperties)
96 {
97  XiEqModel::read(XiEqProperties);
98 
99  XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef_;
100  XiEqModelCoeffs_.lookup("uPrimeCoef") >> uPrimeCoef_;
101  XiEqModelCoeffs_.lookup("subGridSchelkin") >> subGridSchelkin_;
102 
103  return true;
104 }
105 
106 
107 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
virtual tmp< volScalarField > muu() const =0
Dynamic viscosity of unburnt gas [kg/ms].
const double e
Elementary charge.
Definition: doubleFloat.H:78
const compressible::RASModel & turbulence_
Definition: XiModel.H:118
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > rhou() const
Unburnt gas density [kg/m^3].
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
bool readBool(Istream &)
Definition: boolIO.C:60
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
const volScalarField & Su_
Definition: XiModel.H:119
virtual ~Gulder()
Destructor.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar epsilon
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
const psiuReactionThermo & thermo_
Definition: XiModel.H:117