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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "Gulder.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace XiEqModels
36 {
37  defineTypeNameAndDebug(Gulder, 0);
38  addToRunTimeSelectionTable(XiEqModel, Gulder, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::XiEqModels::Gulder::Gulder
46 (
47  const dictionary& XiEqProperties,
48  const psiuReactionThermo& thermo,
50  const volScalarField& Su
51 )
52 :
53  XiEqModel(XiEqProperties, thermo, turbulence, Su),
54  XiEqCoef_(XiEqModelCoeffs_.get<scalar>("XiEqCoef")),
55  SuMin_(0.01*Su.average()),
56  uPrimeCoef_(XiEqModelCoeffs_.get<scalar>("uPrimeCoef")),
57  subGridSchelkin_(XiEqModelCoeffs_.get<bool>("subGridSchelkin"))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
70 {
71  volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
72  const volScalarField& epsilon = turbulence_.epsilon();
73 
74  if (subGridSchelkin_)
75  {
76  up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
77  }
78 
79  volScalarField tauEta(sqrt(mag(thermo_.muu()/(thermo_.rhou()*epsilon))));
80 
81  volScalarField Reta
82  (
83  up
84  / (
85  sqrt(epsilon*tauEta)
86  + dimensionedScalar("1e-8", up.dimensions(), 1e-8)
87  )
88  );
89 
90  return (1.0 + XiEqCoef_*sqrt(up/(Su_ + SuMin_))*Reta);
91 }
92 
93 
94 bool Foam::XiEqModels::Gulder::read(const dictionary& XiEqProperties)
95 {
96  XiEqModel::read(XiEqProperties);
97 
98  XiEqModelCoeffs_.readEntry("XiEqCoef", XiEqCoef_);
99  XiEqModelCoeffs_.readEntry("uPrimeCoef", uPrimeCoef_);
100  XiEqModelCoeffs_.readEntry("subGridSchelkin", subGridSchelkin_);
101 
102  return true;
103 }
104 
105 
106 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
zeroField Su
Definition: alphaSuSp.H:1
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:40
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual ~Gulder()
Destructor.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
defineTypeNameAndDebug(combustionModel, 0)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar epsilon
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)