SCOPEXiEq.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-2018 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 "SCOPEXiEq.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiEqModels
34 {
35  defineTypeNameAndDebug(SCOPEXiEq, 0);
36  addToRunTimeSelectionTable(XiEqModel, SCOPEXiEq, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
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  XiEqExp_(readScalar(XiEqModelCoeffs_.lookup("XiEqExp"))),
54  lCoef_(readScalar(XiEqModelCoeffs_.lookup("lCoef"))),
55  SuMin_(0.01*Su.average()),
56  uPrimeCoef_(readScalar(XiEqModelCoeffs_.lookup("uPrimeCoef"))),
57  subGridSchelkin_
58  (
59  readBool(XiEqModelCoeffs_.lookup("subGridSchelkin"))
60  ),
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 
78 {
79  const volScalarField& k = turbulence_.k();
80  const volScalarField& epsilon = turbulence_.epsilon();
81 
82  volScalarField up(sqrt((2.0/3.0)*k));
83  if (subGridSchelkin_)
84  {
85  up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
86  }
87 
88  volScalarField l(lCoef_*sqrt(3.0/2.0)*up*k/epsilon);
89  volScalarField Rl(up*l*thermo_.rhou()/thermo_.muu());
90 
91  volScalarField upBySu(up/(Su_ + SuMin_));
92  volScalarField K(0.157*upBySu/sqrt(Rl));
93  volScalarField Ma(MaModel.Ma());
94 
95  tmp<volScalarField> tXiEq
96  (
97  new volScalarField
98  (
99  IOobject
100  (
101  "XiEq",
102  epsilon.time().timeName(),
103  epsilon.db(),
106  ),
107  epsilon.mesh(),
108  dimensionedScalar("XiEq", dimless, 0.0)
109  )
110  );
111  volScalarField& xieq = tXiEq.ref();
112 
113  forAll(xieq, celli)
114  {
115  if (Ma[celli] > 0.01)
116  {
117  xieq[celli] =
118  XiEqCoef_*pow(K[celli]*Ma[celli], -XiEqExp_)*upBySu[celli];
119  }
120  }
121 
122  volScalarField::Boundary& xieqBf = xieq.boundaryFieldRef();
123 
124  forAll(xieq.boundaryField(), patchi)
125  {
126  scalarField& xieqp = xieqBf[patchi];
127  const scalarField& Kp = K.boundaryField()[patchi];
128  const scalarField& Map = Ma.boundaryField()[patchi];
129  const scalarField& upBySup = upBySu.boundaryField()[patchi];
130 
131  forAll(xieqp, facei)
132  {
133  if (Ma[facei] > 0.01)
134  {
135  xieqp[facei] =
136  XiEqCoef_*pow(Kp[facei]*Map[facei], -XiEqExp_)
137  *upBySup[facei];
138  }
139  }
140  }
141 
142  return tXiEq;
143 }
144 
145 
146 bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
147 {
148  XiEqModel::read(XiEqProperties);
149 
150  XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef_;
151  XiEqModelCoeffs_.lookup("XiEqExp") >> XiEqExp_;
152  XiEqModelCoeffs_.lookup("lCoef") >> lCoef_;
153  XiEqModelCoeffs_.lookup("uPrimeCoef") >> uPrimeCoef_;
154  XiEqModelCoeffs_.lookup("subGridSchelkin") >> subGridSchelkin_;
155 
156  return true;
157 }
158 
159 
160 // ************************************************************************* //
dictionary XiEqModelCoeffs_
Dictionary.
Definition: XiEqModel.H:65
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
zeroField Su
Definition: alphaSuSp.H:1
virtual tmp< volScalarField > muu() const =0
Dynamic viscosity of unburnt gas [kg/ms].
tmp< volScalarField > calculateSchelkinEffect(const scalar) const
Return the sub-grid Schelkin effect.
virtual tmp< volScalarField > rhou() const
Unburnt gas density [kg/m^3].
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
bool readBool(Istream &)
Definition: boolIO.C:60
rhoReactionThermo & thermo
Definition: createFields.H:28
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
virtual ~SCOPEXiEq()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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 successful.
Definition: doubleScalar.H:68
const psiuReactionThermo & thermo_
Thermo.
Definition: XiEqModel.H:68
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
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
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
const compressible::RASModel & turbulence_
Turbulence.
Definition: XiEqModel.H:71