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