diskInjection.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) 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::diskInjection
26 
27 Description
28  Disk injection model. This injects particles continuously over a disk with
29  a given number rate. The disk is characterised by a centre, and axis and an
30  inner and an outer diameter. The centre and axis and number rate are all
31  Function1-s and can vary with time.
32 
33  Note that this model only controls the number and position of injected
34  Lagrangian particles. All physical properties are specified by
35  corresponding source conditions. So the velocity/direction/angle/etc..., is
36  controlled by the velocity source condition, the size distribution by the
37  diameter source condition, and the flow rate by the number source
38  condition.
39 
40 Usage
41  \table
42  Property | Description | Required? | Default
43  centre | The centre of the disk | yes |
44  axis | The axis normal to the disk | yes |
45  diameter | The diameter of the circle | if innerDiameter \
46  and outerDiameter \
47  are not specified |
48  innerDiameter | The inner diameter of the annulus | if diameter is not \
49  specified |
50  outerDiameter | The outer diameter of the annulus | if diameter is not \
51  specified |
52  numberRate | The number of particles to \
53  inject per unit time | yes |
54  \endtable
55 
56  Example specification:
57  \verbatim
58  <LagrangianModelName>
59  {
60  type pointInjection;
61  point (1 2 3) [cm];
62  numberRate 10000;
63  }
64  \endverbatim
65 
66 See also
67  Foam::Function1s
68 
69 SourceFiles
70  diskInjection.C
71 
72 \*---------------------------------------------------------------------------*/
73 
74 #ifndef diskInjection_H
75 #define diskInjection_H
76 
77 #include "LagrangianInjection.H"
78 #include "Function1.H"
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 namespace Foam
84 {
85 namespace Lagrangian
86 {
87 
88 /*---------------------------------------------------------------------------*\
89  Class diskInjection Declaration
90 \*---------------------------------------------------------------------------*/
91 
92 class diskInjection
93 :
94  public LagrangianInjection
95 {
96 private:
97 
98  // Private Data
99 
100  //- The centre of the disk
101  autoPtr<Function1<point>> centre_;
102 
103  //- The axis of the disk
104  autoPtr<Function1<vector>> axis_;
105 
106  //- The inner diameter
107  scalar innerDiameter_;
108 
109  //- The outer diameter
110  scalar outerDiameter_;
111 
112  //- The rate at which to inject
113  autoPtr<Function1<scalar>> numberRate_;
114 
115  //- The number deferred to the next injection step
116  mutable scalar numberDeferred_;
117 
118  //- A random number generator
119  mutable restartableRandomGenerator rndGen_;
120 
121  //- The time index
122  mutable label timeIndex_;
123 
124  //- Cached radius fractions
125  mutable autoPtr<scalarField> rFracPtr_;
126 
127  //- Cached axes
128  mutable autoPtr<vectorField> axisPtr_;
129 
130  //- Cached radial vectors
131  mutable autoPtr<vectorField> radialPtr_;
132 
133 
134  // Private Member Functions
135 
136  //- Non-virtual read
137  void readCoeffs(const dictionary& modelDict);
138 
139 
140 public:
141 
142  //- Runtime type information
143  TypeName("diskInjection");
144 
145 
146  // Constructors
147 
148  //- Construct from components
150  (
151  const word& name,
152  const LagrangianMesh& mesh,
153  const dictionary& modelDict,
154  const dictionary& stateDict
155  );
156 
157 
158  // Member Functions
159 
160  // Access
161 
162  //- Access the cached radius fractions. Only valid during injection.
163  const scalarField& rFrac() const;
164 
165  //- Access the cached axes. Only valid during injection.
166  const vectorField& axis() const;
167 
168  //- Access the cached radial vectors. Only valid during injection.
169  const vectorField& radial() const;
170 
171 
172  //- Create new elements in the Lagrangian mesh
174  (
176  const LagrangianSubMesh&
177  ) const;
178 
179 
180  // IO
181 
182  //- Read dictionary
183  virtual bool read(const dictionary& modelDict);
184 
185  //- Write state
186  virtual void writeState(Ostream& os) const;
187 };
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 } // End namespace Lagrangian
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
Class containing Lagrangian geometry and topology.
const LagrangianMesh & mesh() const
The mesh.
const word & name() const
The source name.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
virtual void writeState(Ostream &os) const
Write state.
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &) const
Create new elements in the Lagrangian mesh.
diskInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
const scalarField & rFrac() const
Access the cached radius fractions. Only valid during injection.
const vectorField & radial() const
Access the cached radial vectors. Only valid during injection.
virtual bool read(const dictionary &modelDict)
Read dictionary.
TypeName("diskInjection")
Runtime type information.
const vectorField & axis() const
Access the cached axes. Only valid during injection.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Disk injection model. This injects particles continuously over a disk with a given number rate....
static dictionary stateDict(const word &name, const objectRegistry &db)
Construct and return the state dictionary for reading.
Definition: stateModel.C:137
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59