sphericalXiCorr.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::XiCorrModels::spherical
26 
27 Description
28  Spherical ignition kernel flame-wrinkling correction model
29 
30 Usage
31  Example usage:
32  \verbatim
33  XiCorr
34  {
35  type spherical;
36  cellZone all;
37  sphereFraction 1;
38  }
39  \endverbatim
40 
41  Where:
42  \table
43  Property | Description | Required | Default value
44  cellZone | Correction cellZone | yes |
45  sphereFraction | Kernel sphere fraction | no | 1
46  bMin | Min b below which no correction | no | 0.01
47  \endtable
48 
49 SourceFiles
50  spherical.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef sphericalXiCorr_H
55 #define sphericalXiCorr_H
56 
57 #include "XiCorrModel.H"
58 #include "dimensionedTypes.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace XiCorrModels
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class spherical Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class spherical
72 :
73  public XiCorrModel
74 {
75  // Private Data
76 
77  //- Kernel sphere fraction
78  // If the ignition is adjacent to n symmetry planes: 1/2^n
79  dimensionedScalar sphereFraction_;
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("spherical");
92 
93 
94  // Constructors
95 
96  //- Construct from mesh and dictionary
97  spherical(const fvMesh& mesh, const dictionary& dict);
98 
99  //- Disallow default bitwise copy construction
100  spherical(const spherical&) = delete;
101 
102 
103  //- Destructor
104  virtual ~spherical();
105 
106 
107  // Member Functions
108 
109  //- Return the area of the ignition kernel
110  // calculated from the given kernel volume
111  virtual dimensionedScalar Ak(const dimensionedScalar& Vk) const;
112 
113 
114  // Member Operators
115 
116  //- Disallow default bitwise assignment
117  void operator=(const spherical&) = delete;
118 };
119 
120 
121 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122 
123 } // End namespace XiCorrModels
124 } // End namespace Foam
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 #endif
129 
130 // ************************************************************************* //
Spherical ignition kernel flame-wrinkling correction model.
virtual dimensionedScalar Ak(const dimensionedScalar &Vk) const
Return the area of the ignition kernel.
TypeName("spherical")
Runtime type information.
spherical(const fvMesh &mesh, const dictionary &dict)
Construct from mesh and dictionary.
virtual bool readCoeffs(const dictionary &dict)
Update coefficients from given dictionary.
void operator=(const spherical &)=delete
Disallow default bitwise assignment.
virtual ~spherical()
Destructor.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict