Loading...
Searching...
No Matches
SCOPEXiEq.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 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "SCOPEXiEq.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace XiEqModels
37{
40}
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
47(
48 const dictionary& XiEqProperties,
49 const psiuReactionThermo& thermo,
50 const compressible::RASModel& turbulence,
51 const volScalarField& Su
52)
53:
54 XiEqModel(XiEqProperties, thermo, turbulence, Su),
55 XiEqCoef_(XiEqModelCoeffs_.get<scalar>("XiEqCoef")),
56 XiEqExp_(XiEqModelCoeffs_.get<scalar>("XiEqExp")),
57 lCoef_(XiEqModelCoeffs_.get<scalar>("lCoef")),
58 SuMin_(0.01*Su.average()),
59 uPrimeCoef_(XiEqModelCoeffs_.get<scalar>("uPrimeCoef")),
60 subGridSchelkin_(XiEqModelCoeffs_.get<bool>("subGridSchelkin")),
61 MaModel
62 (
63 Su.mesh().lookupObject<IOdictionary>("combustionProperties"),
64 thermo
65 )
66{}
67
68
69// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70
72{}
73
74
75// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
76
77Foam::tmp<Foam::volScalarField> Foam::XiEqModels::SCOPEXiEq::XiEq() const
78{
79 const tmp<volScalarField> tk(turbulence_.k());
80 const volScalarField& k = tk();
81 const tmp<volScalarField> tepsilon(turbulence_.epsilon());
82 const volScalarField& epsilon = tepsilon();
83
84 volScalarField up(sqrt((2.0/3.0)*k));
85 if (subGridSchelkin_)
86 {
87 up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
88 }
89
90 volScalarField l(lCoef_*sqrt(3.0/2.0)*up*k/epsilon);
91 volScalarField Rl(up*l*thermo_.rhou()/thermo_.muu());
92
93 volScalarField upBySu(up/(Su_ + SuMin_));
94 volScalarField K(0.157*upBySu/sqrt(Rl));
95 volScalarField Ma(MaModel.Ma());
96
97 auto tXiEq = volScalarField::New
98 (
99 "XiEq",
101 epsilon.mesh(),
103 );
104 auto& xieq = tXiEq.ref();
105
106 forAll(xieq, celli)
107 {
108 if (Ma[celli] > 0.01)
109 {
110 xieq[celli] =
111 XiEqCoef_*pow(K[celli]*Ma[celli], -XiEqExp_)*upBySu[celli];
112 }
113 }
114
115 volScalarField::Boundary& xieqBf = xieq.boundaryFieldRef();
116
117 forAll(xieq.boundaryField(), patchi)
118 {
119 scalarField& xieqp = xieqBf[patchi];
120 const scalarField& Kp = K.boundaryField()[patchi];
121 const scalarField& Map = Ma.boundaryField()[patchi];
122 const scalarField& upBySup = upBySu.boundaryField()[patchi];
123
124 forAll(xieqp, facei)
125 {
126 if (Ma[facei] > 0.01)
127 {
128 xieqp[facei] =
129 XiEqCoef_*pow(Kp[facei]*Map[facei], -XiEqExp_)
130 *upBySup[facei];
131 }
132 }
133 }
134
135 return tXiEq;
136}
137
138
139bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
140{
141 XiEqModel::read(XiEqProperties);
142
143 XiEqModelCoeffs_.readEntry("XiEqCoef", XiEqCoef_);
144 XiEqModelCoeffs_.readEntry("XiEqExp", XiEqExp_);
145 XiEqModelCoeffs_.readEntry("lCoef", lCoef_);
146 XiEqModelCoeffs_.readEntry("uPrimeCoef", uPrimeCoef_);
147 XiEqModelCoeffs_.readEntry("subGridSchelkin", subGridSchelkin_);
148
149 return true;
150}
151
152
153// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Simple SCOPEXiEq model for XiEq based on SCOPEXiEqs correlation with a linear correction function to ...
Definition SCOPEXiEq.H:55
virtual ~SCOPEXiEq()
Destructor.
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
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
dynamicFvMesh & mesh
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.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299