Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
33namespace Foam
34{
35namespace XiEqModels
36{
39}
40}
41
42
43// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44
45Foam::XiEqModels::Gulder::Gulder
46(
47 const dictionary& XiEqProperties,
48 const psiuReactionThermo& thermo,
49 const compressible::RASModel& turbulence,
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
69Foam::tmp<Foam::volScalarField> Foam::XiEqModels::Gulder::XiEq() const
70{
71 volScalarField up(sqrt((2.0/3.0)*turbulence_.k()));
72 const tmp<volScalarField> tepsilon(turbulence_.epsilon());
73 const volScalarField& epsilon = tepsilon();
74
75 if (subGridSchelkin_)
76 {
77 up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
78 }
79
80 volScalarField tauEta(sqrt(mag(thermo_.muu()/(thermo_.rhou()*epsilon))));
81
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
95bool Foam::XiEqModels::Gulder::read(const dictionary& XiEqProperties)
96{
97 XiEqModel::read(XiEqProperties);
98
99 XiEqModelCoeffs_.readEntry("XiEqCoef", XiEqCoef_);
100 XiEqModelCoeffs_.readEntry("uPrimeCoef", uPrimeCoef_);
101 XiEqModelCoeffs_.readEntry("subGridSchelkin", subGridSchelkin_);
102
103 return true;
104}
105
106
107// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Simple Gulder model for XiEq based on Gulders correlation with a linear correction function to give a...
Definition Gulder.H:52
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
virtual ~Gulder()
Destructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
scalar epsilon
compressible::turbulenceModel & turbulence
zeroField Su
Definition alphaSuSp.H:1
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition facAverage.C:40
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & e