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