basicXiSubXiEq.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-2012 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::basicSubGrid
28 
29 Description
30  Basic sub-grid obstacle flame-wrinkling enhancement factor model.
31  Details supplied by J Puttock 2/7/06.
32 
33  <b> Sub-grid flame area generation </b>
34 
35  \f$ n = N - \hat{\dwea{\vec{U}}}.n_{s}.\hat{\dwea{\vec{U}}} \f$
36  \f$ n_{r} = \sqrt{n} \f$
37 
38  where:
39 
40  \f$ \hat{\dwea{\vec{U}}} = \dwea{\vec{U}} / \vert \dwea{\vec{U}}
41  \vert \f$
42 
43  \f$ b = \hat{\dwea{\vec{U}}}.B.\hat{\dwea{\vec{U}}} / n_{r} \f$
44 
45  where:
46 
47  \f$ B \f$ is the file "B".
48 
49  \f$ N \f$ is the file "N".
50 
51  \f$ n_{s} \f$ is the file "ns".
52 
53  The flame area enhancement factor \f$ \Xi_{sub} \f$ is expected to
54  approach:
55 
56  \f[
57  \Xi_{{sub}_{eq}} =
58  1 + max(2.2 \sqrt{b}, min(0.34 \frac{\vert \dwea{\vec{U}}
59  \vert}{{\vec{U}}^{'}}, 1.6)) \times min(\frac{n}{4}, 1)
60  \f]
61 
62 
63 SourceFiles
64  basicSubGrid.C
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #ifndef basicSubGrid_H
69 #define basicSubGrid_H
70 
71 #include "XiEqModel.H"
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 namespace Foam
76 {
77 namespace XiEqModels
78 {
79 
80 /*---------------------------------------------------------------------------*\
81  Class basicSubGrid Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class basicSubGrid
85 :
86  public XiEqModel
87 {
88  // Private data
89 
90  //- tblock
92 
93  //- Equilibrium Xi model due to turbulence
94  autoPtr<XiEqModel> XiEqModel_;
95 
96 
97  // Private Member Functions
98 
99  //- No copy construct
100  basicSubGrid(const basicSubGrid&) = delete;
101 
102  //- No copy assignment
103  void operator=(const basicSubGrid&) = delete;
104 
105 
106 public:
107 
108  //- Runtime type information
109  TypeName("basicSubGrid");
110 
111 
112  // Constructors
113 
114  //- Construct from components
116  (
117  const dictionary& XiEqProperties,
118  const psiuReactionThermo& thermo,
120  const volScalarField& Su
121  );
122 
123 
124  //- Destructor
125  virtual ~basicSubGrid();
126 
127 
128  // Member Functions
129 
130  //- Return the flame-wrinkling XiEq
131  virtual tmp<volScalarField> XiEq() const;
132 
133  //- Update properties from given dictionary
134  virtual bool read(const dictionary& XiEqProperties);
135 };
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace XiEqModels
141 } // End namespace Foam
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 #endif
146 
147 // ************************************************************************* //
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
XiEqModel.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::XiEqModels::basicSubGrid::~basicSubGrid
virtual ~basicSubGrid()
Destructor.
Su
zeroField Su
Definition: alphaSuSp.H:1
Foam::XiEqModels::basicSubGrid
Basic sub-grid obstacle flame-wrinkling enhancement factor model. Details supplied by J Puttock 2/7/0...
Definition: basicXiSubXiEq.H:83
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:66
Foam::XiEqModels::basicSubGrid::XiEq
virtual tmp< volScalarField > XiEq() const
Return the flame-wrinkling XiEq.
Foam::XiEqModels::basicSubGrid::TypeName
TypeName("basicSubGrid")
Runtime type information.
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::XiEqModels::basicSubGrid::read
virtual bool read(const dictionary &XiEqProperties)
Update properties from given dictionary.