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-2020 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 cyclicAMIPolyPatch;
60 class cyclicACMIPolyPatch;
61 class cyclicRepeatAMIPolyPatch;
62 class processorPolyPatch;
63 class symmetryPlanePolyPatch;
64 class symmetryPolyPatch;
65 class wallPolyPatch;
66 class wedgePolyPatch;
67 
68 // Forward declaration of friend functions and operators
69 
70 Ostream& operator<<
71 (
72  Ostream&,
73  const particle&
74 );
75 
76 bool operator==(const particle&, const particle&);
77 
78 bool operator!=(const particle&, const particle&);
79 
80 /*---------------------------------------------------------------------------*\
81  Class Particle Declaration
82 \*---------------------------------------------------------------------------*/
83 
84 class particle
85 :
86  public IDLList<particle>::link
87 {
88  // Private member data
89 
90  //- Size in bytes of the position data
91  static const std::size_t sizeofPosition_;
92 
93  //- Size in bytes of the fields
94  static const std::size_t sizeofFields_;
95 
96  //- The value of nBehind_ at which tracking is abandoned. See the
97  // description of nBehind_.
98  static const label maxNBehind_;
99 
100 
101 public:
103  class trackingData
104  {
105  public:
106 
107  // Public data
108 
109  //- Flag to switch processor
110  bool switchProcessor;
111 
112  //- Flag to indicate whether to keep particle (false = delete)
113  bool keepParticle;
114 
115 
116  // Constructor
117  template <class TrackCloudType>
118  trackingData(const TrackCloudType& cloud)
119  {}
120  };
121 
122 
123 private:
124 
125  // Private Data
126 
127  //- Reference to the polyMesh database
128  const polyMesh& mesh_;
129 
130  //- Coordinates of particle
131  barycentric coordinates_;
132 
133  //- Index of the cell it is in
134  label celli_;
135 
136  //- Index of the face that owns the decomposed tet that the
137  // particle is in
138  label tetFacei_;
139 
140  //- Index of the point on the face that defines the decomposed
141  // tet that the particle is in. Relative to the face base
142  // point.
143  label tetPti_;
144 
145  //- Face index if the particle is on a face otherwise -1
146  label facei_;
147 
148  //- Fraction of time-step completed
149  scalar stepFraction_;
150 
151  //- The distance behind the maximum distance reached so far
152  scalar behind_;
153 
154  //- The number of tracks carried out that ended in a distance behind the
155  // maximum distance reached so far. Once this reaches maxNBehind_,
156  // tracking is abandoned for the current step. This is needed because
157  // when tetrahedra are inverted a straight trajectory can form a closed
158  // loop through regions of overlapping positive and negative space.
159  // Without this break clause, such loops can result in a valid track
160  // which never ends.
161  label nBehind_;
162 
163  //- Originating processor id
164  label origProc_;
165 
166  //- Local particle id on originating processor
167  label origId_;
168 
169 
170 private:
171 
172  // Private Member Functions
173 
174  // Tetrahedra functions
175 
176  //- Get the vertices of the current tet
177  inline void stationaryTetGeometry
178  (
179  vector& centre,
180  vector& base,
181  vector& vertex1,
182  vector& vertex2
183  ) const;
184 
185  //- Get the transformation associated with the current tet. This
186  // will convert a barycentric position within the tet to a
187  // cartesian position in the global coordinate system. The
188  // conversion is x = A & y, where x is the cartesian position, y is
189  // the barycentric position and A is the transformation tensor.
190  inline barycentricTensor stationaryTetTransform() const;
191 
192  //- Get the reverse transform associated with the current tet. The
193  // conversion is detA*y = (x - centre) & T. The variables x, y and
194  // centre have the same meaning as for the forward transform. T is
195  // the transposed inverse of the forward transform tensor, A,
196  // multiplied by its determinant, detA. This separation allows
197  // the barycentric tracking algorithm to function on inverted or
198  // degenerate tetrahedra.
199  void stationaryTetReverseTransform
200  (
201  vector& centre,
202  scalar& detA,
204  ) const;
205 
206  //- Get the vertices of the current moving tet. Two values are
207  // returned for each vertex. The first is a constant, and the
208  // second is a linear coefficient of the track fraction.
209  inline void movingTetGeometry
210  (
211  const scalar endStepFraction,
212  Pair<vector>& centre,
213  Pair<vector>& base,
214  Pair<vector>& vertex1,
215  Pair<vector>& vertex2
216  ) const;
217 
218  //- Get the transformation associated with the current, moving, tet.
219  // This is of the same form as for the static case. As with the
220  // moving geometry, a linear function of the tracking fraction is
221  // returned for each component.
222  inline Pair<barycentricTensor> movingTetTransform
223  (
224  const scalar endStepFraction
225  ) const;
226 
227  //- Get the reverse transformation associated with the current,
228  // moving, tet. This is of the same form as for the static case. As
229  // with the moving geometry, a function of the tracking fraction is
230  // returned for each component. The functions are higher order than
231  // for the forward transform; the determinant is cubic, and the
232  // tensor is quadratic.
233  void movingTetReverseTransform
234  (
235  const scalar endStepFraction,
236  Pair<vector>& centre,
237  FixedList<scalar, 4>& detA,
239  ) const;
240 
241 
242  // Transformations
243 
244  //- Reflection transform. Corrects the coordinates when the particle
245  // moves between two tets which share a base vertex, but for which
246  // the other two non cell-centre vertices are reversed. All hits
247  // which retain the same face behave this way, as do face hits.
248  void reflect();
249 
250  //- Rotation transform. Corrects the coordinates when the particle
251  // moves between two tets with different base vertices, but are
252  // otherwise similarly oriented. Hits which change the face within
253  // the cell make use of both this and the reflect transform.
254  void rotate(const bool direction);
255 
256 
257  // Topology changes
258 
259  //- Change tet within a cell. Called after a triangle is hit.
260  void changeTet(const label tetTriI);
261 
262  //- Change tet face within a cell. Called by changeTet.
263  void changeFace(const label tetTriI);
264 
265  //- Change cell. Called when the particle hits an internal face.
266  void changeCell();
267 
268  //- Put the particle on the lowest indexed patch for the current set
269  // of coincident faces. In the case of an ACMI-wall pair, this will
270  // move the particle from the wall face to the ACMI face, because
271  // ACMI patches are always listed before their associated non-
272  // overlapping patch.
273  void changeToMasterPatch();
274 
275 
276  // Geometry changes
277 
278  //- Locate the particle at the given position
279  void locate
280  (
281  const vector& position,
282  label celli,
283  const bool boundaryFail,
284  const string boundaryMsg
285  );
286 
287 
288 protected:
289 
290  // Patch interactions
291 
292  //- Overridable function to handle the particle hitting a patch.
293  // Executed before other patch-hitting functions.
294  template<class TrackCloudType>
295  bool hitPatch(TrackCloudType&, trackingData&);
296 
297  //- Overridable function to handle the particle hitting a wedgePatch
298  template<class TrackCloudType>
299  void hitWedgePatch(TrackCloudType&, trackingData&);
300 
301  //- Overridable function to handle the particle hitting a
302  // symmetryPlanePatch
303  template<class TrackCloudType>
304  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
305 
306  //- Overridable function to handle the particle hitting a symmetryPatch
307  template<class TrackCloudType>
308  void hitSymmetryPatch(TrackCloudType&, trackingData&);
309 
310  //- Overridable function to handle the particle hitting a cyclicPatch
311  template<class TrackCloudType>
312  void hitCyclicPatch(TrackCloudType&, trackingData&);
313 
314  //- Overridable function to handle the particle hitting a cyclicAMIPatch
315  template<class TrackCloudType>
316  void hitCyclicAMIPatch
317  (
318  const vector& displacement,
319  const scalar fraction,
320  TrackCloudType& cloud,
321  trackingData& td
322  );
323 
324  //- Overridable function to handle the particle hitting a
325  // cyclicACMIPatch
326  template<class TrackCloudType>
327  void hitCyclicACMIPatch
328  (
329  const vector& displacement,
330  const scalar fraction,
331  TrackCloudType& cloud,
332  trackingData& td
333  );
334 
335  //- Overridable function to handle the particle hitting an
336  // cyclicRepeatAMIPolyPatch
337  template<class TrackCloudType>
339  (
340  const vector& displacement,
341  const scalar fraction,
342  TrackCloudType& cloud,
343  trackingData& td
344  );
345 
346  //- Overridable function to handle the particle hitting a processorPatch
347  template<class TrackCloudType>
348  void hitProcessorPatch(TrackCloudType&, trackingData&);
349 
350  //- Overridable function to handle the particle hitting a wallPatch
351  template<class TrackCloudType>
352  void hitWallPatch(TrackCloudType&, trackingData&);
353 
354 
355 public:
356 
357  // Static Data Members
358 
359  //- Runtime type information
360  TypeName("particle");
361 
362  //- String representation of properties
364  (
365  "(coordinatesa coordinatesb coordinatesc coordinatesd) "
366  "celli tetFacei tetPti facei stepFraction "
367  "behind nBehind origProc origId"
368  );
369 
370  //- Cumulative particle counter - used to provide unique ID
371  static label particleCount_;
372 
373 
374  // Constructors
375 
376  //- Construct from components
377  particle
378  (
379  const polyMesh& mesh,
380  const barycentric& coordinates,
381  const label celli,
382  const label tetFacei,
383  const label tetPti
384  );
385 
386  //- Construct from a position and a cell, searching for the rest of the
387  // required topology
388  particle
389  (
390  const polyMesh& mesh,
391  const vector& position,
392  const label celli
393  );
394 
395  //- Construct from Istream
396  particle(const polyMesh& mesh, Istream&, bool readFields = true);
397 
398  //- Construct as a copy
399  particle(const particle& p);
400 
401  //- Construct as a copy with references to a new mesh
402  particle(const particle& p, const polyMesh& mesh);
403 
404  //- Construct a clone
405  virtual autoPtr<particle> clone() const
406  {
407  return autoPtr<particle>(new particle(*this));
408  }
409 
410  //- Factory class to read-construct particles used for
411  // parallel transfer
412  class iNew
413  {
414  const polyMesh& mesh_;
415 
416  public:
418  iNew(const polyMesh& mesh)
419  :
420  mesh_(mesh)
421  {}
423  autoPtr<particle> operator()(Istream& is) const
424  {
425  return autoPtr<particle>(new particle(mesh_, is, true));
426  }
427  };
428 
429 
430  //- Destructor
431  virtual ~particle()
432  {}
433 
434 
435  // Member Functions
436 
437  // Access
438 
439  //- Get unique particle creation id
440  inline label getNewParticleID() const;
441 
442  //- Return the mesh database
443  inline const polyMesh& mesh() const;
444 
445  //- Return current particle coordinates
446  inline const barycentric& coordinates() const;
447 
448  //- Return current cell particle is in
449  inline label cell() const;
450 
451  //- Return current tet face particle is in
452  inline label tetFace() const;
453 
454  //- Return current tet face particle is in
455  inline label tetPt() const;
456 
457  //- Return current face particle is on otherwise -1
458  inline label face() const;
459 
460  //- Return the fraction of time-step completed
461  inline scalar stepFraction() const;
462 
463  //- Return the fraction of time-step completed
464  inline scalar& stepFraction();
465 
466  //- Return the originating processor ID
467  inline label origProc() const;
468 
469  //- Return the originating processor ID
470  inline label& origProc();
471 
472  //- Return the particle ID on the originating processor
473  inline label origId() const;
474 
475  //- Return the particle ID on the originating processor
476  inline label& origId();
477 
478 
479  // Check
480 
481  //- Return the step fraction change within the overall time-step.
482  // Returns the start value and the change as a scalar pair. Always
483  // return Pair<scalar>(0, 1), unless sub-cycling is in effect, in
484  // which case the values will reflect the span of the sub-cycle
485  // within the time-step.
486  inline Pair<scalar> stepFractionSpan() const;
487 
488  //- Return the current fraction within the timestep. This differs
489  // from the stored step fraction due to sub-cycling.
490  inline scalar currentTimeFraction() const;
491 
492  //- Return the indices of the current tet that the
493  // particle occupies.
494  inline tetIndices currentTetIndices() const;
495 
496  //- Return the current tet transformation tensor
497  inline barycentricTensor currentTetTransform() const;
498 
499  //- Return the normal of the tri on tetFacei_ for the
500  // current tet.
501  inline vector normal() const;
502 
503  //- Is the particle on a face?
504  inline bool onFace() const;
505 
506  //- Is the particle on an internal face?
507  inline bool onInternalFace() const;
508 
509  //- Is the particle on a boundary face?
510  inline bool onBoundaryFace() const;
511 
512  //- Return the index of patch that the particle is on
513  inline label patch() const;
514 
515  //- Return current particle position
516  inline vector position() const;
517 
518 
519  // Track
520 
521  //- Set step fraction and behind data to zero in preparation for a new
522  // track
523  inline void reset();
524 
525  //- Track along the displacement for a given fraction of the overall
526  // step. End when the track is complete, or when a boundary is hit.
527  // On exit, stepFraction_ will have been incremented to the current
528  // position, and facei_ will be set to the index of the boundary
529  // face that was hit, or -1 if the track completed within a cell.
530  // The proportion of the displacement still to be completed is
531  // returned.
532  scalar track
533  (
534  const vector& displacement,
535  const scalar fraction
536  );
537 
538  //- As particle::track, but stops when a new cell is reached.
539  scalar trackToCell
540  (
541  const vector& displacement,
542  const scalar fraction
543  );
544 
545  //- As particle::track, but stops when a face is hit.
546  scalar trackToFace
547  (
548  const vector& displacement,
549  const scalar fraction
550  );
551 
552  //- As particle::trackToFace, but stops when a tet triangle is hit. On
553  // exit, tetTriI is set to the index of the tet triangle that was hit,
554  // or -1 if the end position was reached.
555  scalar trackToTri
556  (
557  const vector& displacement,
558  const scalar fraction,
559  label& tetTriI
560  );
561 
562  //- As particle::trackToTri, but for stationary meshes
563  scalar trackToStationaryTri
564  (
565  const vector& displacement,
566  const scalar fraction,
567  label& tetTriI
568  );
569 
570  //- As particle::trackToTri, but for moving meshes
571  scalar trackToMovingTri
572  (
573  const vector& displacement,
574  const scalar fraction,
575  label& tetTriI
576  );
577 
578  //- Hit the current face. If the current face is internal than this
579  // crosses into the next cell. If it is a boundary face then this will
580  // interact the particle with the relevant patch.
581  template<class TrackCloudType>
582  void hitFace
583  (
584  const vector& displacement,
585  const scalar fraction,
586  TrackCloudType& cloud,
587  trackingData& td
588  );
589 
590  //- As above, but does not change the master patch. Needed in order for
591  // ACMI to be able to delegate a hit to the non-overlap patch.
592  template<class TrackCloudType>
594  (
595  const vector& displacement,
596  const scalar fraction,
597  TrackCloudType& cloud,
598  trackingData& td
599  );
600 
601  //- Convenience function. Combines trackToFace and hitFace
602  template<class TrackCloudType>
603  scalar trackToAndHitFace
604  (
605  const vector& displacement,
606  const scalar fraction,
607  TrackCloudType& cloud,
608  trackingData& td
609  );
610 
611  //- Get the displacement from the mesh centre. Used to correct the
612  // particle position in cases with reduced dimensionality. Returns a
613  // zero vector for three-dimensional cases.
615 
616 
617  // Patch data
618 
619  //- Get the normal and displacement of the current patch location
620  inline void patchData(vector& normal, vector& displacement) const;
621 
622 
623  // Transformations
624 
625  //- Transform the physical properties of the particle
626  // according to the given transformation tensor
627  virtual void transformProperties(const transformer&);
628 
629 
630  // Parallel transfer
631 
632  //- Convert global addressing to the processor patch local equivalents
634 
635  //- Convert processor patch addressing to the global equivalents
636  // and set the celli to the face-neighbour
638 
639 
640  // Interaction list referral
641 
642  //- Break the topology and store the particle position so that the
643  // particle can be referred.
645  (
646  const transformer& transform
647  );
648 
649  //- Correct the topology after referral. The particle may still be
650  // outside the stored tet and therefore not track-able.
651  void correctAfterInteractionListReferral(const label celli);
652 
653 
654  // Decompose and reconstruct
655 
656  //- Return the tet point appropriate for decomposition or reconstruction
657  // to or from the given mesh.
659  (
660  const polyMesh& procMesh,
661  const label procCell,
662  const label procTetFace
663  ) const;
664 
665 
666  // Mapping
667 
668  //- Map after a topology change
669  void autoMap(const vector& position, const mapPolyMesh& mapper);
670 
671 
672  // I-O
673 
674  //- Read the fields associated with the owner cloud
675  template<class TrackCloudType>
676  static void readFields(TrackCloudType& c);
677 
678  //- Write the fields associated with the owner cloud
679  template<class TrackCloudType>
680  static void writeFields(const TrackCloudType& c);
681 
682  //- Write the particle position and cell
683  void writePosition(Ostream&) const;
684 
685 
686  // Friend Operators
687 
688  friend Ostream& operator<<(Ostream&, const particle&);
689 
690  friend bool operator==(const particle& pA, const particle& pB);
691 
692  friend bool operator!=(const particle& pA, const particle& pB);
693 };
694 
695 
696 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
697 
698 } // End namespace Foam
699 
700 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
701 
702 #include "particleI.H"
703 
704 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
705 
706 #ifdef NoRepository
707  #include "particleTemplates.C"
708 #endif
709 
710 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
711 
712 #endif
713 
714 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
void hitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
label origProc() const
Return the originating processor ID.
Definition: particleI.H:173
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:1060
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:125
trackingData(const TrackCloudType &cloud)
Definition: particle.H:117
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
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space...
Definition: transformer.H:83
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:1033
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:404
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:266
friend bool operator!=(const particle &pA, const particle &pB)
Macros for adding to particle property lists.
label procTetPt(const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or reconstruction.
Definition: particle.C:1156
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:645
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but stops when a tet triangle is hit. On.
Definition: particle.C:1016
scalar trackToAndHitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
Pair< scalar > stepFractionSpan() const
Return the step fraction change within the overall time-step.
Definition: particleI.H:197
scalar trackToCell(const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
Definition: particle.C:623
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
Definition: particle.C:1048
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
const dimensionedScalar c
Speed of light in a vacuum.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Base particle class.
Definition: particle.H:83
label tetPt() const
Return current tet face particle is in.
Definition: particleI.H:149
friend Ostream & operator<<(Ostream &, const particle &)
vector normal() const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:248
void patchData(vector &normal, vector &displacement) const
Get the normal and displacement of the current patch location.
Definition: particleI.H:292
label cell() const
Return current cell particle is in.
Definition: particleI.H:137
void hitFaceNoChangeToMasterPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
As above, but does not change the master patch. Needed in order for.
void prepareForInteractionListReferral(const transformer &transform)
Break the topology and store the particle position so that the.
Definition: particle.C:1094
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
TypeName("particle")
Runtime type information.
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:155
void reset()
Set step fraction and behind data to zero in preparation for a new.
Definition: particleI.H:284
Intrusive doubly-linked list.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
bool switchProcessor
Flag to switch processor.
Definition: particle.H:109
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:513
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
void hitCyclicACMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition: particle.C:1118
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition: particle.C:1052
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:235
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:185
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1182
bool onFace() const
Is the particle on a face?
Definition: particleI.H:254
label tetFace() const
Return current tet face particle is in.
Definition: particleI.H:143
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:111
Factory class to read-construct particles used for.
Definition: particle.H:411
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:161
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:229
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:830
label patchi
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:112
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:87
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
friend bool operator==(const particle &pA, const particle &pB)
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:131
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:367
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:260
void hitCyclicAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a cyclicAMIPatch.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:706
volScalarField & p
#define DefinePropertyList(str)
virtual ~particle()
Destructor.
Definition: particle.H:430
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:599
bool operator!=(const particle &, const particle &)
Definition: particle.C:1205
vector position() const
Return current particle position.
Definition: particleI.H:278
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:272
void hitCyclicRepeatAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:221