SCOPEXiEq.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2013 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::XiEqModels::SCOPEXiEq
28 
29 Description
30  Simple SCOPEXiEq model for XiEq based on SCOPEXiEqs correlation
31  with a linear correction function to give a plausible profile for XiEq.
32  See \link SCOPELaminarFlameSpeed.H \endlink for details on the SCOPE laminar
33  flame speed model.
34 
35 SourceFiles
36  SCOPEXiEq.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef SCOPEXiEq_H
41 #define SCOPEXiEq_H
42 
43 #include "XiEqModel.H"
44 #include "SCOPELaminarFlameSpeed.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 namespace XiEqModels
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class SCOPEXiEq Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class SCOPEXiEq
58 :
59  public XiEqModel
60 {
61  // Private data
62 
63  // Model constant
64  scalar XiEqCoef_;
65 
66  // Model constant
67  scalar XiEqExp_;
68 
69  // Model constant
70  scalar lCoef_;
71 
72  //- Minimum Su
73  dimensionedScalar SuMin_;
74 
75  //- Schelkin effect Model constant
76  scalar uPrimeCoef_;
77 
78  //- Use sub-grid Schelkin effect
79  bool subGridSchelkin_;
80 
81  //- The SCOPE laminar flame speed model used to obtain the
82  // Marstein number. Note: the laminar flame speed need not be
83  // obtained form the same model.
85 
86 
87  // Private Member Functions
88 
89  //- No copy construct
90  SCOPEXiEq(const SCOPEXiEq&) = delete;
91 
92  //- No copy assignment
93  void operator=(const SCOPEXiEq&) = delete;
94 
95 
96 public:
97 
98  //- Runtime type information
99  TypeName("SCOPEXiEq");
100 
101 
102  // Constructors
103 
104  //- Construct from components
105  SCOPEXiEq
106  (
107  const dictionary& XiEqProperties,
108  const psiuReactionThermo& thermo,
110  const volScalarField& Su
111  );
112 
113 
114  //- Destructor
115  virtual ~SCOPEXiEq();
116 
117 
118  // Member Functions
119 
120  //- Return the flame-wrinkling XiEq
121  virtual tmp<volScalarField> XiEq() const;
122 
123  //- Update properties from given dictionary
124  virtual bool read(const dictionary& XiEqProperties);
125 };
126 
127 
128 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
129 
130 } // End namespace XiEqModels
131 } // End namespace Foam
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 #endif
136 
137 // ************************************************************************* //
Foam::XiEqModel
Base-class for all XiEq models used by the b-XiEq combustion model. The available models are : basicX...
Definition: XiEqModel.H:59
turbulence
Info<< "Reading field U\n"<< 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\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::psiuReactionThermo
Foam::psiuReactionThermo.
Definition: psiuReactionThermo.H:55
SCOPELaminarFlameSpeed.H
XiEqModel.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::XiEqModels::SCOPEXiEq::XiEq
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Foam::laminarFlameSpeedModels::SCOPE
Laminar flame speed obtained from the SCOPE correlation.
Definition: SCOPELaminarFlameSpeed.H:79
Foam::XiEqModels::SCOPEXiEq::read
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::XiEqModels::SCOPEXiEq::TypeName
TypeName("SCOPEXiEq")
Runtime type information.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:66
Foam::XiEqModels::SCOPEXiEq
Simple SCOPEXiEq model for XiEq based on SCOPEXiEqs correlation with a linear correction function to ...
Definition: SCOPEXiEq.H:56
Foam::XiEqModels::SCOPEXiEq::~SCOPEXiEq
virtual ~SCOPEXiEq()
Destructor.
Foam::GeometricField< scalar, fvPatchField, volMesh >