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