basicXiSubG.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 "basicXiSubG.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiGModels
34 {
35  defineTypeNameAndDebug(basicSubGrid, 0);
36  addToRunTimeSelectionTable(XiGModel, basicSubGrid, dictionary);
37 };
38 };
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiGModels::basicSubGrid::basicSubGrid
44 (
45  const dictionary& XiGProperties,
46  const psiuReactionThermo& thermo,
48  const volScalarField& Su
49 )
50 :
51  XiGModel(XiGProperties, thermo, turbulence, Su),
52  k1(readScalar(XiGModelCoeffs_.lookup("k1"))),
53  XiGModel_(XiGModel::New(XiGModelCoeffs_, thermo, turbulence, Su))
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
64 
66 {
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");
71 
72  tmp<volScalarField> tGtot = XiGModel_->G();
73  volScalarField& Gtot = tGtot.ref();
74 
75  const scalarField Cw = pow(Su_.mesh().V(), 2.0/3.0);
76  scalarField N(Nv.primitiveField()*Cw);
77 
78  forAll(N, celli)
79  {
80  if (N[celli] > 1e-3)
81  {
82  Gtot[celli] += k1*mag(U[celli])/Lobs[celli];
83  }
84  }
85 
86  return tGtot;
87 }
88 
89 
91 {
92  const objectRegistry& db = Su_.db();
93  const volScalarField& Xi = db.lookupObject<volScalarField>("Xi");
94  const volScalarField& rho = db.lookupObject<volScalarField>("rho");
95  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
96  const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
97 
98  return XiGModel_->Db()
99  + rho*Su_*(Xi - 1.0)*mgb*(0.5*Lobs)*Lobs/(mgb*Lobs + 1.0);
100 }
101 
102 
103 bool Foam::XiGModels::basicSubGrid::read(const dictionary& XiGProperties)
104 {
105  XiGModel::read(XiGProperties);
106 
107  XiGModelCoeffs_.lookup("k1") >> k1;
108 
109  return true;
110 }
111 
112 
113 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
zeroField Su
Definition: alphaSuSp.H:1
virtual tmp< volScalarField > G() const
Return the flame-wrinking generation rate.
virtual bool read(const dictionary &XiGProperties)
Update properties from given dictionary.
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
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.
virtual bool read(const dictionary &XiGProperties)=0
Update properties from given dictionary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual ~basicSubGrid()
Destructor.
label N
Definition: createFields.H:22
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:361
Namespace for OpenFOAM.