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 -------------------------------------------------------------------------------
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 "SCOPEXiEq.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace XiEqModels
36 {
37  defineTypeNameAndDebug(SCOPEXiEq, 0);
38  addToRunTimeSelectionTable(XiEqModel, SCOPEXiEq, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
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  XiEqExp_(XiEqModelCoeffs_.get<scalar>("XiEqExp")),
56  lCoef_(XiEqModelCoeffs_.get<scalar>("lCoef")),
57  SuMin_(0.01*Su.average()),
58  uPrimeCoef_(XiEqModelCoeffs_.get<scalar>("uPrimeCoef")),
59  subGridSchelkin_(XiEqModelCoeffs_.get<bool>("subGridSchelkin")),
60  MaModel
61  (
62  Su.mesh().lookupObject<IOdictionary>("combustionProperties"),
63  thermo
64  )
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
77 {
78  const volScalarField& k = turbulence_.k();
79  const volScalarField& epsilon = turbulence_.epsilon();
80 
81  volScalarField up(sqrt((2.0/3.0)*k));
82  if (subGridSchelkin_)
83  {
84  up.primitiveFieldRef() += calculateSchelkinEffect(uPrimeCoef_);
85  }
86 
87  volScalarField l(lCoef_*sqrt(3.0/2.0)*up*k/epsilon);
88  volScalarField Rl(up*l*thermo_.rhou()/thermo_.muu());
89 
90  volScalarField upBySu(up/(Su_ + SuMin_));
91  volScalarField K(0.157*upBySu/sqrt(Rl));
92  volScalarField Ma(MaModel.Ma());
93 
94  tmp<volScalarField> tXiEq
95  (
96  new volScalarField
97  (
98  IOobject
99  (
100  "XiEq",
101  epsilon.time().timeName(),
102  epsilon.db(),
105  ),
106  epsilon.mesh(),
108  )
109  );
110  volScalarField& xieq = tXiEq.ref();
111 
112  forAll(xieq, celli)
113  {
114  if (Ma[celli] > 0.01)
115  {
116  xieq[celli] =
117  XiEqCoef_*pow(K[celli]*Ma[celli], -XiEqExp_)*upBySu[celli];
118  }
119  }
120 
121  volScalarField::Boundary& xieqBf = xieq.boundaryFieldRef();
122 
123  forAll(xieq.boundaryField(), patchi)
124  {
125  scalarField& xieqp = xieqBf[patchi];
126  const scalarField& Kp = K.boundaryField()[patchi];
127  const scalarField& Map = Ma.boundaryField()[patchi];
128  const scalarField& upBySup = upBySu.boundaryField()[patchi];
129 
130  forAll(xieqp, facei)
131  {
132  if (Ma[facei] > 0.01)
133  {
134  xieqp[facei] =
135  XiEqCoef_*pow(Kp[facei]*Map[facei], -XiEqExp_)
136  *upBySup[facei];
137  }
138  }
139  }
140 
141  return tXiEq;
142 }
143 
144 
145 bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
146 {
147  XiEqModel::read(XiEqProperties);
148 
149  XiEqModelCoeffs_.readEntry("XiEqCoef", XiEqCoef_);
150  XiEqModelCoeffs_.readEntry("XiEqExp", XiEqExp_);
151  XiEqModelCoeffs_.readEntry("lCoef", lCoef_);
152  XiEqModelCoeffs_.readEntry("uPrimeCoef", uPrimeCoef_);
153  XiEqModelCoeffs_.readEntry("subGridSchelkin", subGridSchelkin_);
154 
155  return true;
156 }
157 
158 
159 // ************************************************************************* //
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
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 bool read(const dictionary &XiEqProperties)=0
Update properties from given dictionary.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
label k
Boltzmann constant.
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
virtual ~SCOPEXiEq()
Destructor.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
scalar epsilon
Nothing to be read.
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133