algebraic.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-2023 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::XiModels::algebraic
26 
27 Description
28  Simple algebraic model for Xi based on Gulders correlation
29  with a linear correction function to give a plausible profile for Xi.
30  See report TR/HGW/10 for details on the Weller two equations model.
31  See \link XiModel.H \endlink for more details on flame wrinkling modelling.
32 
33 SourceFiles
34  algebraic.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef algebraic_H
39 #define algebraic_H
40 
41 #include "XiModel.H"
42 #include "XiEqModel.H"
43 #include "XiGModel.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace XiModels
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class algebraic Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class algebraic
57 :
58  public XiModel
59 {
60  // Private Data
61 
62  scalar XiShapeCoef;
63 
64  autoPtr<XiEqModel> XiEqModel_;
65  autoPtr<XiGModel> XiGModel_;
66 
67 
68 public:
69 
70  //- Runtime type information
71  TypeName("algebraic");
72 
73 
74  // Constructors
75 
76  //- Construct from components
77  algebraic
78  (
79  const dictionary& XiProperties,
82  const volScalarField& Su,
83  const volScalarField& rho,
84  const volScalarField& b,
85  const surfaceScalarField& phi
86  );
87 
88  //- Disallow default bitwise copy construction
89  algebraic(const algebraic&) = delete;
90 
91 
92  //- Destructor
93  virtual ~algebraic();
94 
95 
96  // Member Functions
97 
98  //- Return the flame diffusivity
99  virtual tmp<volScalarField> Db() const;
100 
101  //- Correct the flame-wrinkling Xi
102  virtual void correct();
103 
104  //- Update properties from given dictionary
105  virtual bool read(const dictionary& XiProperties);
106 
107  //- Write fields of the XiEq model
108  virtual void writeFields()
109  {
110  XiEqModel_().writeFields();
111  }
112 
113 
114  // Member Operators
115 
116  //- Disallow default bitwise assignment
117  void operator=(const algebraic&) = delete;
118 };
119 
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 } // End namespace XiModels
124 } // End namespace Foam
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 #endif
129 
130 // ************************************************************************* //
Generic GeometricField class.
Base-class for all Xi models used by the b-Xi combustion model. See Technical Report SH/RE/01R for de...
Definition: XiModel.H:109
Simple algebraic model for Xi based on Gulders correlation with a linear correction function to give ...
Definition: algebraic.H:58
virtual void writeFields()
Write fields of the XiEq model.
Definition: algebraic.H:107
virtual ~algebraic()
Destructor.
Definition: algebraic.C:63
void operator=(const algebraic &)=delete
Disallow default bitwise assignment.
virtual void correct()
Correct the flame-wrinkling Xi.
Definition: algebraic.C:75
virtual bool read(const dictionary &XiProperties)
Update properties from given dictionary.
Definition: algebraic.C:88
TypeName("algebraic")
Runtime type information.
virtual tmp< volScalarField > Db() const
Return the flame diffusivity.
Definition: algebraic.C:69
algebraic(const dictionary &XiProperties, const psiuMulticomponentThermo &thermo, const compressible::RASModel &turbulence, const volScalarField &Su, const volScalarField &rho, const volScalarField &b, const surfaceScalarField &phi)
Construct from components.
Definition: algebraic.C:44
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Base-class for combustion fluid thermodynamic properties based on compressibility.
A class for managing temporary objects.
Definition: tmp.H:55
volScalarField & b
Definition: createFields.H:27
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
fluidMulticomponentThermo & thermo
Definition: createFields.H:31