SCOPEXiEq.H
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 Class
25  Foam::XiEqModels::SCOPEXiEq
26 
27 Description
28  Simple SCOPEXiEq model for XiEq based on SCOPEXiEqs correlation
29  with a linear correction function to give a plausible profile for XiEq.
30  See \link SCOPELaminarFlameSpeed.H \endlink for details on the SCOPE laminar
31  flame speed model.
32 
33 SourceFiles
34  SCOPEXiEq.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef SCOPEXiEq_H
39 #define SCOPEXiEq_H
40 
41 #include "XiEqModel.H"
42 #include "SCOPELaminarFlameSpeed.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace XiEqModels
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class SCOPEXiEq Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class SCOPEXiEq
56 :
57  public XiEqModel
58 {
59  // Private data
60 
61  // Model constant
62  scalar XiEqCoef_;
63 
64  // Model constant
65  scalar XiEqExp_;
66 
67  // Model constant
68  scalar lCoef_;
69 
70  //- Minimum Su
71  dimensionedScalar SuMin_;
72 
73  //- Schelkin effect Model constant
74  scalar uPrimeCoef_;
75 
76  //- Use sub-grid Schelkin effect
77  bool subGridSchelkin_;
78 
79  //- The SCOPE laminar flame speed model used to obtain the
80  // Marstein number. Note: the laminar flame speed need not be
81  // obtained form the same model.
83 
84 
85  // Private Member Functions
86 
87  //- Disallow copy construct
88  SCOPEXiEq(const SCOPEXiEq&);
89 
90  //- Disallow default bitwise assignment
91  void operator=(const SCOPEXiEq&);
92 
93 
94 public:
95 
96  //- Runtime type information
97  TypeName("SCOPEXiEq");
98 
99 
100  // Constructors
101 
102  //- Construct from components
103  SCOPEXiEq
104  (
105  const dictionary& XiEqProperties,
106  const psiuReactionThermo& thermo,
108  const volScalarField& Su
109  );
110 
111 
112  //- Destructor
113  virtual ~SCOPEXiEq();
114 
115 
116  // Member Functions
117 
118  //- Return the flame-wrinking XiEq
119  virtual tmp<volScalarField> XiEq() const;
120 
121  //- Update properties from given dictionary
122  virtual bool read(const dictionary& XiEqProperties);
123 };
124 
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 } // End namespace XiEqModels
129 } // End namespace Foam
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 #endif
134 
135 // ************************************************************************* //
zeroField Su
Definition: alphaSuSp.H:1
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinking XiEq.
rhoReactionThermo & thermo
Definition: createFields.H:28
Simple SCOPEXiEq model for XiEq based on SCOPEXiEqs correlation with a linear correction function to ...
Definition: SCOPEXiEq.H:54
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:57
virtual ~SCOPEXiEq()
Destructor.
Laminar flame speed obtained from the SCOPE correlation.
Foam::psiuReactionThermo.
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
TypeName("SCOPEXiEq")
Runtime type information.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.