particle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017, 2020 OpenFOAM Foundation
9  Copyright (C) 2017-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::particle
29 
30 Description
31  Base particle class
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef Foam_particle_H
36 #define Foam_particle_H
37 
38 #include "vector.H"
39 #include "barycentric.H"
40 #include "barycentricTensor.H"
41 #include "Cloud.H"
42 #include "IDLList.H"
43 #include "pointField.H"
44 #include "faceList.H"
45 #include "OFstream.H"
46 #include "FixedList.H"
48 #include "particleMacros.H"
49 #include "vectorTensorTransform.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward Declarations
57 class particle;
58 class polyPatch;
59 class cyclicPolyPatch;
60 class cyclicAMIPolyPatch;
61 class cyclicACMIPolyPatch;
62 class processorPolyPatch;
63 class symmetryPlanePolyPatch;
64 class symmetryPolyPatch;
65 class wallPolyPatch;
66 class wedgePolyPatch;
67 
68 Ostream& operator<<(Ostream&, const particle&);
69 bool operator==(const particle&, const particle&);
70 bool operator!=(const particle&, const particle&);
71 
72 /*---------------------------------------------------------------------------*\
73  Class Particle Declaration
74 \*---------------------------------------------------------------------------*/
75 
76 class particle
77 :
78  public IDLList<particle>::link
79 {
80  // Private Data
81 
82  //- Size in bytes of the position data
83  static const std::size_t sizeofPosition;
84 
85  //- Size in bytes of the fields
86  static const std::size_t sizeofFields;
87 
88  //- The value of nBehind_ at which tracking is abandoned. See the
89  // description of nBehind_.
90  static const label maxNBehind_;
91 
92 
93 public:
94 
96  {
97  public:
98 
99  // Public Data
100 
101  //- Flag to switch processor
102  bool switchProcessor;
103 
104  //- Flag to indicate whether to keep particle (false = delete)
105  bool keepParticle;
106 
107 
108  // Constructor
109  template<class TrackCloudType>
110  trackingData(const TrackCloudType& cloud)
111  {}
112  };
113 
115  //- Old particle positions content for OpenFOAM-1706 and earlier
116  struct positionsCompat1706
117  {
119  label celli;
120  label facei;
121  scalar stepFraction;
122  label tetFacei;
123  label tetPti;
124  label origProc;
125  label origId;
126  };
129 private:
131  // Private Data
132 
133  //- Reference to the polyMesh database
134  const polyMesh& mesh_;
135 
136  //- Coordinates of particle
137  barycentric coordinates_;
138 
139  //- Index of the cell it is in
140  label celli_;
141 
142  //- Index of the face that owns the decomposed tet that the
143  //- particle is in
144  label tetFacei_;
145 
146  //- Index of the point on the face that defines the decomposed
147  //- tet that the particle is in. Relative to the face base
148  //- point.
149  label tetPti_;
150 
151  //- Face index if the particle is on a face otherwise -1
152  label facei_;
153 
154  //- Fraction of time-step completed
155  scalar stepFraction_;
156 
157  //- The distance behind the maximum distance reached so far
158  scalar behind_;
159 
160  //- The number of tracks carried out that ended in a distance behind the
161  //- maximum distance reached so far. Once this reaches maxNBehind_,
162  // tracking is abandoned for the current step. This is needed because
163  // when tetrahedra are inverted a straight trajectory can form a closed
164  // loop through regions of overlapping positive and negative space.
165  // Without this break clause, such loops can result in a valid track
166  // which never ends.
167  label nBehind_;
168 
169  //- Originating processor id
170  label origProc_;
171 
172  //- Local particle id on originating processor
173  label origId_;
174 
175 
176  // Private Member Functions
177 
178  // Tetrahedra functions
179 
180  //- Get the transformation associated with the current tet. This
181  // will convert a barycentric position within the tet to a
182  // Cartesian position in the global coordinate system. The
183  // conversion is x = A & y, where x is the Cartesian position, y is
184  // the barycentric position and A is the transformation tensor.
185  inline barycentricTensor stationaryTetTransform() const;
186 
187  //- Get the reverse transform associated with the current tet. The
188  // conversion is detA*y = (x - centre) & T. The variables x, y and
189  // centre have the same meaning as for the forward transform. T is
190  // the transposed inverse of the forward transform tensor, A,
191  // multiplied by its determinant, detA. This separation allows
192  // the barycentric tracking algorithm to function on inverted or
193  // degenerate tetrahedra.
194  void stationaryTetReverseTransform
195  (
196  vector& centre,
197  scalar& detA,
199  ) const;
200 
201  //- Get the vertices of the current moving tet. Two values are
202  // returned for each vertex. The first is a constant, and the
203  // second is a linear coefficient of the track fraction.
204  inline void movingTetGeometry
205  (
206  const scalar endStepFraction,
207  Pair<vector>& centre,
208  Pair<vector>& base,
209  Pair<vector>& vertex1,
210  Pair<vector>& vertex2
211  ) const;
212 
213  //- Get the transformation associated with the current, moving, tet.
214  // This is of the same form as for the static case. As with the
215  // moving geometry, a linear function of the tracking fraction is
216  // returned for each component.
217  inline Pair<barycentricTensor> movingTetTransform
218  (
219  const scalar endStepFraction
220  ) const;
221 
222  //- Get the reverse transformation associated with the current,
223  // moving, tet. This is of the same form as for the static case. As
224  // with the moving geometry, a function of the tracking fraction is
225  // returned for each component. The functions are higher order than
226  // for the forward transform; the determinant is cubic, and the
227  // tensor is quadratic.
228  void movingTetReverseTransform
229  (
230  const scalar endStepFraction,
231  Pair<vector>& centre,
232  FixedList<scalar, 4>& detA,
234  ) const;
235 
236 
237  // Transformations
238 
239  //- Reflection transform. Corrects the coordinates when the particle
240  // moves between two tets which share a base vertex, but for which
241  // the other two non cell-centre vertices are reversed. All hits
242  // which retain the same face behave this way, as do face hits.
243  void reflect();
244 
245  //- Rotation transform. Corrects the coordinates when the particle
246  // moves between two tets with different base vertices, but are
247  // otherwise similarly oriented. Hits which change the face within
248  // the cell make use of both this and the reflect transform.
249  void rotate(const bool direction);
250 
251 
252  // Topology changes
253 
254  //- Change tet within a cell. Called after a triangle is hit.
255  void changeTet(const label tetTriI);
256 
257  //- Change tet face within a cell. Called by changeTet.
258  void changeFace(const label tetTriI);
259 
260  //- Change cell. Called when the particle hits an internal face.
261  void changeCell();
262 
263  //- Put the particle on the lowest indexed patch for the current set
264  // of coincident faces. In the case of an ACMI-wall pair, this will
265  // move the particle from the wall face to the ACMI face, because
266  // ACMI patches are always listed before their associated non-
267  // overlapping patch.
268  void changeToMasterPatch();
269 
270 
271  // Geometry changes
272 
273  //- Locate the particle at the given position
274  void locate
275  (
276  const vector& position,
277  const vector* direction,
278  const label celli,
279  const bool boundaryFail,
280  const string& boundaryMsg
281  );
282 
283 
284 protected:
285 
286  // Patch Interactions
287 
288  //- Read particle from stream. Optionally (for old format) return
289  // read position. Used by construct-from-Istream
290  void readData
291  (
292  Istream& is,
293  point& position,
294  const bool readFields,
295  const bool newFormat,
296  const bool doLocate
297  );
298 
299  //- Overridable function to handle the particle hitting a patch.
300  // Executed before other patch-hitting functions.
301  template<class TrackCloudType>
302  bool hitPatch(TrackCloudType&, trackingData&);
303 
304  //- Overridable function to handle the particle hitting a wedgePatch
305  template<class TrackCloudType>
306  void hitWedgePatch(TrackCloudType&, trackingData&);
307 
308  //- Overridable function to handle the particle hitting a
309  // symmetryPlanePatch
310  template<class TrackCloudType>
311  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
312 
313  //- Overridable function to handle the particle hitting a symmetryPatch
314  template<class TrackCloudType>
315  void hitSymmetryPatch(TrackCloudType&, trackingData&);
316 
317  //- Overridable function to handle the particle hitting a cyclicPatch
318  template<class TrackCloudType>
319  void hitCyclicPatch(TrackCloudType&, trackingData&);
320 
321  //- Overridable function to handle the particle hitting a cyclicAMIPatch
322  template<class TrackCloudType>
323  void hitCyclicAMIPatch(TrackCloudType&, trackingData&, const vector&);
324 
325  //- Overridable function to handle the particle hitting a
326  // cyclicACMIPatch
327  template<class TrackCloudType>
328  void hitCyclicACMIPatch(TrackCloudType&, trackingData&, const vector&);
329 
330  //- Overridable function to handle the particle hitting a processorPatch
331  template<class TrackCloudType>
332  void hitProcessorPatch(TrackCloudType&, trackingData&);
333 
334  //- Overridable function to handle the particle hitting a wallPatch
335  template<class TrackCloudType>
336  void hitWallPatch(TrackCloudType&, trackingData&);
337 
338 
339  //- Dispatch function for boundary face interaction. Calls one of
340  // the above (hitWedgePatch, hitCyclicPatch etc) depending on the
341  // patch type
342  template<class TrackCloudType>
343  void hitBoundaryFace
344  (
345  const vector& direction,
346  TrackCloudType& cloud,
348  );
349 
350 
351 public:
352 
353  // Static Data Members
354 
355  //- Runtime type information
356  TypeName("particle");
357 
358  //- String representation of properties
360  (
361  "(coordinatesa coordinatesb coordinatesc coordinatesd) "
362  "celli tetFacei tetPti facei stepFraction origProc origId"
363  );
364 
365  //- Cumulative particle counter - used to provide unique ID
366  static label particleCount_;
367 
368  //- Write particle coordinates file (v1712 and later)
369  //- Default is true
370  static bool writeLagrangianCoordinates;
371 
372  //- Write particle positions file (v1706 format and earlier)
373  //- Default is true (disable in etc/controlDict)
374  static bool writeLagrangianPositions;
375 
376 
377  // Constructors
378 
379  //- Construct from components
380  particle
381  (
382  const polyMesh& mesh,
383  const barycentric& coordinates,
384  const label celli,
385  const label tetFacei,
386  const label tetPti
387  );
388 
389  //- Construct from a position and a cell.
390  // Searches for the rest of the required topology.
391  particle
392  (
393  const polyMesh& mesh,
394  const vector& position,
395  const label celli = -1
396  );
397 
398  //- Construct from position components
399  particle
400  (
401  const polyMesh& mesh,
402  const vector& position,
403  const label celli,
404  const label tetFacei,
405  const label tetPti,
406  const bool doLocate = true
407  );
408 
409 
410  //- Construct from Istream
411  particle
412  (
413  const polyMesh& mesh,
414  Istream&,
415  const bool readFields = true,
416  const bool newFormat = true,
417  const bool doLocate = true
418  );
419 
420  //- Construct as a copy with reference to a mesh
421  particle(const particle& p, const polyMesh& mesh);
422 
423  //- Copy construct
424  particle(const particle& p);
425 
426  //- Construct a clone
427  virtual autoPtr<particle> clone() const
428  {
429  return autoPtr<particle>::New(*this);
430  }
431 
432 
433  // Factory Methods
434 
435  //- Clone a particle
436  template<class Derived>
437  static autoPtr<particle> Clone(const Derived& p)
438  {
439  return autoPtr<particle>(new Derived(p));
440  }
441 
442  //- Clone a particle with a mesh reference
443  template<class Derived>
444  static autoPtr<particle> Clone(const Derived& p, const polyMesh& mesh)
445  {
446  return autoPtr<particle>(new Derived(p, mesh));
447  }
448 
449  //- Factory class to read-construct particles (for parallel transfer)
450  class iNew
451  {
452  const polyMesh& mesh_;
453 
454  public:
456  iNew(const polyMesh& mesh)
457  :
458  mesh_(mesh)
459  {}
460 
462  {
463  return autoPtr<particle>::New(mesh_, is, true);
464  }
465  };
467 
468  //- Destructor
469  virtual ~particle() = default;
470 
471 
472  // Member Functions
473 
474  // Access
475 
476  //- Get unique particle creation id
477  inline label getNewParticleID() const;
478 
479  //- Return the mesh database
480  inline const polyMesh& mesh() const noexcept;
481 
482  //- Return current particle coordinates
483  inline const barycentric& coordinates() const noexcept;
484 
485  //- Return current cell particle is in
486  inline label cell() const noexcept;
487 
488  //- Return current cell particle is in for manipulation
489  inline label& cell() noexcept;
490 
491  //- Return current tet face particle is in
492  inline label tetFace() const noexcept;
493 
494  //- Return current tet face particle is in for manipulation
495  inline label& tetFace() noexcept;
496 
497  //- Return current tet face particle is in
498  inline label tetPt() const noexcept;
499 
500  //- Return current tet face particle is in for manipulation
501  inline label& tetPt() noexcept;
502 
503  //- Return current face particle is on otherwise -1
504  inline label face() const noexcept;
505 
506  //- Return current face particle is on for manipulation
507  inline label& face() noexcept;
508 
509  //- Return the fraction of time-step completed
510  inline scalar stepFraction() const noexcept;
511 
512  //- Return the fraction of time-step completed
513  inline scalar& stepFraction() noexcept;
514 
515  //- Return the originating processor ID
516  inline label origProc() const noexcept;
517 
518  //- Return the originating processor ID
519  inline label& origProc() noexcept;
520 
521  //- Return the particle ID on the originating processor
522  inline label origId() const noexcept;
523 
524  //- Return the particle ID on the originating processor
525  inline label& origId() noexcept;
526 
527 
528  // Check
529 
530  //- Return the step fraction change within the overall time-step.
531  // Returns the start value and the change as a scalar pair. Always
532  // return Pair<scalar>(0, 1), unless sub-cycling is in effect, in
533  // which case the values will reflect the span of the sub-cycle
534  // within the time-step.
535  inline Pair<scalar> stepFractionSpan() const;
536 
537  //- Return the current fraction within the timestep. This differs
538  // from the stored step fraction due to sub-cycling.
539  inline scalar currentTimeFraction() const;
541  //- Return indices of the current tet that the particle occupies.
542  inline tetIndices currentTetIndices() const noexcept;
543 
544  //- Return the current tet transformation tensor
545  inline barycentricTensor currentTetTransform() const;
546 
547  //- The (unit) normal of the tri on tetFacei_ for the current tet.
548  inline vector normal() const;
549 
550  //- Is the particle on a face?
551  inline bool onFace() const noexcept;
553  //- Is the particle on an internal face?
554  inline bool onInternalFace() const noexcept;
555 
556  //- Is the particle on a boundary face?
557  inline bool onBoundaryFace() const noexcept;
558 
559  //- Return the index of patch that the particle is on
560  inline label patch() const;
562  //- Return current particle position
563  inline vector position() const;
564 
565  //- Reset particle data
566  inline void reset();
567 
568 
569  // Track
570 
571  //- Track along the displacement for a given fraction of the overall
572  // step. End when the track is complete, or when a boundary is hit.
573  // On exit, stepFraction_ will have been incremented to the current
574  // position, and facei_ will be set to the index of the boundary
575  // face that was hit, or -1 if the track completed within a cell.
576  // The proportion of the displacement still to be completed is
577  // returned.
578  scalar track
579  (
580  const vector& displacement,
581  const scalar fraction
582  );
583 
584  //- As particle::track, but also stops on internal faces.
585  scalar trackToFace
586  (
587  const vector& displacement,
588  const scalar fraction
589  );
590 
591  //- As particle::trackToFace, but also stops on tet triangles. On
592  // exit, tetTriI is set to the index of the tet triangle that was
593  // hit, or -1 if the end position was reached.
594  scalar trackToTri
595  (
596  const vector& displacement,
597  const scalar fraction,
598  label& tetTriI
599  );
600 
601  //- As particle::trackToTri, but for stationary meshes
602  scalar trackToStationaryTri
603  (
604  const vector& displacement,
605  const scalar fraction,
606  label& tetTriI
607  );
608 
609  //- As particle::trackToTri, but for moving meshes
610  scalar trackToMovingTri
611  (
612  const vector& displacement,
613  const scalar fraction,
614  label& tetTriI
615  );
616 
617  //- Hit the current face. If the current face is internal than this
618  // crosses into the next cell. If it is a boundary face then this will
619  // interact the particle with the relevant patch.
620  template<class TrackCloudType>
621  void hitFace
622  (
623  const vector& direction,
624  TrackCloudType& cloud,
626  );
627 
628  //- Convenience function. Combines trackToFace and hitFace
629  template<class TrackCloudType>
630  void trackToAndHitFace
631  (
632  const vector& direction,
633  const scalar fraction,
634  TrackCloudType& cloud,
636  );
637 
638  //- Get the displacement from the mesh centre. Used to correct the
639  // particle position in cases with reduced dimensionality. Returns a
640  // zero vector for three-dimensional cases.
642 
643 
644  // Patch data
645 
646  //- Get the normal and velocity of the current patch location
647  inline void patchData(vector& n, vector& U) const;
648 
649 
650  // Transformations
651 
652  //- Transform the physical properties of the particle
653  // according to the given transformation tensor
654  virtual void transformProperties(const tensor& T);
655 
656  //- Transform the physical properties of the particle
657  // according to the given separation vector
658  virtual void transformProperties(const vector& separation);
659 
660 
661  // Parallel transfer
662 
663  //- Convert global addressing to the processor patch local equivalents
665 
666  //- Convert processor patch addressing to the global equivalents
667  // and set the celli to the face-neighbour
668  void correctAfterParallelTransfer(const label patchi, trackingData& td);
669 
670 
671  // Interaction list referral
672 
673  //- Break the topology and store the particle position so that the
674  // particle can be referred.
676  (
678  );
679 
680  //- Correct the topology after referral. The particle may still be
681  // outside the stored tet and therefore not track-able.
682  void correctAfterInteractionListReferral(const label celli);
683 
684 
685  // Decompose and reconstruct
686 
687  //- Return the tet point appropriate for decomposition or reconstruction
688  // to or from the given mesh.
689  label procTetPt
690  (
691  const polyMesh& procMesh,
692  const label procCell,
693  const label procTetFace
694  ) const;
695 
696 
697  // Mapping
698 
699  //- Map after a topology change
700  void autoMap(const vector& position, const mapPolyMesh& mapper);
701 
702  //- Set the addressing based on the provided position
703  void relocate(const point& position, const label celli = -1);
704 
705 
706 
707  // I-O
708 
709  //- Write the name representation to stream
710  template<class Type>
711  static void writePropertyName
712  (
713  Ostream& os,
714  const word& name,
715  const word& delim
716  );
717 
718  //- Write a named particle property to stream,
719  //- optionally filtered based on its name
720  template<class Type>
721  static void writeProperty
722  (
723  Ostream& os,
724  const word& name,
725  const Type& value,
726  const bool nameOnly,
727  const word& delim,
728  const wordRes& filters = wordRes::null()
729  );
730 
731  //- Write a named particle property list to stream,
732  //- optionally filtered based on its name
733  template<class Type>
734  static void writeProperty
735  (
736  Ostream& os,
737  const word& name,
738  const Field<Type>& values,
739  const bool nameOnly,
740  const word& delim,
741  const wordRes& filters = wordRes::null()
742  );
743 
744  //- Read the fields associated with the owner cloud
745  template<class TrackCloudType>
746  static void readFields(TrackCloudType& c);
747 
748  //- Write the fields associated with the owner cloud
749  template<class TrackCloudType>
750  static void writeFields(const TrackCloudType& c);
751 
752  //- Write individual particle properties to stream
753  void writeProperties
754  (
755  Ostream& os,
756  const wordRes& filters,
757  const word& delim,
758  const bool namesOnly
759  ) const;
760 
761  //- Read particle fields as objects from the obr registry
762  template<class CloudType>
763  static void readObjects(CloudType& c, const objectRegistry& obr);
764 
765  //- Write particle fields as objects into the obr registry
766  // Always writes "position", not "coordinate"
767  template<class CloudType>
768  static void writeObjects(const CloudType& c, objectRegistry& obr);
769 
770  //- Write the particle barycentric coordinates and cell info
771  void writeCoordinates(Ostream& os) const;
772 
773  //- Write the particle position and cell id
774  virtual void writePosition(Ostream& os) const;
775 
776 
777  // Friend Operators
778 
779  friend Ostream& operator<<(Ostream&, const particle&);
780 
781  friend bool operator==(const particle& pA, const particle& pB);
782 
783  friend bool operator!=(const particle& pA, const particle& pB);
784 };
785 
786 
787 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
788 
789 } // End namespace Foam
790 
791 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
792 
793 #include "particleI.H"
794 
795 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
796 
797 #ifdef NoRepository
798  #include "particleTemplates.C"
799 #endif
800 
801 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
802 
803 #endif
804 
805 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:45
autoPtr< particle > operator()(Istream &is) const
Definition: particle.H:580
virtual void writePosition(Ostream &os) const
Write the particle position and cell id.
Definition: particleIO.C:257
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:1066
trackingData(const TrackCloudType &cloud)
Definition: particle.H:114
uint8_t direction
Definition: direction.H:46
label tetFace() const noexcept
Return current tet face particle is in.
Definition: particleI.H:134
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:1035
label tetPt() const noexcept
Return current tet face particle is in.
Definition: particleI.H:146
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:540
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition: particleI.H:116
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:1178
Vector-tensor class used to perform translations and rotations in 3D space.
label origId() const noexcept
Return the particle ID on the originating processor.
Definition: particleI.H:194
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
static void writeProperty(Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
Write a named particle property to stream, optionally filtered based on its name. ...
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 also stops on internal faces.
Definition: particle.C:639
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but also stops on tet triangles. On.
Definition: particle.C:1018
bool onFace() const noexcept
Is the particle on a face?
Definition: particleI.H:259
static bool writeLagrangianCoordinates
Write particle coordinates file (v1712 and later) Default is true.
Definition: particle.H:466
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:206
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:1050
static void writePropertyName(Ostream &os, const word &name, const word &delim)
Write the name representation to stream.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
void prepareForInteractionListReferral(const vectorTensorTransform &transform)
Break the topology and store the particle position so that the.
Definition: particle.C:1115
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Old particle positions content for OpenFOAM-1706 and earlier.
Definition: particle.H:122
Base particle class.
Definition: particle.H:69
vector normal() const
The (unit) normal of the tri on tetFacei_ for the current tet.
Definition: particleI.H:253
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
void hitBoundaryFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Dispatch function for boundary face interaction. Calls one of.
void readData(Istream &is, point &position, const bool readFields, const bool newFormat, const bool doLocate)
Read particle from stream. Optionally (for old format) return.
Definition: particleIO.C:69
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:297
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
TypeName("particle")
Runtime type information.
void hitFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
void reset()
Reset particle data.
Definition: particleI.H:289
Generic templated field type.
Definition: Field.H:62
Intrusive doubly-linked list.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
bool onInternalFace() const noexcept
Is the particle on an internal face?
Definition: particleI.H:265
A class for handling words, derived from Foam::string.
Definition: word.H:63
bool switchProcessor
Flag to switch processor.
Definition: particle.H:104
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:507
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
tetIndices currentTetIndices() const noexcept
Return indices of the current tet that the particle occupies.
Definition: particleI.H:236
void writeCoordinates(Ostream &os) const
Write the particle barycentric coordinates and cell info.
Definition: particleIO.C:238
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition: particle.C:1140
virtual ~particle()=default
Destructor.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition: particle.C:1058
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:78
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:242
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a cyclicAMIPatch.
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:170
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1204
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition: particle.H:552
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:96
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
label origProc() const noexcept
Return the originating processor ID.
Definition: particleI.H:182
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:77
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:824
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:110
U
Definition: pEqn.H:72
static bool writeLagrangianPositions
Write particle positions file (v1706 format and earlier) Default is true (disable in etc/controlDict)...
Definition: particle.H:472
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:109
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
iNew(const polyMesh &mesh)
Definition: particle.H:575
const dimensionedScalar c
Speed of light in a vacuum.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
void relocate(const point &position, const label celli=-1)
Set the addressing based on the provided position.
Definition: particle.C:1220
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
label n
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:455
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:695
bool onBoundaryFace() const noexcept
Is the particle on a boundary face?
Definition: particleI.H:271
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
#define DefinePropertyList(str)
Define a static &#39;propertyList&#39; for particle properties.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition: particleIO.C:213
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:620
Tensor of scalars, i.e. Tensor<scalar>.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
vector position() const
Return current particle position.
Definition: particleI.H:283
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:277
Namespace for OpenFOAM.
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:228