injectedParticle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2024 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::injectedParticle
28 
29 Description
30  Primarily stores particle properties so that it can be injected at a later
31  time. Note that this stores its own local position as opposed to the
32  base particle class barycentric coordinates since the particle is not
33  (usually) attached to a mesh, and instead used for post-processing.
34 
35 SourceFiles
36  injectedParticle.C
37  injectedParticleIO.C
38 
39 SeeAlso
40  Foam::functionObjects::extractEulerianParticles
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef Foam_injectedParticle_H
45 #define Foam_injectedParticle_H
46 
47 #include "particle.H"
48 #include "IOstream.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward Declarations
56 class injectedParticle;
57 Ostream& operator<<(Ostream&, const injectedParticle&);
58 
59 /*---------------------------------------------------------------------------*\
60  Class injectedParticle Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 class injectedParticle
64 :
65  public particle
66 {
67 protected:
68 
69  // Protected Data
70 
71  // Particle properties
72 
73  //- Position
75 
76  //- Tag
77  label tag_;
78 
79  //- Start of injection [s]
80  scalar soi_;
81 
82  //- Diameter [m]
83  scalar d_;
84 
85  //- Velocity [m/s]
87 
88 
89 public:
90 
91  // Static Data Members
92 
93  //- Size in bytes of the fields
94  static const std::size_t sizeofFields;
95 
96  //- Runtime type information
97  TypeName("injectedParticle");
98 
99  //- String representation of properties
101  (
102  particle,
103  " tag"
104  + " soi"
105  + " d"
106  + " (Ux Uy Uz)";
107  );
108 
109 
110  // Constructors
111 
112  //- Construct from a position and a cell.
113  // Searches for the rest of the required topology.
114  // Other properties are zero initialised.
115  inline injectedParticle
116  (
117  const polyMesh& mesh,
118  const vector& position,
119  const label celli = -1
120  );
121 
122  //- Construct from components
123  inline injectedParticle
124  (
125  const polyMesh& mesh,
126  const vector& position,
127  const label tag,
128  const scalar soi,
129  const scalar d,
130  const vector& U,
131  const bool doLocate = true
132  );
133 
134  //- Construct from Istream
136  (
137  const polyMesh& mesh,
138  Istream& is,
139  bool readFields = true,
140  bool newFormat = true
141  );
142 
143  //- Construct as a copy
145 
146  //- Construct as a copy
148 
149  //- Return a (basic particle) clone
150  virtual autoPtr<particle> clone() const
151  {
152  return particle::Clone(*this);
153  }
154 
155  //- Construct and return a (basic particle) clone
156  virtual autoPtr<particle> clone(const polyMesh& mesh) const
157  {
158  return particle::Clone(*this, mesh);
159  }
160 
161  //- Factory class to read-construct particles (for parallel transfer)
162  class iNew
163  {
164  const polyMesh& mesh_;
165 
166  public:
167 
168  iNew(const polyMesh& mesh)
169  :
170  mesh_(mesh)
171  {}
172 
173  autoPtr<injectedParticle> operator()(Istream& is) const
174  {
175  return autoPtr<injectedParticle>::New(mesh_, is, true);
176  }
177  };
178 
179 
180  // Member Functions
181 
182  // Access
183 
184  //- Return const access to the tag
185  label tag() const noexcept { return tag_; }
186 
187  //- Return const access to the start of injection
188  scalar soi() const noexcept { return soi_; }
189 
190  //- Return const access to diameter
191  scalar d() const noexcept { return d_; }
192 
193  //- Return const access to velocity
194  const vector& U() const noexcept { return U_; }
195 
197  // Edit
198 
199  //- Return the tag
200  label& tag() noexcept { return tag_; }
202  //- Return the start of injection
203  scalar& soi() noexcept { return soi_; }
204 
205  //- Return access to diameter
206  scalar& d() noexcept { return d_; }
207 
208  //- Return access to velocity
209  vector& U() noexcept { return U_; }
210 
211 
212  // I-O
213 
214  //- Read fields
216 
217  //- Write individual parcel properties to stream
218  void writeProperties
219  (
221  const wordRes& filters,
222  const word& delim,
223  const bool namesOnly
224  ) const;
226  //- Write fields
227  static void writeFields(const Cloud<injectedParticle>& c);
228 
229  //- Read particle fields as objects from the obr registry
230  static void readObjects
231  (
233  const objectRegistry& obr
234  );
235 
236  //- Write particle fields as objects into the obr registry
237  static void writeObjects
238  (
239  const Cloud<injectedParticle>& c,
240  objectRegistry& obr
241  );
242 
243  //- Write the particle position and cell
244  // Note: This uses the local particle position, and bypasses the
245  // barycentric description
246  virtual void writePosition(Ostream&) const;
247 
249  // Ostream Operator
250 
251  friend Ostream& operator<<
252  (
254  const injectedParticle&
255  );
256 };
257 
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 } // End namespace Foam
262 
263 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264 
265 #include "injectedParticleI.H"
266 
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268 
269 #endif
270 
271 // ************************************************************************* //
TypeName("injectedParticle")
Runtime type information.
injectedParticle(const polyMesh &mesh, const vector &position, const label celli=-1)
Construct from a position and a cell.
autoPtr< injectedParticle > operator()(Istream &is) const
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual parcel properties to stream.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const vector & U() const noexcept
Return const access to velocity.
Primarily stores particle properties so that it can be injected at a later time. Note that this store...
AddToPropertyList(particle, " tag"+" soi"+" d"+" (Ux Uy Uz)";)
String representation of properties.
Base particle class.
Definition: particle.H:69
static const std::size_t sizeofFields
Size in bytes of the fields.
virtual void writePosition(Ostream &) const
Write the particle position and cell.
vector U_
Velocity [m/s].
static void writeObjects(const Cloud< injectedParticle > &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
virtual autoPtr< particle > clone() const
Return a (basic particle) clone.
OBJstream os(runTime.globalPath()/outputName)
point position_
Position.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
scalar d_
Diameter [m].
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:110
const dimensionedScalar c
Speed of light in a vacuum.
label tag() const noexcept
Return const access to the tag.
static void readObjects(Cloud< injectedParticle > &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
volScalarField & p
scalar d() const noexcept
Return const access to diameter.
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
scalar soi() const noexcept
Return const access to the start of injection.
iNew(const polyMesh &mesh)
vector position() const
Return current particle position.
Definition: particleI.H:283
static void writeFields(const Cloud< injectedParticle > &c)
Write fields.
Namespace for OpenFOAM.
scalar soi_
Start of injection [s].
static void readFields(Cloud< injectedParticle > &c)
Read fields.