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-2023 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  //- The value of nTracksBehind_ at which tracking is abandoned. See the
94  // description of nTracksBehind_.
95  static const label maxNTracksBehind_;
96 
97 
98 public:
99 
100  class trackingData
101  {
102  public:
103 
104  // Public data
105 
106  //- Reference to the mesh
107  const polyMesh& mesh;
108 
109  //- Flag to indicate whether to keep particle (false = delete)
110  bool keepParticle;
111 
112  //- Processor to send the particle to. -1 indicates that this
113  // particle is not to be transferred.
115 
116  //- Patch from which to send the particle
118 
119  //- Patch to which to send the particle
121 
122  //- Patch face to which to send the particle
124 
125  //- Position to which to send
127 
128  //- Number of boundary hits that occurred during locate executions
129  // following (non-conformal) patch transfers. For reporting.
131 
132 
133  // Constructor
134  template <class TrackCloudType>
135  trackingData(const TrackCloudType& cloud)
136  :
137  mesh(cloud.pMesh()),
138  keepParticle(false),
139  sendToProc(-1),
140  sendFromPatch(-1),
141  sendToPatch(-1),
142  sendToPatchFace(-1),
143  sendToPosition(vector::uniform(NaN)),
145  (
146  mesh.boundaryMesh().size()
147  - mesh.globalData().processorPatches().size(),
148  0
149  )
150  {}
151  };
152 
153 
154 private:
155 
156  // Private Data
157 
158  //- Coordinates of particle
159  barycentric coordinates_;
160 
161  //- Index of the cell it is in
162  label celli_;
163 
164  //- Index of the face that owns the decomposed tet that the
165  // particle is in
166  label tetFacei_;
167 
168  //- Index of the point on the face that defines the decomposed
169  // tet that the particle is in. Relative to the face base
170  // point.
171  label tetPti_;
172 
173  //- Face index if the particle is on a face otherwise -1
174  label facei_;
175 
176  //- Fraction of time-step completed
177  scalar stepFraction_;
178 
179  //- The step fraction less than the maximum reached so far
180  scalar stepFractionBehind_;
181 
182  //- The number of tracks carried out that ended in a step fraction less
183  // than the maximum reached so far. Once this reaches
184  // maxNTracksBehind_, tracking is abandoned for the current step.
185  //
186  // This is needed because when tetrahedra are inverted a straight
187  // trajectory can form a closed loop through regions of overlapping
188  // positive and negative space. Without this break clause, such loops
189  // can result in a valid track which never ends.
190  //
191  // Because the test is susceptible to round off error, a track of zero
192  // length will also increment the counter. As such, it is important
193  // that maxNTracksBehind_ is set large enough so that valid small
194  // tracks do not result in the track being abandoned. The largest
195  // number of valid small tracks that are likely to be performed
196  // sequentially is equal to the number of tetrahedra that can meet at
197  // a point. An estimate of this number is therefore used to set
198  // maxNTracksBehind_.
199  label nTracksBehind_;
200 
201  //- Originating processor id
202  label origProc_;
203 
204  //- Local particle id on originating processor
205  label origId_;
206 
207 
208 private:
209 
210  // Private Member Functions
211 
212  // Tetrahedra functions
213 
214  //- Get the vertices of the current tet
215  inline void stationaryTetGeometry
216  (
217  const polyMesh& mesh,
218  vector& centre,
219  vector& base,
220  vector& vertex1,
221  vector& vertex2
222  ) const;
223 
224  //- Get the transformation associated with the current tet. This
225  // will convert a barycentric position within the tet to a
226  // cartesian position in the global coordinate system. The
227  // conversion is x = A & y, where x is the cartesian position, y is
228  // the barycentric position and A is the transformation tensor.
229  inline barycentricTensor stationaryTetTransform
230  (
231  const polyMesh& mesh
232  ) const;
233 
234  //- Get the reverse transform associated with the current tet. The
235  // conversion is detA*y = (x - centre) & T. The variables x, y and
236  // centre have the same meaning as for the forward transform. T is
237  // the transposed inverse of the forward transform tensor, A,
238  // multiplied by its determinant, detA. This separation allows
239  // the barycentric tracking algorithm to function on inverted or
240  // degenerate tetrahedra.
241  void stationaryTetReverseTransform
242  (
243  const polyMesh& mesh,
244  vector& centre,
245  scalar& detA,
247  ) const;
248 
249  //- Get the vertices of the current moving tet. Two values are
250  // returned for each vertex. The first is a constant, and the
251  // second is a linear coefficient of the track fraction.
252  inline void movingTetGeometry
253  (
254  const polyMesh& mesh,
255  const scalar endStepFraction,
256  Pair<vector>& centre,
257  Pair<vector>& base,
258  Pair<vector>& vertex1,
259  Pair<vector>& vertex2
260  ) const;
261 
262  //- Get the transformation associated with the current, moving, tet.
263  // This is of the same form as for the static case. As with the
264  // moving geometry, a linear function of the tracking fraction is
265  // returned for each component.
266  inline Pair<barycentricTensor> movingTetTransform
267  (
268  const polyMesh& mesh,
269  const scalar endStepFraction
270  ) const;
271 
272  //- Get the reverse transformation associated with the current,
273  // moving, tet. This is of the same form as for the static case. As
274  // with the moving geometry, a function of the tracking fraction is
275  // returned for each component. The functions are higher order than
276  // for the forward transform; the determinant is cubic, and the
277  // tensor is quadratic.
278  void movingTetReverseTransform
279  (
280  const polyMesh& mesh,
281  const scalar endStepFraction,
282  Pair<vector>& centre,
283  FixedList<scalar, 4>& detA,
285  ) const;
286 
287 
288  // Transformations
289 
290  //- Reflection transform. Corrects the coordinates when the particle
291  // moves between two tets which share a base vertex, but for which
292  // the other two non cell-centre vertices are reversed. All hits
293  // which retain the same face behave this way, as do face hits.
294  void reflect();
295 
296  //- Rotation transform. Corrects the coordinates when the particle
297  // moves between two tets with different base vertices, but are
298  // otherwise similarly oriented. Hits which change the face within
299  // the cell make use of both this and the reflect transform.
300  void rotate(const bool direction);
301 
302 
303  // Topology changes
304 
305  //- Change tet within a cell. Called after a triangle is hit.
306  void changeTet(const polyMesh& mesh, const label tetTriI);
307 
308  //- Change tet face within a cell. Called by changeTet.
309  void changeFace(const polyMesh& mesh, const label tetTriI);
310 
311  //- Change cell. Called when the particle hits an internal face.
312  void changeCell(const polyMesh& mesh);
313 
314 
315 public:
316 
317  // Static Data Members
318 
319  //- Runtime type information
320  TypeName("particle");
321 
322  //- String representation of properties
324  (
325  "(coordinatesa coordinatesb coordinatesc coordinatesd) "
326  "celli tetFacei tetPti facei stepFraction "
327  "behind nBehind origProc origId"
328  );
329 
330  //- Cumulative particle counter - used to provide unique ID
331  static label particleCount_;
332 
333 
334  // Constructors
335 
336  //- Construct from components
337  particle
338  (
339  const polyMesh& mesh,
340  const barycentric& coordinates,
341  const label celli,
342  const label tetFacei,
343  const label tetPti,
344  const label facei
345  );
346 
347  //- Construct from a position and a cell, searching for the rest of the
348  // required topology
349  particle
350  (
351  const polyMesh& mesh,
352  const vector& position,
353  const label celli,
354  label& nLocateBoundaryHits
355  );
356 
357  //- Construct from Istream
358  particle(Istream&, bool readFields = true);
359 
360  //- Construct as a copy
361  particle(const particle& p);
362 
363  //- Construct and return a clone
364  virtual autoPtr<particle> clone() const
365  {
366  return autoPtr<particle>(new particle(*this));
367  }
368 
369  //- Construct from Istream and return
370  static autoPtr<particle> New(Istream& is)
371  {
372  return autoPtr<particle>(new particle(is));
373  }
374 
375 
376  //- Destructor
377  virtual ~particle()
378  {}
379 
380 
381  // Member Functions
382 
383  // Access
384 
385  //- Get unique particle creation id
386  inline label getNewParticleIndex() const;
387 
388  //- Return current particle coordinates
389  inline const barycentric& coordinates() const;
390 
391  //- Return current cell particle is in
392  inline label cell() const;
393 
394  //- Return current tet face particle is in
395  inline label tetFace() const;
396 
397  //- Return current tet face particle is in
398  inline label tetPt() const;
399 
400  //- Return current face particle is on otherwise -1
401  inline label face() const;
402 
403  //- Return current face particle is on otherwise -1
404  inline label& face();
405 
406  //- Return the fraction of time-step completed
407  inline scalar stepFraction() const;
408 
409  //- Return the fraction of time-step completed
410  inline scalar& stepFraction();
411 
412  //- Return the originating processor ID
413  inline label origProc() const;
414 
415  //- Return the originating processor ID
416  inline label& origProc();
417 
418  //- Return the particle ID on the originating processor
419  inline label origId() const;
420 
421  //- Return the particle ID on the originating processor
422  inline label& origId();
423 
424 
425  // Check
426 
427  //- Return the indices of the current tet that the
428  // particle occupies.
429  inline tetIndices currentTetIndices(const polyMesh& mesh) const;
430 
431  //- Return the current tet transformation tensor
433  (
434  const polyMesh& mesh
435  ) const;
436 
437  //- Return the normal of the tri on tetFacei_ for the
438  // current tet.
439  inline vector normal(const polyMesh& mesh) const;
440 
441  //- Is the particle on a face?
442  inline bool onFace() const;
443 
444  //- Is the particle on an internal face?
445  inline bool onInternalFace(const polyMesh& mesh) const;
446 
447  //- Is the particle on a boundary face?
448  inline bool onBoundaryFace(const polyMesh& mesh) const;
449 
450  //- Return the index of patch that the particle is on
451  inline label patch(const polyMesh& mesh) const;
452 
453  //- Return current particle position
454  inline vector position(const polyMesh& mesh) const;
455 
456 
457  // Track
458 
459  //- Set the step fraction and clear the behind data in preparation
460  // for a new track
461  inline void reset(const scalar stepFraction);
462 
463  //- Locate the particle at the given position
464  bool locate
465  (
466  const polyMesh& mesh,
467  const vector& position,
468  label celli
469  );
470 
471  //- Track along the displacement for a given fraction of the overall
472  // step. End when the track is complete, or when a boundary is hit.
473  // On exit, stepFraction_ will have been incremented to the current
474  // position, and facei_ will be set to the index of the boundary
475  // face that was hit, or -1 if the track completed within a cell.
476  // The proportion of the displacement still to be completed is
477  // returned.
478  scalar track
479  (
480  const polyMesh& mesh,
481  const vector& displacement,
482  const scalar fraction
483  );
484 
485  //- As particle::track, but stops when a new cell is reached.
486  scalar trackToCell
487  (
488  const polyMesh& mesh,
489  const vector& displacement,
490  const scalar fraction
491  );
492 
493  //- As particle::track, but stops when a face is hit.
494  scalar trackToFace
495  (
496  const polyMesh& mesh,
497  const vector& displacement,
498  const scalar fraction
499  );
500 
501  //- As particle::trackToFace, but stops when a tet triangle is hit.
502  // On exit, tetTriI is set to the index of the tet triangle that
503  // was hit, or -1 if the end position was reached.
504  scalar trackToTri
505  (
506  const polyMesh& mesh,
507  const vector& displacement,
508  const scalar fraction,
509  label& tetTriI
510  );
511 
512  //- As particle::trackToTri, but for stationary meshes
513  scalar trackToStationaryTri
514  (
515  const polyMesh& mesh,
516  const vector& displacement,
517  const scalar fraction,
518  label& tetTriI
519  );
520 
521  //- As particle::trackToTri, but for moving meshes
522  scalar trackToMovingTri
523  (
524  const polyMesh& mesh,
525  const vector& displacement,
526  const scalar fraction,
527  label& tetTriI
528  );
529 
530  //- Hit the current face. If the current face is internal than this
531  // crosses into the next cell. If it is a boundary face then this
532  // will interact the particle with the relevant patch.
533  template<class TrackCloudType>
534  void hitFace
535  (
536  const vector& displacement,
537  const scalar fraction,
538  TrackCloudType& cloud,
539  trackingData& td
540  );
541 
542  //- Convenience function. Combines trackToFace and hitFace
543  template<class TrackCloudType>
544  scalar trackToAndHitFace
545  (
546  const vector& displacement,
547  const scalar fraction,
548  TrackCloudType& cloud,
549  trackingData& td
550  );
551 
552  //- Get the displacement from the mesh centre. Used to correct the
553  // particle position in cases with reduced dimensionality. Returns
554  // a zero vector for three-dimensional cases.
555  vector deviationFromMeshCentre(const polyMesh& mesh) const;
556 
557 
558  // Patch data
559 
560  //- Get the normal and displacement of the current patch location
561  inline void patchData
562  (
563  const polyMesh& mesh,
564  vector& normal,
565  vector& displacement
566  ) const;
567 
568 
569  // Transformations
570 
571  //- Transform the physical properties of the particle
572  // according to the given transformation tensor
573  virtual void transformProperties(const transformer&);
574 
575 
576  // Transfers
577 
578  //- Make changes prior to a parallel transfer. Runs either
579  // processor or nonConformalCyclic variant below.
580  template<class TrackCloudType>
581  void prepareForParallelTransfer(TrackCloudType&, trackingData&);
582 
583  //- Make changes following a parallel transfer. Runs either
584  // processor or nonConformalCyclic variant below.
585  template<class TrackCloudType>
586  void correctAfterParallelTransfer(TrackCloudType&, trackingData&);
587 
588  //- Make changes prior to a transfer across a processor boundary.
589  // Stores the local patch face index (in facei_) so that the mesh
590  // face index can be determined on the other side.
592 
593  //- Make changes following a transfer across a processor boundary.
594  // Converts the stored patch index to a mesh index. Accounts for
595  // the receiving face being reversed relative to the sending face.
597 
598  //- Make changes prior to a transfer across a non conformal cyclic
599  // boundary. Stores the receiving patch face (in facei_). Breaks
600  // the topology and stores the cartesian position.
602  (
603  const polyMesh& mesh,
604  const label sendToPatch,
605  const label sendToPatchFace,
606  const vector& sendToPosition
607  );
608 
609  //- Make changes following a transfer across a non conformal cyclic
610  // boundary. Locates the particle using the stored face index and
611  // cartesian position.
613  (
614  const polyMesh& mesh,
615  const label sendToPatch,
616  labelList& patchNLocateBoundaryHits
617  );
618 
619 
620  // Patch interactions
621 
622  //- Overridable function to handle the particle hitting a patch.
623  // Executed before other patch-hitting functions.
624  template<class TrackCloudType>
625  bool hitPatch(TrackCloudType&, trackingData&);
626 
627  //- Overridable function to handle the particle hitting a wedgePatch
628  template<class TrackCloudType>
629  void hitWedgePatch(TrackCloudType&, trackingData&);
630 
631  //- Overridable function to handle the particle hitting a
632  // symmetryPlanePatch
633  template<class TrackCloudType>
634  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
635 
636  //- Overridable function to handle the particle hitting a
637  // symmetryPatch
638  template<class TrackCloudType>
639  void hitSymmetryPatch(TrackCloudType&, trackingData&);
640 
641  //- Overridable function to handle the particle hitting a
642  // cyclicPatch
643  template<class TrackCloudType>
644  void hitCyclicPatch(TrackCloudType&, trackingData&);
645 
646  //- Overridable function to handle the particle hitting an
647  // nonConformalCyclicPolyPatch
648  template<class TrackCloudType>
650  (
651  const vector& displacement,
652  const scalar fraction,
653  const label patchi,
654  TrackCloudType& cloud,
655  trackingData& td
656  );
657 
658  //- Overridable function to handle the particle hitting a
659  // processorPatch
660  template<class TrackCloudType>
661  void hitProcessorPatch(TrackCloudType&, trackingData&);
662 
663  //- Overridable function to handle the particle hitting a wallPatch
664  template<class TrackCloudType>
665  void hitWallPatch(TrackCloudType&, trackingData&);
666 
667  //- Overridable function to handle the particle hitting a basic
668  // patch. Fall-through for the above.
669  template<class TrackCloudType>
670  void hitBasicPatch(TrackCloudType&, trackingData&);
671 
672 
673  // Interaction list referral
674 
675  //- Break the topology and store the cartesian position so that the
676  // particle can be referred.
678  (
679  const polyMesh& mesh,
680  const transformer& transform
681  );
682 
683  //- Correct the topology after referral. Locates the particle
684  // relative to a nearby cell/tet. The particle may end up outside
685  // this cell/tet and cannot therefore be tracked.
687  (
688  const polyMesh& mesh,
689  const label celli
690  );
691 
692 
693  // Decompose and reconstruct
694 
695  //- Return the tet point appropriate for decomposition or
696  // reconstruction to or from the given mesh.
698  (
699  const polyMesh& mesh,
700  const polyMesh& procMesh,
701  const label procCell,
702  const label procTetFace
703  ) const;
704 
705 
706  // I-O
707 
708  //- Read the fields associated with the owner cloud
709  template<class TrackCloudType>
710  static void readFields(TrackCloudType& c);
711 
712  //- Write the fields associated with the owner cloud
713  template<class TrackCloudType>
714  static void writeFields(const TrackCloudType& c);
715 
716  //- Write the particle position and cell
717  void writePosition(Ostream&) const;
718 
719 
720  // Friend Operators
721 
722  friend Ostream& operator<<(Ostream&, const particle&);
723 
724  friend bool operator==(const particle& pA, const particle& pB);
725 
726  friend bool operator!=(const particle& pA, const particle& pB);
727 };
728 
729 
730 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
731 
732 } // End namespace Foam
733 
734 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
735 
736 #include "particleI.H"
737 
738 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
739 
740 #ifdef NoRepository
741  #include "particleTemplates.C"
742 #endif
743 
744 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
745 
746 #endif
747 
748 // ************************************************************************* //
Intrusive doubly-linked list.
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
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 ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:65
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
labelList patchNLocateBoundaryHits
Number of boundary hits that occurred during locate executions.
Definition: particle.H:129
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:113
trackingData(const TrackCloudType &cloud)
Definition: particle.H:134
label sendFromPatch
Patch from which to send the particle.
Definition: particle.H:116
vector sendToPosition
Position to which to send.
Definition: particle.H:125
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:109
label sendToPatch
Patch to which to send the particle.
Definition: particle.H:119
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:106
label sendToPatchFace
Patch face to which to send the particle.
Definition: particle.H:122
Base particle class.
Definition: particle.H:83
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:152
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:225
label tetPt() const
Return current tet face particle is in.
Definition: particleI.H:146
label origProc() const
Return the originating processor ID.
Definition: particleI.H:176
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:128
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:1036
virtual ~particle()
Destructor.
Definition: particle.H:376
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:542
label cell() const
Return current cell particle is in.
Definition: particleI.H:134
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:994
void reset(const scalar stepFraction)
Set the step fraction and clear the behind data in preparation.
Definition: particleI.H:261
label patch(const polyMesh &mesh) const
Return the index of patch that the particle is on.
Definition: particleI.H:249
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:375
void correctAfterInteractionListReferral(const polyMesh &mesh, const label celli)
Correct the topology after referral. Locates the particle.
Definition: particle.C:1135
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:1076
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:517
scalar trackToFace(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:565
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:164
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:1110
bool onFace() const
Is the particle on a face?
Definition: particleI.H:231
scalar trackToTri(const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but stops when a tet triangle is hit.
Definition: particle.C:954
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:369
friend Ostream & operator<<(Ostream &, const particle &)
virtual autoPtr< particle > clone() const
Construct and return a clone.
Definition: particle.H:363
void patchData(const polyMesh &mesh, vector &normal, vector &displacement) const
Get the normal and displacement of the current patch location.
Definition: particleI.H:270
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:441
scalar trackToMovingTri(const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:756
vector deviationFromMeshCentre(const polyMesh &mesh) const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:973
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:327
void correctAfterProcessorTransfer(trackingData &td)
Make changes following a transfer across a processor boundary.
Definition: particle.C:1001
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:243
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:201
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:1175
bool onInternalFace(const polyMesh &mesh) const
Is the particle on an internal face?
Definition: particleI.H:237
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:255
void correctAfterParallelTransfer(TrackCloudType &, trackingData &)
Make changes following a parallel transfer. Runs either.
scalar trackToStationaryTri(const polyMesh &mesh, const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:620
label getNewParticleIndex() const
Get unique particle creation id.
Definition: particleI.H:114
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:990
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:188
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
barycentricTensor currentTetTransform(const polyMesh &mesh) const
Return the current tet transformation tensor.
Definition: particleI.H:210
label tetFace() const
Return current tet face particle is in.
Definition: particleI.H:140
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
label patchi
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1210
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:504
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
uint8_t direction
Definition: direction.H:45
Macros for adding to particle property lists.
#define DefinePropertyList(str)
volScalarField & p