constantbXiIgnition.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::fv::constantbXiIgnition
26 
27 Description
28  Simple constant rate ignition model for the Weller b-Xi combustion models
29 
30 Usage
31  Example usage:
32  \verbatim
33  constantbXiIgnition1
34  {
35  type constantbXiIgnition;
36 
37  cellZone ignition;
38 
39  start 0;
40  duration 0.003;
41  strength 2;
42 
43  XiCorr
44  {
45  type spherical;
46  cellZone all;
47  }
48  }
49  \endverbatim
50 
51  Where:
52  \table
53  Property | Description | Required | Default value
54  cellZone | Correction cellZone | yes |
55  start | Ignition start time | yes |
56  duration | Ignition duration | yes |
57  strength | Ignition strength [1/s] | yes |
58  XiCorr | Flame-wrinkling correction | yes |
59  \endtable
60 
61 See also
62  Foam::XiCorrModel
63  Foam::XiCorrModels::planar
64  Foam::XiCorrModels::cylindrical
65  Foam::XiCorrModels::spherical
66 
67 SourceFiles
68  constantbXiIgnition.C
69 
70 \*---------------------------------------------------------------------------*/
71 
72 #ifndef constantbXiIgnition_H
73 #define constantbXiIgnition_H
74 
75 #include "bXiIgnition.H"
76 #include "XiCorrModel.H"
77 
78 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 namespace fv
83 {
84 
85 /*---------------------------------------------------------------------------*\
86  Class constantbXiIgnition Declaration
87 \*---------------------------------------------------------------------------*/
88 
89 class constantbXiIgnition
90 :
91  public bXiIgnition
92 {
93 
94 protected:
95 
96  // Protected Data
97 
98  //- The ignition cellZone
99  fvCellZone zone_;
100 
101  autoPtr<XiCorrModel> XiCorrModel_;
102 
103  //- Ignition start time
105 
106  //- Ignition duration
108 
109  //- Ignition strength
111 
112 
113 private:
114 
115  // Private Member Functions
116 
117  //- Non-virtual read
118  void readCoeffs(const dictionary& dict);
119 
120 
121 public:
122 
123  //- Runtime type information
124  TypeName("constantbXiIgnition");
125 
126 
127  // Constructors
128 
129  //- Construct from explicit source name and mesh
131  (
132  const word& name,
133  const word& modelType,
134  const fvMesh& mesh,
135  const dictionary& dict
136  );
137 
138  //- Disallow default bitwise copy construction
140  (
141  const constantbXiIgnition&
142  ) = delete;
143 
144 
145  // Member Functions
146 
147  // Checks
148 
149  //- Return true during the ignition duration
150  virtual bool igniting() const;
151 
152  //- Return true during the combustion duration
153  virtual bool ignited() const;
154 
155 
156  // Add explicit and implicit contributions to compressible equation
157 
158  //- Add ignition contribution to b equation
159  virtual void addSup
160  (
161  const volScalarField& rho,
162  const volScalarField& b,
163  fvMatrix<scalar>& eqn
164  ) const;
165 
166  virtual void XiCorr
167  (
169  const volScalarField& b,
170  const volScalarField& mgb
171  ) const;
172 
173 
174  // Mesh motion
175 
176  //- Update topology using the given map
177  virtual void topoChange(const polyTopoChangeMap&);
178 
179  //- Update from another mesh using the given map
180  virtual void mapMesh(const polyMeshMap&);
181 
182  //- Redistribute or update using the given distribution map
183  virtual void distribute(const polyDistributionMap&);
184 
185  //- Update for mesh motion
186  virtual bool movePoints();
187 
188 
189  // IO
190 
191  //- Read source dictionary
192  virtual bool read(const dictionary& dict);
193 
194 
195  // Member Operators
196 
197  //- Disallow default bitwise assignment
198  void operator=(const constantbXiIgnition&) = delete;
199 };
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace fv
205 } // End namespace Foam
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
Generic GeometricField class.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Simple constant rate ignition model for the Weller b-Xi combustion models.
virtual bool movePoints()
Update for mesh motion.
dimensionedScalar strength_
Ignition strength.
virtual void XiCorr(volScalarField &Xi, const volScalarField &b, const volScalarField &mgb) const
TypeName("constantbXiIgnition")
Runtime type information.
autoPtr< XiCorrModel > XiCorrModel_
fvCellZone zone_
The ignition cellZone.
dimensionedScalar duration_
Ignition duration.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
constantbXiIgnition(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual bool ignited() const
Return true during the combustion duration.
dimensionedScalar start_
Ignition start time.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual void addSup(const volScalarField &rho, const volScalarField &b, fvMatrix< scalar > &eqn) const
Add ignition contribution to b equation.
void operator=(const constantbXiIgnition &)=delete
Disallow default bitwise assignment.
virtual bool igniting() const
Return true during the ignition duration.
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.
A class for handling words, derived from string.
Definition: word.H:62
volScalarField & b
Definition: createFields.H:27
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict