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 | Copyright (C) 2011-2015 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 
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.internalField() += 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();
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  forAll(xieq.boundaryField(), patchi)
123  {
124  scalarField& xieqp = xieq.boundaryField()[patchi];
125  const scalarField& Kp = K.boundaryField()[patchi];
126  const scalarField& Map = Ma.boundaryField()[patchi];
127  const scalarField& upBySup = upBySu.boundaryField()[patchi];
128 
129  forAll(xieqp, facei)
130  {
131  if (Ma[facei] > 0.01)
132  {
133  xieqp[facei] =
134  XiEqCoef_*pow(Kp[facei]*Map[facei], -XiEqExp_)
135  *upBySup[facei];
136  }
137  }
138  }
139 
140  return tXiEq;
141 }
142 
143 
144 bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
145 {
146  XiEqModel::read(XiEqProperties);
147 
148  XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef_;
149  XiEqModelCoeffs_.lookup("XiEqExp") >> XiEqExp_;
150  XiEqModelCoeffs_.lookup("lCoef") >> lCoef_;
151  XiEqModelCoeffs_.lookup("uPrimeCoef") >> uPrimeCoef_;
152  XiEqModelCoeffs_.lookup("subGridSchelkin") >> subGridSchelkin_;
153 
154  return true;
155 }
156 
157 
158 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::fvc::average
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
Foam::XiEqModels::SCOPEXiEq::XiEq
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
Foam::XiEqModel::read
virtual bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::XiEqModels::SCOPEXiEq::read
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
SCOPEXiEq(const SCOPEXiEq &)
Disallow copy construct.
SCOPEXiEq.H
epsilon
epsilon
Definition: createTDFields.H:56
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
RASModel
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);singlePhaseTransportModel laminarTransport(U, phi);autoPtr< incompressible::RASModel > RASModel(incompressible::New< incompressible::RASModel >(U, phi, laminarTransport))
patchi
label patchi
Definition: getPatchFieldScalar.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::XiEqModels::SCOPEXiEq::~SCOPEXiEq
virtual ~SCOPEXiEq()
Destructor.
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress