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 OpenFOAM Foundation
9  Copyright (C) 2017-2019 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 member 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 factor by which the displacement is increased when passing
90  // through negative space. This should be slightly bigger than one.
91  // This is needed as a straight trajectory can form a closed loop
92  // through regions of overlapping positive and negative space, leading
93  // to a track which never ends. By increasing the rate of displacement
94  // through negative regions, the change in track fraction over this
95  // part of the loop no longer exactly cancels the change over the
96  // positive part, and the track therefore terminates.
97  static const scalar negativeSpaceDisplacementFactor;
98 
99 
100 public:
101 
102  class trackingData
103  {
104  public:
105 
106  // Public data
107 
108  //- Flag to switch processor
109  bool switchProcessor;
110 
111  //- Flag to indicate whether to keep particle (false = delete)
112  bool keepParticle;
113 
114 
115  // Constructor
116  template <class TrackCloudType>
117  trackingData(const TrackCloudType& cloud)
118  {}
119  };
120 
121 
122  //- Old particle positions content for OpenFOAM-1706 and earlier
123  struct positionsCompat1706
124  {
126  label celli;
127  label facei;
128  scalar stepFraction;
129  label tetFacei;
130  label tetPti;
131  label origProc;
132  label origId;
133  };
134 
135 
136 private:
137 
138  // Private data
139 
140  //- Reference to the polyMesh database
141  const polyMesh& mesh_;
142 
143  //- Coordinates of particle
144  barycentric coordinates_;
145 
146  //- Index of the cell it is in
147  label celli_;
148 
149  //- Index of the face that owns the decomposed tet that the
150  // particle is in
151  label tetFacei_;
152 
153  //- Index of the point on the face that defines the decomposed
154  // tet that the particle is in. Relative to the face base
155  // point.
156  label tetPti_;
157 
158  //- Face index if the particle is on a face otherwise -1
159  label facei_;
160 
161  //- Fraction of time-step completed
162  scalar stepFraction_;
163 
164  //- Originating processor id
165  label origProc_;
166 
167  //- Local particle id on originating processor
168  label origId_;
169 
170 
171  // Private Member Functions
172 
173  // Tetrahedra functions
174 
175  //- Get the vertices of the current tet
176  inline void stationaryTetGeometry
177  (
178  vector& centre,
179  vector& base,
180  vector& vertex1,
181  vector& vertex2
182  ) const;
183 
184  //- Get the transformation associated with the current tet. This
185  // will convert a barycentric position within the tet to a
186  // Cartesian position in the global coordinate system. The
187  // conversion is x = A & y, where x is the Cartesian position, y is
188  // the barycentric position and A is the transformation tensor.
189  inline barycentricTensor stationaryTetTransform() const;
190 
191  //- Get the reverse transform associated with the current tet. The
192  // conversion is detA*y = (x - centre) & T. The variables x, y and
193  // centre have the same meaning as for the forward transform. T is
194  // the transposed inverse of the forward transform tensor, A,
195  // multiplied by its determinant, detA. This separation allows
196  // the barycentric tracking algorithm to function on inverted or
197  // degenerate tetrahedra.
198  void stationaryTetReverseTransform
199  (
200  vector& centre,
201  scalar& detA,
203  ) const;
204 
205  //- Get the vertices of the current moving tet. Two values are
206  // returned for each vertex. The first is a constant, and the
207  // second is a linear coefficient of the track fraction.
208  inline void movingTetGeometry
209  (
210  const scalar endStepFraction,
211  Pair<vector>& centre,
212  Pair<vector>& base,
213  Pair<vector>& vertex1,
214  Pair<vector>& vertex2
215  ) const;
216 
217  //- Get the transformation associated with the current, moving, tet.
218  // This is of the same form as for the static case. As with the
219  // moving geometry, a linear function of the tracking fraction is
220  // returned for each component.
221  inline Pair<barycentricTensor> movingTetTransform
222  (
223  const scalar endStepFraction
224  ) const;
225 
226  //- Get the reverse transformation associated with the current,
227  // moving, tet. This is of the same form as for the static case. As
228  // with the moving geometry, a function of the tracking fraction is
229  // returned for each component. The functions are higher order than
230  // for the forward transform; the determinant is cubic, and the
231  // tensor is quadratic.
232  void movingTetReverseTransform
233  (
234  const scalar endStepFraction,
235  Pair<vector>& centre,
236  FixedList<scalar, 4>& detA,
238  ) const;
239 
240 
241  // Transformations
242 
243  //- Reflection transform. Corrects the coordinates when the particle
244  // moves between two tets which share a base vertex, but for which
245  // the other two non cell-centre vertices are reversed. All hits
246  // which retain the same face behave this way, as do face hits.
247  void reflect();
248 
249  //- Rotation transform. Corrects the coordinates when the particle
250  // moves between two tets with different base vertices, but are
251  // otherwise similarly oriented. Hits which change the face within
252  // the cell make use of both this and the reflect transform.
253  void rotate(const bool direction);
254 
255 
256  // Topology changes
257 
258  //- Change tet within a cell. Called after a triangle is hit.
259  void changeTet(const label tetTriI);
260 
261  //- Change tet face within a cell. Called by changeTet.
262  void changeFace(const label tetTriI);
263 
264  //- Change cell. Called when the particle hits an internal face.
265  void changeCell();
266 
267  //- Put the particle on the lowest indexed patch for the current set
268  // of coincident faces. In the case of an ACMI-wall pair, this will
269  // move the particle from the wall face to the ACMI face, because
270  // ACMI patches are always listed before their associated non-
271  // overlapping patch.
272  void changeToMasterPatch();
273 
274 
275  // Geometry changes
276 
277  //- Locate the particle at the given position
278  void locate
279  (
280  const vector& position,
281  const vector* direction,
282  const label celli,
283  const bool boundaryFail,
284  const string& boundaryMsg
285  );
286 
287 
288 protected:
289 
290  // Patch interactions
291 
292  //- Overridable function to handle the particle hitting a patch.
293  // Executed before other patch-hitting functions.
294  template<class TrackCloudType>
295  bool hitPatch(TrackCloudType&, trackingData&);
296 
297  //- Overridable function to handle the particle hitting a wedgePatch
298  template<class TrackCloudType>
299  void hitWedgePatch(TrackCloudType&, trackingData&);
300 
301  //- Overridable function to handle the particle hitting a
302  // symmetryPlanePatch
303  template<class TrackCloudType>
304  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
305 
306  //- Overridable function to handle the particle hitting a symmetryPatch
307  template<class TrackCloudType>
308  void hitSymmetryPatch(TrackCloudType&, trackingData&);
309 
310  //- Overridable function to handle the particle hitting a cyclicPatch
311  template<class TrackCloudType>
312  void hitCyclicPatch(TrackCloudType&, trackingData&);
313 
314  //- Overridable function to handle the particle hitting a cyclicAMIPatch
315  template<class TrackCloudType>
316  void hitCyclicAMIPatch(TrackCloudType&, trackingData&, const vector&);
317 
318  //- Overridable function to handle the particle hitting a
319  // cyclicACMIPatch
320  template<class TrackCloudType>
321  void hitCyclicACMIPatch(TrackCloudType&, trackingData&, const vector&);
322 
323  //- Overridable function to handle the particle hitting a processorPatch
324  template<class TrackCloudType>
325  void hitProcessorPatch(TrackCloudType&, trackingData&);
326 
327  //- Overridable function to handle the particle hitting a wallPatch
328  template<class TrackCloudType>
329  void hitWallPatch(TrackCloudType&, trackingData&);
330 
331 
332 public:
333 
334  // Static data members
335 
336  //- Runtime type information
337  TypeName("particle");
338 
339  //- String representation of properties
341  (
342  "(coordinatesa coordinatesb coordinatesc coordinatesd) "
343  "celli tetFacei tetPti facei stepFraction origProc origId"
344  );
345 
346  //- Cumulative particle counter - used to provide unique ID
347  static label particleCount_;
348 
349  //- Write particle coordinates file (v1712 and later)
350  //- Default is true
351  static bool writeLagrangianCoordinates;
352 
353  //- Write particle positions file (v1706 format and earlier)
354  //- Default is true (disable in etc/controlDict)
355  static bool writeLagrangianPositions;
356 
357 
358  // Constructors
359 
360  //- Construct from components
361  particle
362  (
363  const polyMesh& mesh,
364  const barycentric& coordinates,
365  const label celli,
366  const label tetFacei,
367  const label tetPti
368  );
369 
370  //- Construct from a position and a cell.
371  // Searches for the rest of the required topology.
372  particle
373  (
374  const polyMesh& mesh,
375  const vector& position,
376  const label celli = -1
377  );
378 
379  //- Construct from position components
380  particle
381  (
382  const polyMesh& mesh,
383  const vector& position,
384  const label celli,
385  const label tetFacei,
386  const label tetPti,
387  const bool doLocate = true
388  );
389 
390 
391  //- Construct from Istream
392  particle
393  (
394  const polyMesh& mesh,
395  Istream&,
396  bool readFields = true,
397  bool newFormat = true
398  );
399 
400  //- Construct as a copy
401  particle(const particle& p);
402 
403  //- Construct as a copy with reference to a new mesh
404  particle(const particle& p, const polyMesh& mesh);
405 
406  //- Construct a clone
407  virtual autoPtr<particle> clone() const
408  {
409  return autoPtr<particle>::New(*this);
410  }
411 
412  //- Factory class to read-construct particles (for parallel transfer)
413  class iNew
414  {
415  const polyMesh& mesh_;
416 
417  public:
418 
419  iNew(const polyMesh& mesh)
420  :
421  mesh_(mesh)
422  {}
423 
425  {
426  return autoPtr<particle>::New(mesh_, is, true);
427  }
428  };
429 
430 
431  //- Destructor
432  virtual ~particle() = default;
433 
434 
435  // Member Functions
436 
437  // Access
438 
439  //- Get unique particle creation id
440  inline label getNewParticleID() const;
441 
442  //- Return the mesh database
443  inline const polyMesh& mesh() const;
444 
445  //- Return current particle coordinates
446  inline const barycentric& coordinates() const;
447 
448  //- Return current cell particle is in
449  inline label cell() const;
450 
451  //- Return current cell particle is in for manipulation
452  inline label& cell();
453 
454  //- Return current tet face particle is in
455  inline label tetFace() const;
456 
457  //- Return current tet face particle is in for manipulation
458  inline label& tetFace();
459 
460  //- Return current tet face particle is in
461  inline label tetPt() const;
462 
463  //- Return current tet face particle is in for manipulation
464  inline label& tetPt();
465 
466  //- Return current face particle is on otherwise -1
467  inline label face() const;
468 
469  //- Return current face particle is on for manipulation
470  inline label& face();
471 
472  //- Return the fraction of time-step completed
473  inline scalar stepFraction() const;
474 
475  //- Return the fraction of time-step completed
476  inline scalar& stepFraction();
477 
478  //- Return the originating processor ID
479  inline label origProc() const;
480 
481  //- Return the originating processor ID
482  inline label& origProc();
483 
484  //- Return the particle ID on the originating processor
485  inline label origId() const;
486 
487  //- Return the particle ID on the originating processor
488  inline label& origId();
489 
490 
491  // Check
492 
493  //- Return the step fraction change within the overall time-step.
494  // Returns the start value and the change as a scalar pair. Always
495  // return Pair<scalar>(0, 1), unless sub-cycling is in effect, in
496  // which case the values will reflect the span of the sub-cycle
497  // within the time-step.
498  inline Pair<scalar> stepFractionSpan() const;
499 
500  //- Return the current fraction within the timestep. This differs
501  // from the stored step fraction due to sub-cycling.
502  inline scalar currentTimeFraction() const;
503 
504  //- Return the indices of the current tet that the
505  // particle occupies.
506  inline tetIndices currentTetIndices() const;
507 
508  //- Return the current tet transformation tensor
509  inline barycentricTensor currentTetTransform() const;
510 
511  //- The (unit) normal of the tri on tetFacei_ for the current tet.
512  inline vector normal() const;
513 
514  //- Is the particle on a face?
515  inline bool onFace() const;
516 
517  //- Is the particle on an internal face?
518  inline bool onInternalFace() const;
519 
520  //- Is the particle on a boundary face?
521  inline bool onBoundaryFace() const;
522 
523  //- Return the index of patch that the particle is on
524  inline label patch() const;
525 
526  //- Return current particle position
527  inline vector position() const;
528 
529 
530  // Track
531 
532  //- Track along the displacement for a given fraction of the overall
533  // step. End when the track is complete, or when a boundary is hit.
534  // On exit, stepFraction_ will have been incremented to the current
535  // position, and facei_ will be set to the index of the boundary
536  // face that was hit, or -1 if the track completed within a cell.
537  // The proportion of the displacement still to be completed is
538  // returned.
539  scalar track
540  (
541  const vector& displacement,
542  const scalar fraction
543  );
544 
545  //- As particle::track, but also stops on internal faces.
546  scalar trackToFace
547  (
548  const vector& displacement,
549  const scalar fraction
550  );
551 
552  //- As particle::trackToFace, but also stops on tet triangles. On
553  // exit, tetTriI is set to the index of the tet triangle that was
554  // hit, or -1 if the end position was reached.
555  scalar trackToTri
556  (
557  const vector& displacement,
558  const scalar fraction,
559  label& tetTriI
560  );
561 
562  //- As particle::trackToTri, but for stationary meshes
563  scalar trackToStationaryTri
564  (
565  const vector& displacement,
566  const scalar fraction,
567  label& tetTriI
568  );
569 
570  //- As particle::trackToTri, but for moving meshes
571  scalar trackToMovingTri
572  (
573  const vector& displacement,
574  const scalar fraction,
575  label& tetTriI
576  );
577 
578  //- Hit the current face. If the current face is internal than this
579  // crosses into the next cell. If it is a boundary face then this will
580  // interact the particle with the relevant patch.
581  template<class TrackCloudType>
582  void hitFace
583  (
584  const vector& direction,
585  TrackCloudType& cloud,
586  trackingData& td
587  );
588 
589  //- Convenience function. Combines trackToFace and hitFace
590  template<class TrackCloudType>
591  void trackToAndHitFace
592  (
593  const vector& direction,
594  const scalar fraction,
595  TrackCloudType& cloud,
596  trackingData& td
597  );
598 
599  //- Get the displacement from the mesh centre. Used to correct the
600  // particle position in cases with reduced dimensionality. Returns a
601  // zero vector for three-dimensional cases.
603 
604 
605  // Patch data
606 
607  //- Get the normal and velocity of the current patch location
608  inline void patchData(vector& n, vector& U) const;
609 
610 
611  // Transformations
612 
613  //- Transform the physical properties of the particle
614  // according to the given transformation tensor
615  virtual void transformProperties(const tensor& T);
616 
617  //- Transform the physical properties of the particle
618  // according to the given separation vector
619  virtual void transformProperties(const vector& separation);
620 
621 
622  // Parallel transfer
623 
624  //- Convert global addressing to the processor patch local equivalents
626 
627  //- Convert processor patch addressing to the global equivalents
628  // and set the celli to the face-neighbour
629  void correctAfterParallelTransfer(const label patchi, trackingData& td);
630 
631 
632  // Interaction list referral
633 
634  //- Break the topology and store the particle position so that the
635  // particle can be referred.
637  (
639  );
640 
641  //- Correct the topology after referral. The particle may still be
642  // outside the stored tet and therefore not track-able.
643  void correctAfterInteractionListReferral(const label celli);
644 
645 
646  // Decompose and reconstruct
647 
648  //- Return the tet point appropriate for decomposition or reconstruction
649  // to or from the given mesh.
650  label procTetPt
651  (
652  const polyMesh& procMesh,
653  const label procCell,
654  const label procTetFace
655  ) const;
656 
657 
658  // Mapping
659 
660  //- Map after a topology change
661  void autoMap(const vector& position, const mapPolyMesh& mapper);
662 
663  //- Set the addressing based on the provided position
664  void relocate(const point& position, const label celli = -1);
665 
666 
667  // I-O
668 
669  //- Write the name representation to stream
670  template<class Type>
671  static void writePropertyName
672  (
673  Ostream& os,
674  const word& name,
675  const word& delim
676  );
677 
678  //- Write a named particle property to stream,
679  //- optionally filtered based on its name
680  template<class Type>
681  static void writeProperty
682  (
683  Ostream& os,
684  const word& name,
685  const Type& value,
686  const bool nameOnly,
687  const word& delim,
688  const wordRes& filters = wordRes::null()
689  );
690 
691  //- Write a named particle property list to stream,
692  //- optionally filtered based on its name
693  template<class Type>
694  static void writeProperty
695  (
696  Ostream& os,
697  const word& name,
698  const Field<Type>& values,
699  const bool nameOnly,
700  const word& delim,
701  const wordRes& filters = wordRes::null()
702  );
703 
704  //- Read the fields associated with the owner cloud
705  template<class TrackCloudType>
706  static void readFields(TrackCloudType& c);
707 
708  //- Write the fields associated with the owner cloud
709  template<class TrackCloudType>
710  static void writeFields(const TrackCloudType& c);
711 
712  //- Write individual particle properties to stream
713  void writeProperties
714  (
715  Ostream& os,
716  const wordRes& filters,
717  const word& delim,
718  const bool namesOnly
719  ) const;
720 
721  //- Read particle fields as objects from the obr registry
722  template<class CloudType>
723  static void readObjects(CloudType& c, const objectRegistry& obr);
724 
725  //- Write particle fields as objects into the obr registry
726  // Always writes "position", not "coordinate"
727  template<class CloudType>
728  static void writeObjects(const CloudType& c, objectRegistry& obr);
729 
730  //- Write the particle barycentric coordinates and cell info
731  void writeCoordinates(Ostream& os) const;
732 
733  //- Write the particle position and cell id
734  virtual void writePosition(Ostream& os) const;
735 
736 
737  // Friend Operators
738 
739  friend Ostream& operator<<(Ostream&, const particle&);
740 
741  friend bool operator==(const particle& pA, const particle& pB);
742 
743  friend bool operator!=(const particle& pA, const particle& pB);
744 };
745 
746 
747 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
748 
749 } // End namespace Foam
750 
751 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
752 
753 #include "particleI.H"
754 
755 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
756 
757 #ifdef NoRepository
758  #include "particleTemplates.C"
759 #endif
760 
761 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
762 
763 #endif
764 
765 // ************************************************************************* //
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:573
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:952
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:1034
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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:116
Foam::particle::onFace
bool onFace() const
Is the particle on a face?
Definition: particleI.H:285
Foam::particle::iNew
Factory class to read-construct particles (for parallel transfer)
Definition: particle.H:412
Foam::particle::positionsCompat1706::tetFacei
label tetFacei
Definition: particle.H:128
Foam::particle::origId
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:216
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:192
Foam::particle::autoMap
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1098
Foam::particle::trackToMovingTri
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:784
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:156
Foam::particle::positionsCompat1706::position
vector position
Definition: particle.H:124
Foam::particle::getNewParticleID
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:118
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:260
Foam::particle::relocate
void relocate(const point &position, const label celli=-1)
Set the addressing based on the provided position.
Definition: particle.C:1114
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:519
faceList.H
Foam::particle::positionsCompat1706::tetPti
label tetPti
Definition: particle.H:129
Foam::particle::particleCount_
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:343
Foam::particle::onInternalFace
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:291
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
polyMeshTetDecomposition.H
Foam::particle::writeLagrangianPositions
static bool writeLagrangianPositions
Definition: particle.H:354
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:131
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:138
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:618
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
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:127
Foam::particle::trackToStationaryTri
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:673
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:130
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
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:531
Foam::particle::position
vector position() const
Return current particle position.
Definition: particleI.H:309
Foam::Barycentric< scalar >
Foam::particle::trackingData::keepParticle
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:111
Foam::particle::trackingData::switchProcessor
bool switchProcessor
Flag to switch processor.
Definition: particle.H:108
Foam::particle::patch
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:303
Foam::particle::mesh
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:132
Foam::particle::deviationFromMeshCentre
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:929
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:252
Foam::particle::clone
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:406
Foam::particle::iNew::operator()
autoPtr< particle > operator()(Istream &is) const
Definition: particle.H:423
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:204
Foam::particle::correctAfterParallelTransfer
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:960
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:1072
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:637
particleI.H
Foam::particle::onBoundaryFace
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:297
Foam::particle::operator!=
friend bool operator!=(const particle &pA, const particle &pB)
pointField.H
Foam::particle::writeLagrangianCoordinates
static bool writeLagrangianCoordinates
Definition: particle.H:350
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:122
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:912
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:944
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:1009
Foam::particle::stepFractionSpan
Pair< scalar > stepFractionSpan() const
Return the step fraction change within the overall time-step.
Definition: particleI.H:228
Foam::direction
uint8_t direction
Definition: direction.H:47
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:578
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:168
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
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:279
Foam::particle::cell
label cell() const
Return current cell particle is in.
Definition: particleI.H:144
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:418
vectorTensorTransform.H
Foam::particle::face
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:180
Foam::particle::trackingData
Definition: particle.H:101
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:505
Foam::particle::positionsCompat1706::facei
label facei
Definition: particle.H:126
IDLList.H
Intrusive doubly-linked list.
Foam::particle::positionsCompat1706::celli
label celli
Definition: particle.H:125
Foam::particle::patchData
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:315
Foam::particle::currentTetTransform
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:266