XiCorrModel.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) 2024-2025 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::XiCorrModel
26 
27 Description
28  Base class for ignition kernel flame wrinkling \c Xi correction
29 
30 Usage
31  Example usage:
32  \verbatim
33  XiCorr
34  {
35  type <type>;
36  cellZone all;
37  bMin 0.05;
38  }
39  \endverbatim
40 
41  Where:
42  \table
43  Property | Description | Required | Default value
44  cellZone | Correction cellZone | yes |
45  bMin | Min b below which no correction | no | 0.01
46  \endtable
47 
48 SourceFiles
49  XiCorrModel.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef XiCorrModel_H
54 #define XiCorrModel_H
55 
56 #include "fvCellZone.H"
57 #include "volFieldsFwd.H"
58 #include "dimensionedTypes.H"
59 #include "runTimeSelectionTables.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 namespace Foam
64 {
65 
66 /*---------------------------------------------------------------------------*\
67  Class XiCorrModel Declaration
68 \*---------------------------------------------------------------------------*/
69 
70 class XiCorrModel
71 {
72  // Private Data
73 
74  //- The Xi correction cellZone
75  fvCellZone zone_;
76 
77  //- Minimum b below which the Xi correction is not applied
78  // Defaults to 0.01
79  dimensionedScalar bMin_;
80 
81 
82 protected:
83 
84  //- Update coefficients from given dictionary
85  virtual bool readCoeffs(const dictionary& dict);
86 
87 
88 public:
89 
90  //- Runtime type information
91  TypeName("XiCorrModel");
92 
93 
94  // Declare run-time constructor selection table
95 
97  (
98  autoPtr,
100  dictionary,
101  (const fvMesh& mesh, const dictionary& dict),
102  (mesh, dict)
103  );
104 
105 
106  // Constructors
107 
108  //- Construct from mesh and dictionary
109  XiCorrModel(const fvMesh& mesh, const dictionary& dict);
110 
111  //- Disallow default bitwise copy construction
112  XiCorrModel(const XiCorrModel&) = delete;
113 
114 
115  // Selectors
116 
117  //- Return a reference to the selected XiCorrModel model
119  (
120  const fvMesh& mesh,
121  const dictionary& dict
122  );
123 
124 
125  //- Destructor
126  virtual ~XiCorrModel();
127 
128 
129  // Member Functions
130 
131  //- Return the area of the ignition kernel
132  // calculated from the given kernel volume
133  virtual dimensionedScalar Ak(const dimensionedScalar& Vk) const = 0;
134 
135  virtual void XiCorr
136  (
137  volScalarField& Xi,
138  const volScalarField& b,
139  const volScalarField& mgb
140  ) const;
141 
142 
143  // Mesh motion
144 
145  //- Update topology using the given map
146  virtual void topoChange(const polyTopoChangeMap&);
147 
148  //- Update from another mesh using the given map
149  virtual void mapMesh(const polyMeshMap&);
150 
151  //- Redistribute or update using the given distribution map
152  virtual void distribute(const polyDistributionMap&);
153 
154  //- Update for mesh motion
155  virtual bool movePoints();
156 
157 
158  // IO
159 
160  //- Update properties from the given XiProperties dictionary
161  bool read(const dictionary& XiProperties);
162 
163 
164  // Member Operators
165 
166  //- Disallow default bitwise assignment
167  void operator=(const XiCorrModel&) = delete;
168 };
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #endif
178 
179 // ************************************************************************* //
Generic GeometricField class.
Base class for ignition kernel flame wrinkling Xi correction.
Definition: XiCorrModel.H:85
virtual bool movePoints()
Update for mesh motion.
Definition: XiCorrModel.C:139
virtual ~XiCorrModel()
Destructor.
Definition: XiCorrModel.C:64
virtual void XiCorr(volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
Definition: XiCorrModel.C:71
virtual dimensionedScalar Ak(const dimensionedScalar &Vk) const =0
Return the area of the ignition kernel.
TypeName("XiCorrModel")
Runtime type information.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: XiCorrModel.C:116
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: XiCorrModel.C:131
void operator=(const XiCorrModel &)=delete
Disallow default bitwise assignment.
static autoPtr< XiCorrModel > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected XiCorrModel model.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: XiCorrModel.C:124
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
Definition: XiCorrModel.C:40
declareRunTimeSelectionTable(autoPtr, XiCorrModel, dictionary,(const fvMesh &mesh, const dictionary &dict),(mesh, dict))
bool read(const dictionary &XiProperties)
Update properties from the given XiProperties dictionary.
Definition: XiCorrModel.C:146
XiCorrModel(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
Definition: XiCorrModel.C:50
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField & b
Definition: createFields.H:27
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Macros to ease declaration of run-time selection tables.
dictionary dict