1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 option) any later version.
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
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM; if not, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "SCOPEXiEq.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 defineTypeNameAndDebug(SCOPEXiEq, 0);
37 addToRunTimeSelectionTable(XiEqModel, SCOPEXiEq, dictionary);
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 Foam::XiEqModels::SCOPEXiEq::SCOPEXiEq
46 const dictionary& XiEqProperties,
47 const hhuCombustionThermo& thermo,
48 const compressible::RASModel& turbulence,
49 const volScalarField& Su
52 XiEqModel(XiEqProperties, thermo, turbulence, Su),
53 XiEqCoef(readScalar(XiEqModelCoeffs_.lookup("XiEqCoef"))),
54 XiEqExp(readScalar(XiEqModelCoeffs_.lookup("XiEqExp"))),
55 lCoef(readScalar(XiEqModelCoeffs_.lookup("lCoef"))),
56 SuMin(0.01*Su.average()),
63 "combustionProperties",
64 Su.mesh().time().constant(),
74 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
76 Foam::XiEqModels::SCOPEXiEq::~SCOPEXiEq()
80 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 Foam::tmp<Foam::volScalarField> Foam::XiEqModels::SCOPEXiEq::XiEq() const
84 const volScalarField& k = turbulence_.k();
85 const volScalarField& epsilon = turbulence_.epsilon();
87 volScalarField up = sqrt((2.0/3.0)*k);
88 volScalarField l = (lCoef*sqrt(3.0/2.0))*up*k/epsilon;
89 volScalarField Rl = up*l*thermo_.rhou()/thermo_.muu();
91 volScalarField upBySu = up/(Su_ + SuMin);
92 volScalarField K = 0.157*upBySu/sqrt(Rl);
93 volScalarField Ma = MaModel.Ma();
95 tmp<volScalarField> tXiEq
102 epsilon.time().timeName(),
108 dimensionedScalar("XiEq", dimless, 0.0)
111 volScalarField& xieq = tXiEq();
115 if (Ma[celli] > 0.01)
118 XiEqCoef*pow(K[celli]*Ma[celli], -XiEqExp)*upBySu[celli];
122 forAll(xieq.boundaryField(), patchi)
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];
131 if (Ma[facei] > 0.01)
134 XiEqCoef*pow(Kp[facei]*Map[facei], -XiEqExp)*upBySup[facei];
143 bool Foam::XiEqModels::SCOPEXiEq::read(const dictionary& XiEqProperties)
145 XiEqModel::read(XiEqProperties);
147 XiEqModelCoeffs_.lookup("XiEqCoef") >> XiEqCoef;
148 XiEqModelCoeffs_.lookup("XiEqExp") >> XiEqExp;
149 XiEqModelCoeffs_.lookup("lCoef") >> lCoef;
155 // ************************************************************************* //