1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2011 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
13 the Free Software Foundation, either version 3 of the License, or
14 (at your 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, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "basicXiSubG.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 defineTypeNameAndDebug(basicSubGrid, 0);
36 addToRunTimeSelectionTable(XiGModel, basicSubGrid, dictionary);
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 Foam::XiGModels::basicSubGrid::basicSubGrid
45 const dictionary& XiGProperties,
46 const hhuCombustionThermo& thermo,
47 const compressible::RASModel& turbulence,
48 const volScalarField& Su
51 XiGModel(XiGProperties, thermo, turbulence, Su),
52 k1(readScalar(XiGModelCoeffs_.lookup("k1"))),
53 XiGModel_(XiGModel::New(XiGModelCoeffs_, thermo, turbulence, Su))
57 // * * * * * * * * * * * * * * * * Destructors * * * * * * * * * * * * * * * //
59 Foam::XiGModels::basicSubGrid::~basicSubGrid()
63 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 Foam::tmp<Foam::volScalarField> Foam::XiGModels::basicSubGrid::G() const
67 const objectRegistry& db = Su_.db();
68 const volVectorField& U = db.lookupObject<volVectorField>("U");
69 const volScalarField& Nv = db.lookupObject<volScalarField>("Nv");
70 const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
72 tmp<volScalarField> tGtot = XiGModel_->G();
73 volScalarField& Gtot = tGtot();
75 const scalarField Cw = pow(Su_.mesh().V(), 2.0/3.0);
77 tmp<volScalarField> tN
84 Su_.mesh().time().timeName(),
91 dimensionedScalar("zero", Nv.dimensions(), 0.0),
92 zeroGradientFvPatchVectorField::typeName
96 volScalarField& N = tN();
98 N.internalField() = Nv.internalField()*Cw;
104 Gtot[celli] += k1*mag(U[celli])/Lobs[celli];
112 Foam::tmp<Foam::volScalarField> Foam::XiGModels::basicSubGrid::Db() const
114 const objectRegistry& db = Su_.db();
115 const volScalarField& Xi = db.lookupObject<volScalarField>("Xi");
116 const volScalarField& rho = db.lookupObject<volScalarField>("rho");
117 const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
118 const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
120 return XiGModel_->Db()
121 + rho*Su_*(Xi - 1.0)*mgb*(0.5*Lobs)*Lobs/(mgb*Lobs + 1.0);
125 bool Foam::XiGModels::basicSubGrid::read(const dictionary& XiGProperties)
127 XiGModel::read(XiGProperties);
129 XiGModelCoeffs_.lookup("k1") >> k1;
135 // ************************************************************************* //