particle.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) 2011-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::particle
26 
27 Description
28  Base particle class
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef particle_H
33 #define particle_H
34 
35 #include "vector.H"
36 #include "barycentric.H"
37 #include "barycentricTensor.H"
38 #include "IDLList.H"
39 #include "pointField.H"
40 #include "faceList.H"
41 #include "OFstream.H"
42 #include "tetPointRef.H"
43 #include "FixedList.H"
45 #include "particleMacros.H"
46 #include "transformer.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // Forward declaration of classes
54 class particle;
55 
56 class polyPatch;
57 
58 class cyclicPolyPatch;
59 class processorPolyPatch;
60 class symmetryPlanePolyPatch;
61 class symmetryPolyPatch;
62 class wallPolyPatch;
63 class wedgePolyPatch;
64 
65 // Forward declaration of friend functions and operators
66 
67 Ostream& operator<<
68 (
69  Ostream&,
70  const particle&
71 );
72 
73 bool operator==(const particle&, const particle&);
74 
75 bool operator!=(const particle&, const particle&);
76 
77 /*---------------------------------------------------------------------------*\
78  Class Particle Declaration
79 \*---------------------------------------------------------------------------*/
80 
81 class particle
82 :
83  public IDLList<particle>::link
84 {
85  // Private member data
86 
87  //- Size in bytes of the position data
88  static const std::size_t sizeofPosition_;
89 
90  //- Size in bytes of the fields
91  static const std::size_t sizeofFields_;
92 
93 
94 public:
95 
96  class trackingData
97  {
98  public:
99 
100  // Public data
101 
102  //- Reference to the mesh
103  const polyMesh& mesh;
104 
105  //- Flag to indicate whether to keep particle (false = delete)
106  bool keepParticle;
107 
108  //- Processor to send the particle to. -1 indicates that this
109  // particle is not to be transferred.
111 
112  //- Patch from which to send the particle
114 
115  //- Patch to which to send the particle
117 
118  //- Patch face to which to send the particle
120 
121  //- Position to which to send
123 
124  //- Number of boundary hits that occurred during locate executions
125  // following (non-conformal) patch transfers. For reporting.
127 
128 
129  // Constructor
130  template <class TrackCloudType>
131  trackingData(const TrackCloudType& cloud)
132  :
133  mesh(cloud.pMesh()),
134  keepParticle(false),
135  sendToProc(-1),
136  sendFromPatch(-1),
137  sendToPatch(-1),
138  sendToPatchFace(-1),
139  sendToPosition(vector::uniform(NaN)),
141  (
142  mesh.boundaryMesh().size()
143  - mesh.globalData().processorPatches().size(),
144  0
145  )
146  {}
147  };
148 
149 
150 private:
151 
152  // Private Data
153 
154  //- Coordinates of particle
155  barycentric coordinates_;
156 
157  //- Index of the cell it is in
158  label celli_;
159 
160  //- Index of the face that owns the decomposed tet that the
161  // particle is in
162  label tetFacei_;
163 
164  //- Index of the point on the face that defines the decomposed
165  // tet that the particle is in. Relative to the face base
166  // point.
167  label tetPti_;
168 
169  //- Face index if the particle is on a face otherwise -1
170  label facei_;
171 
172  //- Fraction of time-step completed
173  scalar stepFraction_;
174 
175  //- The step fraction less than the maximum reached so far. See
176  // tracking.H for details.
177  scalar stepFractionBehind_;
178 
179  //- The number of tracks carried out that ended in a step fraction less
180  // than the maximum reached so far. See tracking.H for details.
181  // maxNTracksBehind_.
182  label nTracksBehind_;
183 
184  //- Originating processor id
185  label origProc_;
186 
187  //- Local particle id on originating processor
188  label origId_;
189 
190 
191 public:
192 
193  // Static Data Members
194 
195  //- Runtime type information
196  TypeName("particle");
197 
198  //- String representation of properties
200  (
201  "(coordinatesa coordinatesb coordinatesc coordinatesd) "
202  "celli tetFacei tetPti facei stepFraction "
203  "behind nBehind origProc origId"
204  );
205 
206  //- Cumulative particle counter - used to provide unique ID
207  static label particleCount;
208 
209  //- Particles are not instantaneous by default
210  static const bool instantaneous = false;
211 
212 
213  // Constructors
214 
215  //- Construct from components
216  particle
217  (
218  const polyMesh& mesh,
219  const barycentric& coordinates,
220  const label celli,
221  const label tetFacei,
222  const label tetPti,
223  const label facei
224  );
225 
226  //- Construct from a position and a cell, searching for the rest of the
227  // required topology
228  particle
229  (
230  const polyMesh& mesh,
231  const vector& position,
232  const label celli,
233  label& nLocateBoundaryHits
234  );
235 
236  //- Construct from Istream
237  particle(Istream&, bool readFields = true);
238 
239  //- Construct as a copy
240  particle(const particle& p);
241 
242  //- Construct and return a clone
243  virtual autoPtr<particle> clone() const
244  {
245  return autoPtr<particle>(new particle(*this));
246  }
247 
248  //- Construct from Istream and return
249  static autoPtr<particle> New(Istream& is)
250  {
251  return autoPtr<particle>(new particle(is));
252  }
253 
254 
255  //- Destructor
256  virtual ~particle()
257  {}
258 
259 
260  // Member Functions
261 
262  // Access
263 
264  //- Get unique particle creation id
265  inline label getNewParticleIndex() const;
266 
267  //- Return current particle coordinates
268  inline const barycentric& coordinates() const;
269 
270  //- Return current cell particle is in
271  inline label cell() const;
272 
273  //- Return current tet face particle is in
274  inline label tetFace() const;
275 
276  //- Return current tet face particle is in
277  inline label tetPt() const;
278 
279  //- Return current face particle is on otherwise -1
280  inline label face() const;
281 
282  //- Return current face particle is on otherwise -1
283  inline label& face();
284 
285  //- Return the fraction of time-step completed
286  inline scalar stepFraction() const;
287 
288  //- Return the fraction of time-step completed
289  inline scalar& stepFraction();
290 
291  //- Return the originating processor ID
292  inline label origProc() const;
293 
294  //- Return the originating processor ID
295  inline label& origProc();
296 
297  //- Return the particle ID on the originating processor
298  inline label origId() const;
299 
300  //- Return the particle ID on the originating processor
301  inline label& origId();
302 
303 
304  // Check
305 
306  //- Return the indices of the current tet that the
307  // particle occupies.
308  inline tetIndices currentTetIndices(const polyMesh& mesh) const;
309 
310  //- Return the normal of the tri on tetFacei_ for the
311  // current tet.
312  inline vector normal(const polyMesh& mesh) const;
313 
314  //- Is the particle on a face?
315  inline bool onFace() const;
316 
317  //- Is the particle on an internal face?
318  inline bool onInternalFace(const polyMesh& mesh) const;
319 
320  //- Is the particle on a boundary face?
321  inline bool onBoundaryFace(const polyMesh& mesh) const;
322 
323  //- Return the index of patch that the particle is on
324  inline label patch(const polyMesh& mesh) const;
325 
326  //- Return current particle position
327  inline vector position(const polyMesh& mesh) const;
328 
329 
330  // Track
331 
332  //- Set the step fraction and clear the behind data in preparation
333  // for a new track
334  inline void reset(const scalar stepFraction);
335 
336  //- Locate the particle at the given position
337  bool locate
338  (
339  const polyMesh& mesh,
340  const vector& position,
341  label celli
342  );
343 
344  //- Track along the displacement for a given fraction of the overall
345  // step. End when the track is complete, or when a boundary is hit.
346  // On exit, stepFraction_ will have been incremented to the current
347  // position, and facei_ will be set to the index of the boundary
348  // face that was hit, or -1 if the track completed within a cell.
349  // The proportion of the displacement still to be completed is
350  // returned.
351  scalar track
352  (
353  const polyMesh& mesh,
354  const vector& displacement,
355  const scalar fraction
356  );
357 
358  //- As particle::track, but stops when a new cell is reached.
359  scalar trackToCell
360  (
361  const polyMesh& mesh,
362  const vector& displacement,
363  const scalar fraction
364  );
365 
366  //- As particle::track, but stops when a face is hit.
367  scalar trackToFace
368  (
369  const polyMesh& mesh,
370  const vector& displacement,
371  const scalar fraction
372  );
373 
374  //- Hit the current face. If the current face is internal than this
375  // crosses into the next cell. If it is a boundary face then this
376  // will interact the particle with the relevant patch.
377  template<class TrackCloudType>
378  void hitFace
379  (
380  const vector& displacement,
381  const scalar fraction,
382  TrackCloudType& cloud,
383  trackingData& td
384  );
385 
386  //- Convenience function. Combines trackToFace and hitFace
387  template<class TrackCloudType>
388  scalar trackToAndHitFace
389  (
390  const vector& displacement,
391  const scalar fraction,
392  TrackCloudType& cloud,
393  trackingData& td
394  );
395 
396  //- Get the displacement from the mesh centre. Used to correct the
397  // particle position in cases with reduced dimensionality. Returns
398  // a zero vector for three-dimensional cases.
400 
401 
402  // Patch data
403 
404  //- Get the normal and displacement of the current patch location
405  inline void patchData
406  (
407  const polyMesh& mesh,
408  vector& normal,
409  vector& displacement
410  ) const;
411 
412 
413  // Transformations
414 
415  //- Transform the physical properties of the particle
416  // according to the given transformation tensor
417  virtual void transformProperties(const transformer&);
418 
419 
420  // Transfers
421 
422  //- Make changes prior to a parallel transfer. Runs either
423  // processor or nonConformalCyclic variant below.
424  template<class TrackCloudType>
425  void prepareForParallelTransfer(TrackCloudType&, trackingData&);
426 
427  //- Make changes following a parallel transfer. Runs either
428  // processor or nonConformalCyclic variant below.
429  template<class TrackCloudType>
430  void correctAfterParallelTransfer(TrackCloudType&, trackingData&);
431 
432  //- Make changes prior to a transfer across a processor boundary.
433  // Stores the local patch face index (in facei_) so that the mesh
434  // face index can be determined on the other side.
436 
437  //- Make changes following a transfer across a processor boundary.
438  // Converts the stored patch index to a mesh index. Accounts for
439  // the receiving face being reversed relative to the sending face.
441 
442  //- Make changes prior to a transfer across a non conformal cyclic
443  // boundary. Stores the receiving patch face (in facei_). Breaks
444  // the topology and stores the cartesian position.
446  (
447  const polyMesh& mesh,
448  const label sendToPatch,
449  const label sendToPatchFace,
450  const vector& sendToPosition
451  );
452 
453  //- Make changes following a transfer across a non conformal cyclic
454  // boundary. Locates the particle using the stored face index and
455  // cartesian position.
457  (
458  const polyMesh& mesh,
459  const label sendToPatch,
460  labelList& patchNLocateBoundaryHits
461  );
462 
463 
464  // Patch interactions
465 
466  //- Overridable function to handle the particle hitting a patch.
467  // Executed before other patch-hitting functions.
468  template<class TrackCloudType>
469  bool hitPatch(TrackCloudType&, trackingData&);
470 
471  //- Overridable function to handle the particle hitting a wedgePatch
472  template<class TrackCloudType>
473  void hitWedgePatch(TrackCloudType&, trackingData&);
474 
475  //- Overridable function to handle the particle hitting a
476  // symmetryPlanePatch
477  template<class TrackCloudType>
478  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
479 
480  //- Overridable function to handle the particle hitting a
481  // symmetryPatch
482  template<class TrackCloudType>
483  void hitSymmetryPatch(TrackCloudType&, trackingData&);
484 
485  //- Overridable function to handle the particle hitting a
486  // cyclicPatch
487  template<class TrackCloudType>
488  void hitCyclicPatch(TrackCloudType&, trackingData&);
489 
490  //- Overridable function to handle the particle hitting an
491  // nonConformalCyclicPolyPatch
492  template<class TrackCloudType>
494  (
495  const vector& displacement,
496  const scalar fraction,
497  const label patchi,
498  TrackCloudType& cloud,
499  trackingData& td
500  );
501 
502  //- Overridable function to handle the particle hitting a
503  // processorPatch
504  template<class TrackCloudType>
505  void hitProcessorPatch(TrackCloudType&, trackingData&);
506 
507  //- Overridable function to handle the particle hitting a wallPatch
508  template<class TrackCloudType>
509  void hitWallPatch(TrackCloudType&, trackingData&);
510 
511  //- Overridable function to handle the particle hitting a basic
512  // patch. Fall-through for the above.
513  template<class TrackCloudType>
514  void hitBasicPatch(TrackCloudType&, trackingData&);
515 
516 
517  // Interaction list referral
518 
519  //- Break the topology and store the cartesian position so that the
520  // particle can be referred.
522  (
523  const polyMesh& mesh,
524  const transformer& transform
525  );
526 
527  //- Correct the topology after referral. Locates the particle
528  // relative to a nearby cell/tet. The particle may end up outside
529  // this cell/tet and cannot therefore be tracked.
531  (
532  const polyMesh& mesh,
533  const label celli
534  );
535 
536 
537  // Decompose and reconstruct
538 
539  //- Return the tet point appropriate for decomposition or
540  // reconstruction to or from the given mesh.
542  (
543  const polyMesh& mesh,
544  const polyMesh& procMesh,
545  const label procCell,
546  const label procTetFace
547  ) const;
548 
549 
550  // I-O
551 
552  //- Read the fields associated with the owner cloud
553  template<class TrackCloudType>
554  static void readFields(TrackCloudType& c);
555 
556  //- Write the fields associated with the owner cloud
557  template<class TrackCloudType>
558  static void writeFields(const TrackCloudType& c);
559 
560  //- Write the particle position and cell
561  void writePosition(Ostream&) const;
562 
563 
564  // Friend Operators
565 
566  friend Ostream& operator<<(Ostream&, const particle&);
567 
568  friend bool operator==(const particle& pA, const particle& pB);
569 
570  friend bool operator!=(const particle& pA, const particle& pB);
571 };
572 
573 
574 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
575 
576 } // End namespace Foam
577 
578 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
579 
580 #include "particleI.H"
581 
582 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
583 
584 #ifdef NoRepository
585  #include "particleTemplates.C"
586 #endif
587 
588 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
589 
590 #endif
591 
592 // ************************************************************************* //
Intrusive doubly-linked list.
Template class for intrusive linked lists.
Definition: ILList.H:67
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
labelList patchNLocateBoundaryHits
Number of boundary hits that occurred during locate executions.
Definition: particle.H:125
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:109
trackingData(const TrackCloudType &cloud)
Definition: particle.H:130
label sendFromPatch
Patch from which to send the particle.
Definition: particle.H:112
vector sendToPosition
Position to which to send.
Definition: particle.H:121
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
label sendToPatch
Patch to which to send the particle.
Definition: particle.H:115
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:102
label sendToPatchFace
Patch face to which to send the particle.
Definition: particle.H:118
Base particle class.
Definition: particle.H:83
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:72
void prepareForParallelTransfer(TrackCloudType &, trackingData &)
Make changes prior to a parallel transfer. Runs either.
vector normal(const polyMesh &mesh) const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:129
label tetPt() const
Return current tet face particle is in.
Definition: particleI.H:66
label origProc() const
Return the originating processor ID.
Definition: particleI.H:96
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:48
void prepareForNonConformalCyclicTransfer(const polyMesh &mesh, const label sendToPatch, const label sendToPatchFace, const vector &sendToPosition)
Make changes prior to a transfer across a non conformal cyclic.
Definition: particle.C:261
virtual ~particle()
Destructor.
Definition: particle.H:255
scalar trackToCell(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
Definition: particle.C:158
label cell() const
Return current cell particle is in.
Definition: particleI.H:54
void hitBasicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a basic.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void prepareForProcessorTransfer(trackingData &td)
Make changes prior to a transfer across a processor boundary.
Definition: particle.C:227
void reset(const scalar stepFraction)
Set the step fraction and clear the behind data in preparation.
Definition: particleI.H:174
static const bool instantaneous
Particles are not instantaneous by default.
Definition: particle.H:209
label patch(const polyMesh &mesh) const
Return the index of patch that the particle is on.
Definition: particleI.H:153
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from components.
Definition: particle.C:47
void correctAfterInteractionListReferral(const polyMesh &mesh, const label celli)
Correct the topology after referral. Locates the particle.
Definition: particle.C:360
scalar trackToAndHitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
void correctAfterNonConformalCyclicTransfer(const polyMesh &mesh, const label sendToPatch, labelList &patchNLocateBoundaryHits)
Make changes following a transfer across a non conformal cyclic.
Definition: particle.C:301
friend bool operator!=(const particle &pA, const particle &pB)
scalar track(const polyMesh &mesh, const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:134
scalar trackToFace(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:182
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:84
void hitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
friend bool operator==(const particle &pA, const particle &pB)
void prepareForInteractionListReferral(const polyMesh &mesh, const transformer &transform)
Break the topology and store the cartesian position so that the.
Definition: particle.C:335
bool onFace() const
Is the particle on a face?
Definition: particleI.H:135
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
static autoPtr< particle > New(Istream &is)
Construct from Istream and return.
Definition: particle.H:248
friend Ostream & operator<<(Ostream &, const particle &)
virtual autoPtr< particle > clone() const
Construct and return a clone.
Definition: particle.H:242
void patchData(const polyMesh &mesh, vector &normal, vector &displacement) const
Get the normal and displacement of the current patch location.
Definition: particleI.H:183
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
bool locate(const polyMesh &mesh, const vector &position, label celli)
Locate the particle at the given position.
Definition: particle.C:113
vector deviationFromMeshCentre(const polyMesh &mesh) const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:206
void correctAfterProcessorTransfer(trackingData &td)
Make changes following a transfer across a processor boundary.
Definition: particle.C:239
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:86
TypeName("particle")
Runtime type information.
bool onBoundaryFace(const polyMesh &mesh) const
Is the particle on a boundary face?
Definition: particleI.H:147
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
tetIndices currentTetIndices(const polyMesh &mesh) const
Return the indices of the current tet that the.
Definition: particleI.H:121
label procTetPt(const polyMesh &mesh, const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or.
Definition: particle.C:389
bool onInternalFace(const polyMesh &mesh) const
Is the particle on an internal face?
Definition: particleI.H:141
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:159
void correctAfterParallelTransfer(TrackCloudType &, trackingData &)
Make changes following a parallel transfer. Runs either.
label getNewParticleIndex() const
Get unique particle creation id.
Definition: particleI.H:33
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
Definition: particle.C:223
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:108
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
static label particleCount
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:203
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
label tetFace() const
Return current tet face particle is in.
Definition: particleI.H:60
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:424
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
Macros for adding to particle property lists.
#define DefinePropertyList(str)
volScalarField & p