pointInjection.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::pointInjection
26 
27 Description
28  Point injection model. This injects particles continuously at a point with
29  a given number rate. The point and the number rate are both Function1-s and
30  can both vary with time.
31 
32  Note that this model only controls the number and position of injected
33  Lagrangian particles. All physical properties are specified by
34  corresponding source conditions. So the velocity/direction/angle/etc..., is
35  controlled by the velocity source condition, the size distribution by the
36  diameter source condition, and the flow rate by the number source
37  condition.
38 
39 Usage
40  \table
41  Property | Description | Required? | Default
42  point | The position at which to inject | yes |
43  numberRate | The number of particles to inject per unit time | yes |
44  \endtable
45 
46  Example specification:
47  \verbatim
48  <LagrangianModelName>
49  {
50  type pointInjection;
51  point (1 2 3) [cm];
52  numberRate 10000;
53  }
54  \endverbatim
55 
56 See also
57  Foam::Function1s
58 
59 SourceFiles
60  pointInjection.C
61 
62 \*---------------------------------------------------------------------------*/
63 
64 #ifndef pointInjection_H
65 #define pointInjection_H
66 
67 #include "LagrangianInjection.H"
68 #include "Function1.H"
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 namespace Foam
74 {
75 namespace Lagrangian
76 {
77 
78 /*---------------------------------------------------------------------------*\
79  Class pointInjection Declaration
80 \*---------------------------------------------------------------------------*/
81 
82 class pointInjection
83 :
84  public LagrangianInjection
85 {
86 private:
87 
88  // Private Data
89 
90  //- The point at which to inject
91  autoPtr<Function1<point>> point_;
92 
93  //- The rate at which to inject
94  autoPtr<Function1<scalar>> numberRate_;
95 
96  //- The number deferred to the next injection step
97  mutable scalar numberDeferred_;
98 
99  //- A random number generator
100  mutable restartableRandomGenerator rndGen_;
101 
102  //- The time index
103  mutable label timeIndex_;
104 
105  //- Where is the injection?
106  enum class injectionLocation
107  {
108  unset,
109  fixedPointOnThisProcessor,
110  fixedPointOnAnotherProcessor,
111  multiplePoints
112  } injectionLocation_;
113 
114  //- The coordinates of the point
115  barycentric coordinates_;
116 
117  //- The cell containing the point
118  label celli_;
119 
120  //- The face adjacent to the tet containing the point
121  label facei_;
122 
123  //- The face-tri adjacent to the tet containing the point
124  label faceTrii_;
125 
126 
127  // Private Member Functions
128 
129  //- Non-virtual read
130  void readCoeffs(const dictionary& modelDict);
131 
132 
133 public:
134 
135  //- Runtime type information
136  TypeName("pointInjection");
137 
138 
139  // Constructors
140 
141  //- Construct from components
143  (
144  const word& name,
145  const LagrangianMesh& mesh,
146  const dictionary& modelDict,
147  const dictionary& stateDict
148  );
149 
150 
151  // Member Functions
152 
153  //- Correct the LagrangianModel
154  virtual void correct();
155 
156  //- Create new elements in the Lagrangian mesh
158  (
160  const LagrangianSubMesh&
161  ) const;
162 
163 
164  // Mesh changes
165 
166  //- Update topology using the given map
167  virtual void topoChange(const polyTopoChangeMap&);
168 
169  //- Update from another mesh using the given map
170  virtual void mapMesh(const polyMeshMap&);
171 
172  //- Redistribute or update using the given distribution map
173  virtual void distribute(const polyDistributionMap&);
174 
175 
176  // IO
177 
178  //- Read dictionary
179  virtual bool read(const dictionary& modelDict);
180 
181  //- Write state
182  virtual void writeState(Ostream& os) const;
183 };
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 } // End namespace Lagrangian
189 } // End namespace Foam
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #endif
194 
195 // ************************************************************************* //
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...
pointInjection(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
virtual void correct()
Correct the LagrangianModel.
virtual void writeState(Ostream &os) const
Write state.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
LagrangianSubMesh modify(LagrangianMesh &mesh, const LagrangianSubMesh &) const
Create new elements in the Lagrangian mesh.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &modelDict)
Read dictionary.
TypeName("pointInjection")
Runtime type information.
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
Point injection model. This injects particles continuously at a point with a given number rate....
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.
Random number generator with the additional ability to go back to an earlier stored state....
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