instabilityG.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 "instabilityG.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace XiGModels
34 {
35  defineTypeNameAndDebug(instabilityG, 0);
36  addToRunTimeSelectionTable(XiGModel, instabilityG, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 Foam::XiGModels::instabilityG::instabilityG
44 (
45  const dictionary& XiGProperties,
46  const psiuReactionThermo& thermo,
48  const volScalarField& Su
49 )
50 :
51  XiGModel(XiGProperties, thermo, turbulence, Su),
52  GIn_(XiGModelCoeffs_.lookup("GIn")),
53  lambdaIn_(XiGModelCoeffs_.lookup("lambdaIn")),
54  XiGModel_(XiGModel::New(XiGModelCoeffs_, thermo, turbulence, Su))
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
67 {
68  volScalarField turbXiG(XiGModel_->G());
69  return (GIn_*GIn_/(GIn_ + turbXiG) + turbXiG);
70 }
71 
72 
74 {
75  const objectRegistry& db = Su_.db();
76  const volScalarField& Xi = db.lookupObject<volScalarField>("Xi");
77  const volScalarField& rho = db.lookupObject<volScalarField>("rho");
78  const volScalarField& mgb = db.lookupObject<volScalarField>("mgb");
79 
80  return XiGModel_->Db()
81  + rho*Su_*(Xi - 1.0)*mgb*(0.5*lambdaIn_)/(mgb + 1.0/lambdaIn_);
82 }
83 
84 
85 bool Foam::XiGModels::instabilityG::read(const dictionary& XiGProperties)
86 {
87  XiGModel::read(XiGProperties);
88 
89  XiGModelCoeffs_.lookup("GIn") >> GIn_;
90  XiGModelCoeffs_.lookup("lambdaIn") >> lambdaIn_;
91 
92  return true;
93 }
94 
95 
96 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
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
psiReactionThermo & thermo
Definition: createFields.H:31
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual tmp< volScalarField > G() const
Return the flame-wrinking generation rate.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
virtual bool read(const dictionary &XiGProperties)=0
Update properties from given dictionary.
const volScalarField & Su_
Laminar burning velocity.
Definition: XiEqModel.H:74
A class for managing temporary objects.
Definition: PtrList.H:54
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
virtual bool read(const dictionary &XiGProperties)
Update properties from given dictionary.
virtual ~instabilityG()
Destructor.
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Namespace for OpenFOAM.