basicXiSubG.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 "basicXiSubG.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace XiGModels
35 {
36  defineTypeNameAndDebug(basicSubGrid, 0);
37  addToRunTimeSelectionTable(XiGModel, basicSubGrid, dictionary);
38 };
39 };
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
44 Foam::XiGModels::basicSubGrid::basicSubGrid
45 (
46  const dictionary& XiGProperties,
47  const psiuReactionThermo& thermo,
49  const volScalarField& Su
50 )
51 :
52  XiGModel(XiGProperties, thermo, turbulence, Su),
53  k1(readScalar(XiGModelCoeffs_.lookup("k1"))),
54  XiGModel_(XiGModel::New(XiGModelCoeffs_, thermo, turbulence, Su))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
67 {
68  const objectRegistry& db = Su_.db();
69  const volVectorField& U = db.lookupObject<volVectorField>("U");
70  const volScalarField& Nv = db.lookupObject<volScalarField>("Nv");
71  const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
72 
73  tmp<volScalarField> tGtot = XiGModel_->G();
74  volScalarField& Gtot = tGtot();
75 
76  const scalarField Cw = pow(Su_.mesh().V(), 2.0/3.0);
77 
78  tmp<volScalarField> tN
79  (
80  new volScalarField
81  (
82  IOobject
83  (
84  "tN",
85  Su_.mesh().time().timeName(),
86  Su_.mesh(),
89  false
90  ),
91  Su_.mesh(),
92  dimensionedScalar("zero", Nv.dimensions(), 0.0),
93  zeroGradientFvPatchVectorField::typeName
94  )
95  );
96 
97  volScalarField& N = tN();
98 
99  N.internalField() = Nv.internalField()*Cw;
100 
101  forAll(N, celli)
102  {
103  if (N[celli] > 1e-3)
104  {
105  Gtot[celli] += k1*mag(U[celli])/Lobs[celli];
106  }
107  }
108 
109  return tGtot;
110 }
111 
112 
114 {
115  const objectRegistry& db = Su_.db();
116  const volScalarField& Xi = db.lookupObject<volScalarField>("Xi");
117  const volScalarField& rho = db.lookupObject<volScalarField>("rho");
118  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
119  const volScalarField& Lobs = db.lookupObject<volScalarField>("Lobs");
120 
121  return XiGModel_->Db()
122  + rho*Su_*(Xi - 1.0)*mgb*(0.5*Lobs)*Lobs/(mgb*Lobs + 1.0);
123 }
124 
125 
126 bool Foam::XiGModels::basicSubGrid::read(const dictionary& XiGProperties)
127 {
128  XiGModel::read(XiGProperties);
129 
130  XiGModelCoeffs_.lookup("k1") >> k1;
131 
132  return true;
133 }
134 
135 
136 // ************************************************************************* //
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
const Mesh & mesh() const
Return mesh.
virtual bool read(const dictionary &XiGProperties)
Update properties from given dictionary.
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Namespace for OpenFOAM.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
const double e
Elementary charge.
Definition: doubleFloat.H:78
stressControl lookup("compactNormalStress") >> compactNormalStress
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Macros for easy insertion into run-time selection tables.
virtual bool read(const dictionary &XiGProperties)=0
Update properties from given dictionary.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual tmp< volScalarField > G() const
Return the flame-wrinking generation rate.
psiReactionThermo & thermo
Definition: createFields.H:32
virtual ~basicSubGrid()
Destructor.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:245
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
defineTypeNameAndDebug(combustionModel, 0)