particle.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::particle
26 
27 Description
28  Base particle class
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef particle_H
33 #define particle_H
34 
35 #include "vector.H"
36 #include "barycentric.H"
37 #include "barycentricTensor.H"
38 #include "Cloud.H"
39 #include "IDLList.H"
40 #include "pointField.H"
41 #include "faceList.H"
42 #include "OFstream.H"
43 #include "tetPointRef.H"
44 #include "FixedList.H"
46 #include "particleMacros.H"
47 #include "vectorTensorTransform.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // Forward declaration of classes
55 class particle;
56 
57 class polyPatch;
58 
59 class cyclicPolyPatch;
60 class cyclicAMIPolyPatch;
61 class cyclicACMIPolyPatch;
62 class cyclicRepeatAMIPolyPatch;
63 class processorPolyPatch;
64 class symmetryPlanePolyPatch;
65 class symmetryPolyPatch;
66 class wallPolyPatch;
67 class wedgePolyPatch;
68 
69 // Forward declaration of friend functions and operators
70 
71 Ostream& operator<<
72 (
73  Ostream&,
74  const particle&
75 );
76 
77 bool operator==(const particle&, const particle&);
78 
79 bool operator!=(const particle&, const particle&);
80 
81 /*---------------------------------------------------------------------------*\
82  Class Particle Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 class particle
86 :
87  public IDLList<particle>::link
88 {
89  // Private member data
90 
91  //- Size in bytes of the position data
92  static const std::size_t sizeofPosition_;
93 
94  //- Size in bytes of the fields
95  static const std::size_t sizeofFields_;
96 
97  //- The factor by which the displacement is increased when passing
98  // through negative space. This should be slightly bigger than one.
99  // This is needed as a straight trajectory can form a closed loop
100  // through regions of overlapping positive and negative space, leading
101  // to a track which never ends. By increasing the rate of displacement
102  // through negative regions, the change in track fraction over this
103  // part of the loop no longer exactly cancels the change over the
104  // positive part, and the track therefore terminates.
105  static const scalar negativeSpaceDisplacementFactor;
106 
107 
108 public:
110  class trackingData
111  {
112  public:
113 
114  // Public data
115 
116  //- Flag to switch processor
117  bool switchProcessor;
118 
119  //- Flag to indicate whether to keep particle (false = delete)
120  bool keepParticle;
121 
122 
123  // Constructor
124  template <class TrackCloudType>
125  trackingData(const TrackCloudType& cloud)
126  {}
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  //- Originating processor id
159  label origProc_;
160 
161  //- Local particle id on originating processor
162  label origId_;
163 
164 
165 private:
166 
167  // Private Member Functions
168 
169  // Tetrahedra functions
170 
171  //- Get the vertices of the current tet
172  inline void stationaryTetGeometry
173  (
174  vector& centre,
175  vector& base,
176  vector& vertex1,
177  vector& vertex2
178  ) const;
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  label celli,
279  const bool boundaryFail,
280  const string boundaryMsg
281  );
282 
283 
284 protected:
285 
286  // Patch interactions
287 
288  //- Overridable function to handle the particle hitting a patch.
289  // Executed before other patch-hitting functions.
290  template<class TrackCloudType>
291  bool hitPatch(TrackCloudType&, trackingData&);
292 
293  //- Overridable function to handle the particle hitting a wedgePatch
294  template<class TrackCloudType>
295  void hitWedgePatch(TrackCloudType&, trackingData&);
296 
297  //- Overridable function to handle the particle hitting a
298  // symmetryPlanePatch
299  template<class TrackCloudType>
300  void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
301 
302  //- Overridable function to handle the particle hitting a symmetryPatch
303  template<class TrackCloudType>
304  void hitSymmetryPatch(TrackCloudType&, trackingData&);
305 
306  //- Overridable function to handle the particle hitting a cyclicPatch
307  template<class TrackCloudType>
308  void hitCyclicPatch(TrackCloudType&, trackingData&);
309 
310  //- Overridable function to handle the particle hitting a cyclicAMIPatch
311  template<class TrackCloudType>
312  void hitCyclicAMIPatch(TrackCloudType&, trackingData&, const vector&);
313 
314  //- Overridable function to handle the particle hitting a
315  // cyclicACMIPatch
316  template<class TrackCloudType>
317  void hitCyclicACMIPatch(TrackCloudType&, trackingData&, const vector&);
318 
319  //- Overridable function to handle the particle hitting an
320  // cyclicRepeatAMIPolyPatch
321  template<class TrackCloudType>
323  (
324  TrackCloudType&,
325  trackingData&,
326  const vector&
327  );
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 provode unique ID
353  static label particleCount_;
354 
355 
356  // Constructors
357 
358  //- Construct from components
359  particle
360  (
361  const polyMesh& mesh,
362  const barycentric& coordinates,
363  const label celli,
364  const label tetFacei,
365  const label tetPti
366  );
367 
368  //- Construct from a position and a cell, searching for the rest of the
369  // required topology
370  particle
371  (
372  const polyMesh& mesh,
373  const vector& position,
374  const label celli
375  );
376 
377  //- Construct from Istream
378  particle(const polyMesh& mesh, Istream&, bool readFields = true);
379 
380  //- Construct as a copy
381  particle(const particle& p);
382 
383  //- Construct as a copy with references to a new mesh
384  particle(const particle& p, const polyMesh& mesh);
385 
386  //- Construct a clone
387  virtual autoPtr<particle> clone() const
388  {
389  return autoPtr<particle>(new particle(*this));
390  }
391 
392  //- Factory class to read-construct particles used for
393  // parallel transfer
394  class iNew
395  {
396  const polyMesh& mesh_;
397 
398  public:
400  iNew(const polyMesh& mesh)
401  :
402  mesh_(mesh)
403  {}
405  autoPtr<particle> operator()(Istream& is) const
406  {
407  return autoPtr<particle>(new particle(mesh_, is, true));
408  }
409  };
410 
411 
412  //- Destructor
413  virtual ~particle()
414  {}
415 
416 
417  // Member Functions
418 
419  // Access
420 
421  //- Get unique particle creation id
422  inline label getNewParticleID() const;
423 
424  //- Return the mesh database
425  inline const polyMesh& mesh() const;
426 
427  //- Return current particle coordinates
428  inline const barycentric& coordinates() const;
429 
430  //- Return current cell particle is in
431  inline label cell() const;
432 
433  //- Return current tet face particle is in
434  inline label tetFace() const;
435 
436  //- Return current tet face particle is in
437  inline label tetPt() const;
438 
439  //- Return current face particle is on otherwise -1
440  inline label face() const;
441 
442  //- Return the fraction of time-step completed
443  inline scalar stepFraction() const;
444 
445  //- Return the fraction of time-step completed
446  inline scalar& stepFraction();
447 
448  //- Return the originating processor ID
449  inline label origProc() const;
450 
451  //- Return the originating processor ID
452  inline label& origProc();
453 
454  //- Return the particle ID on the originating processor
455  inline label origId() const;
456 
457  //- Return the particle ID on the originating processor
458  inline label& origId();
459 
460 
461  // Check
462 
463  //- Return the step fraction change within the overall time-step.
464  // Returns the start value and the change as a scalar pair. Always
465  // return Pair<scalar>(0, 1), unless sub-cycling is in effect, in
466  // which case the values will reflect the span of the sub-cycle
467  // within the time-step.
468  inline Pair<scalar> stepFractionSpan() const;
469 
470  //- Return the current fraction within the timestep. This differs
471  // from the stored step fraction due to sub-cycling.
472  inline scalar currentTimeFraction() const;
473 
474  //- Return the indices of the current tet that the
475  // particle occupies.
476  inline tetIndices currentTetIndices() const;
477 
478  //- Return the current tet transformation tensor
479  inline barycentricTensor currentTetTransform() const;
480 
481  //- Return the normal of the tri on tetFacei_ for the
482  // current tet.
483  inline vector normal() const;
484 
485  //- Is the particle on a face?
486  inline bool onFace() const;
487 
488  //- Is the particle on an internal face?
489  inline bool onInternalFace() const;
490 
491  //- Is the particle on a boundary face?
492  inline bool onBoundaryFace() const;
493 
494  //- Return the index of patch that the particle is on
495  inline label patch() const;
496 
497  //- Return current particle position
498  inline vector position() const;
499 
500 
501  // Track
502 
503  //- Track along the displacement for a given fraction of the overall
504  // step. End when the track is complete, or when a boundary is hit.
505  // On exit, stepFraction_ will have been incremented to the current
506  // position, and facei_ will be set to the index of the boundary
507  // face that was hit, or -1 if the track completed within a cell.
508  // The proportion of the displacement still to be completed is
509  // returned.
510  scalar track
511  (
512  const vector& displacement,
513  const scalar fraction
514  );
515 
516  //- As particle::track, but stops when a new cell is reached.
517  scalar trackToCell
518  (
519  const vector& displacement,
520  const scalar fraction
521  );
522 
523  //- As particle::track, but stops when a face is hit.
524  scalar trackToFace
525  (
526  const vector& displacement,
527  const scalar fraction
528  );
529 
530  //- As particle::trackToFace, but stops when a tet triangle is hit. On
531  // exit, tetTriI is set to the index of the tet triangle that was hit,
532  // or -1 if the end position was reached.
533  scalar trackToTri
534  (
535  const vector& displacement,
536  const scalar fraction,
537  label& tetTriI
538  );
539 
540  //- As particle::trackToTri, but for stationary meshes
541  scalar trackToStationaryTri
542  (
543  const vector& displacement,
544  const scalar fraction,
545  label& tetTriI
546  );
547 
548  //- As particle::trackToTri, but for moving meshes
549  scalar trackToMovingTri
550  (
551  const vector& displacement,
552  const scalar fraction,
553  label& tetTriI
554  );
555 
556  //- Hit the current face. If the current face is internal than this
557  // crosses into the next cell. If it is a boundary face then this will
558  // interact the particle with the relevant patch.
559  template<class TrackCloudType>
560  void hitFace
561  (
562  const vector& direction,
563  TrackCloudType& cloud,
564  trackingData& td
565  );
566 
567  //- Convenience function. Cobines trackToFace and hitFace
568  template<class TrackCloudType>
569  void trackToAndHitFace
570  (
571  const vector& direction,
572  const scalar fraction,
573  TrackCloudType& cloud,
574  trackingData& td
575  );
576 
577  //- Get the displacement from the mesh centre. Used to correct the
578  // particle position in cases with reduced dimensionality. Returns a
579  // zero vector for three-dimensional cases.
581 
582 
583  // Patch data
584 
585  //- Get the normal and velocity of the current patch location
586  inline void patchData(vector& n, vector& U) const;
587 
588 
589  // Transformations
590 
591  //- Transform the physical properties of the particle
592  // according to the given transformation tensor
593  virtual void transformProperties(const tensor& T);
594 
595  //- Transform the physical properties of the particle
596  // according to the given separation vector
597  virtual void transformProperties(const vector& separation);
598 
599 
600  // Parallel transfer
601 
602  //- Convert global addressing to the processor patch local equivalents
604 
605  //- Convert processor patch addressing to the global equivalents
606  // and set the celli to the face-neighbour
608 
609 
610  // Interaction list referral
611 
612  //- Break the topology and store the particle position so that the
613  // particle can be referred.
615  (
617  );
618 
619  //- Correct the topology after referral. The particle may still be
620  // outside the stored tet and therefore not track-able.
621  void correctAfterInteractionListReferral(const label celli);
622 
623 
624  // Decompose and reconstruct
625 
626  //- Return the tet point appropriate for decomposition or reconstruction
627  // to or from the given mesh.
629  (
630  const polyMesh& procMesh,
631  const label procCell,
632  const label procTetFace
633  ) const;
634 
635 
636  // Mapping
637 
638  //- Map after a topology change
639  void autoMap(const vector& position, const mapPolyMesh& mapper);
640 
641 
642  // I-O
643 
644  //- Read the fields associated with the owner cloud
645  template<class TrackCloudType>
646  static void readFields(TrackCloudType& c);
647 
648  //- Write the fields associated with the owner cloud
649  template<class TrackCloudType>
650  static void writeFields(const TrackCloudType& c);
651 
652  //- Write the particle position and cell
653  void writePosition(Ostream&) const;
654 
655 
656  // Friend Operators
657 
658  friend Ostream& operator<<(Ostream&, const particle&);
659 
660  friend bool operator==(const particle& pA, const particle& pB);
661 
662  friend bool operator!=(const particle& pA, const particle& pB);
663 };
664 
665 
666 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
667 
668 } // End namespace Foam
669 
670 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
671 
672 #include "particleI.H"
673 
674 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
675 
676 #ifdef NoRepository
677  #include "particleTemplates.C"
678 #endif
679 
680 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
681 
682 #endif
683 
684 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
label origProc() const
Return the originating processor ID.
Definition: particleI.H:178
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:988
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:130
trackingData(const TrackCloudType &cloud)
Definition: particle.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:957
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:386
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:271
friend bool operator!=(const particle &pA, const particle &pB)
Macros for adding to particle property lists.
label procTetPt(const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or reconstruction.
Definition: particle.C:1100
Vector-tensor class used to perform translations and rotations in 3D space.
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:625
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but stops when a tet triangle is hit. On.
Definition: particle.C:940
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:202
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:972
scalar trackToCell(const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
Definition: particle.C:608
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:1037
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Base particle class.
Definition: particle.H:84
label tetPt() const
Return current tet face particle is in.
Definition: particleI.H:154
friend Ostream & operator<<(Ostream &, const particle &)
vector normal() const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:253
label cell() const
Return current cell particle is in.
Definition: particleI.H:142
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:289
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.
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:160
Intrusive doubly-linked list.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
bool switchProcessor
Flag to switch processor.
Definition: particle.H:116
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition: particle.C:510
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition: particle.C:1062
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition: particle.C:980
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:240
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a cyclicAMIPatch.
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:190
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1126
bool onFace() const
Is the particle on a face?
Definition: particleI.H:259
label tetFace() const
Return current tet face particle is in.
Definition: particleI.H:148
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:116
Factory class to read-construct particles used for.
Definition: particle.H:393
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:166
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:234
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:772
label patchi
U
Definition: pEqn.H:72
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:119
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:84
const dimensionedScalar c
Speed of light in a vacuum.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
friend bool operator==(const particle &pA, const particle &pB)
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
label n
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:136
static label particleCount_
Cumulative particle counter - used to provode unique ID.
Definition: particle.H:349
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:265
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Cobines trackToFace and hitFace.
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:661
volScalarField & p
#define DefinePropertyList(str)
virtual ~particle()
Destructor.
Definition: particle.H:412
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:589
bool operator!=(const particle &, const particle &)
Definition: particle.C:1150
vector position() const
Return current particle position.
Definition: particleI.H:283
void hitCyclicRepeatAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting an.
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:277
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:226