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