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 | Copyright (C) 2011-2016 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 "Cloud.H"
37 #include "IDLList.H"
38 #include "pointField.H"
39 #include "faceList.H"
40 #include "OFstream.H"
41 #include "tetrahedron.H"
42 #include "FixedList.H"
44 #include "particleMacros.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward declaration of classes
52 class particle;
53 
54 class polyPatch;
55 
56 class cyclicPolyPatch;
57 class processorPolyPatch;
58 class symmetryPlanePolyPatch;
59 class symmetryPolyPatch;
60 class wallPolyPatch;
61 class wedgePolyPatch;
62 
63 // Forward declaration of friend functions and operators
64 
65 Ostream& operator<<
66 (
67  Ostream&,
68  const particle&
69 );
70 
71 bool operator==(const particle&, const particle&);
72 
73 bool operator!=(const particle&, const particle&);
74 
75 /*---------------------------------------------------------------------------*\
76  Class Particle Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 class particle
80 :
81  public IDLList<particle>::link
82 {
83  // Private member data
84 
85  //- Size in bytes of the position data
86  static const std::size_t sizeofPosition_;
87 
88  //- Size in bytes of the fields
89  static const std::size_t sizeofFields_;
90 
91 
92 public:
93 
94  template<class CloudType>
95  class TrackingData
96  {
97  // Private data
98 
99  //- Reference to the cloud containing (this) particle
100  CloudType& cloud_;
101 
102 
103  public:
104 
105  // Public data
107  typedef CloudType cloudType;
108 
109  //- Flag to switch processor
110  bool switchProcessor;
111 
112  //- Flag to indicate whether to keep particle (false = delete)
113  bool keepParticle;
114 
115 
116  // Constructor
118  :
119  cloud_(cloud)
120  {}
121 
122 
123  // Member functions
124 
125  //- Return a reference to the cloud
126  CloudType& cloud()
127  {
128  return cloud_;
129  }
130  };
131 
132 
133 protected:
134 
135  // Protected data
136 
137  //- Reference to the polyMesh database
138  const polyMesh& mesh_;
139 
140  //- Position of particle
142 
143  //- Index of the cell it is in
144  label celli_;
145 
146  //- Face index if the particle is on a face otherwise -1
147  label facei_;
148 
149  //- Fraction of time-step completed
150  scalar stepFraction_;
151 
152  //- Index of the face that owns the decomposed tet that the
153  // particle is in
155 
156  //- Index of the point on the face that defines the decomposed
157  // tet that the particle is in. Relative to the face base
158  // point.
159  label tetPtI_;
160 
161  //- Originating processor id
163 
164  //- Local particle id on originating processor
165  label origId_;
166 
167 
168  // Private Member Functions
169 
170  //- Find the tet tri faces between position and tet centre
171  void findTris
172  (
173  const vector& position,
175  const tetPointRef& tet,
176  const FixedList<vector, 4>& tetAreas,
177  const FixedList<label, 4>& tetPlaneBasePtIs,
178  const scalar tol
179  ) const;
180 
181  //- Find the lambda value for the line to-from across the
182  // given tri face, where p = from + lambda*(to - from)
183  inline scalar tetLambda
184  (
185  const vector& from,
186  const vector& to,
187  const label triI,
188  const vector& tetArea,
189  const label tetPlaneBasePtI,
190  const label celli,
191  const label tetFacei,
192  const label tetPtI,
193  const scalar tol
194  ) const;
195 
196  //- Find the lambda value for a moving tri face
197  inline scalar movingTetLambda
198  (
199  const vector& from,
200  const vector& to,
201  const label triI,
202  const vector& tetArea,
203  const label tetPlaneBasePtI,
204  const label celli,
205  const label tetFacei,
206  const label tetPtI,
207  const scalar tol
208  ) const;
209 
210  //- Modify the tet owner data by crossing triI
211  inline void tetNeighbour(label triI);
212 
213  //- Cross the from the given face across the given edge of the
214  // given cell to find the resulting face and tetPtI
215  inline void crossEdgeConnectedFace
216  (
217  const label& celli,
218  label& tetFacei,
219  label& tetPtI,
220  const edge& e
221  );
222 
223  //- Hit wall faces in the current cell if the
224  //- wallImpactDistance is non-zero. They may not be in
225  //- Different tets to the current.
226  template<class CloudType>
227  inline void hitWallFaces
228  (
229  const CloudType& td,
230  const vector& from,
231  const vector& to,
232  scalar& lambdaMin,
233  tetIndices& closestTetIs
234  );
235 
236 
237  // Patch interactions
238 
239  //- Overridable function to handle the particle hitting a face
240  template<class TrackData>
241  void hitFace(TrackData& td);
242 
243  //- Overridable function to handle the particle hitting a
244  // patch. Executed before other patch-hitting functions.
245  // trackFraction is passed in to allow mesh motion to
246  // interpolate in time to the correct face state.
247  template<class TrackData>
248  bool hitPatch
249  (
250  const polyPatch&,
251  TrackData& td,
252  const label patchi,
253  const scalar trackFraction,
254  const tetIndices& tetIs
255  );
256 
257  //- Overridable function to handle the particle hitting a wedgePatch
258  template<class TrackData>
259  void hitWedgePatch(const wedgePolyPatch&, TrackData& td);
260 
261  //- Overridable function to handle the particle hitting a
262  // symmetryPlanePatch
263  template<class TrackData>
265  (
266  const symmetryPlanePolyPatch&,
267  TrackData& td
268  );
269 
270  //- Overridable function to handle the particle hitting a
271  // symmetryPatch
272  template<class TrackData>
273  void hitSymmetryPatch(const symmetryPolyPatch&, TrackData& td);
274 
275  //- Overridable function to handle the particle hitting a cyclicPatch
276  template<class TrackData>
277  void hitCyclicPatch(const cyclicPolyPatch&, TrackData& td);
278 
279  //- Overridable function to handle the particle hitting a cyclicAMIPatch
280  template<class TrackData>
281  void hitCyclicAMIPatch
282  (
283  const cyclicAMIPolyPatch&,
284  TrackData& td,
285  const vector& direction
286  );
287 
288  //- Overridable function to handle the particle hitting a
289  // processorPatch
290  template<class TrackData>
291  void hitProcessorPatch(const processorPolyPatch&, TrackData& td);
292 
293  //- Overridable function to handle the particle hitting a wallPatch
294  template<class TrackData>
295  void hitWallPatch
296  (
297  const wallPolyPatch&,
298  TrackData& td,
299  const tetIndices& tetIs
300  );
301 
302  //- Overridable function to handle the particle hitting a
303  // general patch
304  template<class TrackData>
305  void hitPatch(const polyPatch&, TrackData& td);
306 
307 
308 public:
309 
310  // Static data members
311 
312  //- Runtime type information
313  TypeName("particle");
314 
315  //- String representation of properties
317  (
318  "(Px Py Pz) celli facei stepFraction "
319  "tetFacei tetPtI origProc origId"
320  );
321 
322  //- Cumulative particle counter - used to provode unique ID
323  static label particleCount_;
324 
325  //- Fraction of distance to tet centre to move a particle to
326  // 'rescue' it from a tracking problem
327  static const scalar trackingCorrectionTol;
328 
329  //- Fraction of the cell volume to use in determining tolerance values
330  // for the denominator and numerator of lambda
331  static const scalar lambdaDistanceToleranceCoeff;
332 
333  //- Minimum stepFraction tolerance
334  static const scalar minStepFractionTol;
335 
336 
337  // Constructors
338 
339  //- Construct from components
340  particle
341  (
342  const polyMesh& mesh,
343  const vector& position,
344  const label celli,
345  const label tetFacei,
346  const label tetPtI
347  );
348 
349  //- Construct from components, tetFacei_ and tetPtI_ are not
350  // supplied so they will be deduced by a search
351  particle
352  (
353  const polyMesh& mesh,
354  const vector& position,
355  const label celli,
356  bool doCellFacePt = true
357  );
358 
359  //- Construct from Istream
360  particle(const polyMesh& mesh, Istream&, bool readFields = true);
361 
362  //- Construct as a copy
363  particle(const particle& p);
364 
365  //- Construct as a copy with refernce to a new mesh
366  particle(const particle& p, const polyMesh& mesh);
367 
368  //- Construct a clone
369  virtual autoPtr<particle> clone() const
370  {
371  return autoPtr<particle>(new particle(*this));
372  }
373 
374  //- Factory class to read-construct particles used for
375  // parallel transfer
376  class iNew
377  {
378  const polyMesh& mesh_;
379 
380  public:
382  iNew(const polyMesh& mesh)
383  :
384  mesh_(mesh)
385  {}
387  autoPtr<particle> operator()(Istream& is) const
388  {
389  return autoPtr<particle>(new particle(mesh_, is, true));
390  }
391  };
392 
393 
394  //- Destructor
395  virtual ~particle()
396  {}
397 
398 
399  // Member Functions
400 
401  // Access
402 
403  //- Get unique particle creation id
404  inline label getNewParticleID() const;
405 
406  //- Return the mesh database
407  inline const polyMesh& mesh() const;
408 
409  //- Return current particle position
410  inline const vector& position() const;
411 
412  //- Return current particle position
413  inline vector& position();
414 
415  //- Return current cell particle is in
416  inline label& cell();
417 
418  //- Return current cell particle is in
419  inline label cell() const;
420 
421  //- Return current tet face particle is in
422  inline label& tetFace();
423 
424  //- Return current tet face particle is in
425  inline label tetFace() const;
426 
427  //- Return current tet face particle is in
428  inline label& tetPt();
429 
430  //- Return current tet face particle is in
431  inline label tetPt() const;
432 
433  //- Return the indices of the current tet that the
434  // particle occupies.
435  inline tetIndices currentTetIndices() const;
436 
437  //- Return the geometry of the current tet that the
438  // particle occupies.
439  inline tetPointRef currentTet() const;
440 
441  //- Return the normal of the tri on tetFacei_ for the
442  // current tet.
443  inline vector normal() const;
444 
445  //- Return the normal of the tri on tetFacei_ for the
446  // current tet at the start of the timestep, i.e. based
447  // on oldPoints
448  inline vector oldNormal() const;
449 
450  //- Return current face particle is on otherwise -1
451  inline label& face();
452 
453  //- Return current face particle is on otherwise -1
454  inline label face() const;
455 
456  //- Return the impact model to be used, soft or hard (default).
457  inline bool softImpact() const;
458 
459  //- Return the particle current time
460  inline scalar currentTime() const;
461 
462 
463  // Check
464 
465  //- Check the stored cell value (setting if necessary) and
466  // initialise the tetFace and tetPt values
467  inline void initCellFacePt();
468 
469  //- Is the particle on the boundary/(or outside the domain)?
470  inline bool onBoundary() const;
471 
472  //- Is this global face an internal face?
473  inline bool internalFace(const label facei) const;
474 
475  //- Is this global face a boundary face?
476  inline bool boundaryFace(const label facei) const;
477 
478  //- Which patch is particle on
479  inline label patch(const label facei) const;
480 
481  //- Which face of this patch is this particle on
482  inline label patchFace
483  (
484  const label patchi,
485  const label facei
486  ) const;
487 
488  //- Return the fraction of time-step completed
489  inline scalar& stepFraction();
490 
491  //- Return the fraction of time-step completed
492  inline scalar stepFraction() const;
493 
494  //- Return const access to the originating processor id
495  inline label origProc() const;
496 
497  //- Return the originating processor id for manipulation
498  inline label& origProc();
499 
500  //- Return const access to the particle id on originating processor
501  inline label origId() const;
502 
503  //- Return the particle id on originating processor for manipulation
504  inline label& origId();
505 
506 
507  // Track
508 
509  //- Track particle to end of trajectory
510  // or until it hits the boundary.
511  // On entry 'stepFraction()' should be set to the fraction of the
512  // time-step at which the tracking starts and on exit it contains
513  // the fraction of the time-step completed.
514  // Returns the boundary face index if the track stops at the
515  // boundary, -1 otherwise.
516  template<class TrackData>
517  label track(const vector& endPosition, TrackData& td);
518 
519  //- Track particle to a given position and returns 1.0 if the
520  // trajectory is completed without hitting a face otherwise
521  // stops at the face and returns the fraction of the trajectory
522  // completed.
523  // on entry 'stepFraction()' should be set to the fraction of the
524  // time-step at which the tracking starts.
525  template<class TrackData>
526  scalar trackToFace(const vector& endPosition, TrackData& td);
527 
528  //- Return the index of the face to be used in the interpolation
529  // routine
530  inline label faceInterpolation() const;
531 
532 
533  // Transformations
534 
535  //- Transform the physical properties of the particle
536  // according to the given transformation tensor
537  virtual void transformProperties(const tensor& T);
538 
539  //- Transform the physical properties of the particle
540  // according to the given separation vector
541  virtual void transformProperties(const vector& separation);
542 
543  //- The nearest distance to a wall that
544  // the particle can be in the n direction
545  virtual scalar wallImpactDistance(const vector& n) const;
546 
547 
548  // Parallel transfer
549 
550  //- Convert global addressing to the processor patch
551  // local equivalents
552  template<class TrackData>
553  void prepareForParallelTransfer(const label patchi, TrackData& td);
554 
555  //- Convert processor patch addressing to the global equivalents
556  // and set the celli to the face-neighbour
557  template<class TrackData>
558  void correctAfterParallelTransfer(const label patchi, TrackData& td);
559 
560 
561  // I-O
562 
563  //- Read the fields associated with the owner cloud
564  template<class CloudType>
565  static void readFields(CloudType& c);
566 
567  //- Write the fields associated with the owner cloud
568  template<class CloudType>
569  static void writeFields(const CloudType& c);
570 
571  //- Write the particle position and cell
572  void writePosition(Ostream&) const;
573 
574 
575  // Friend Operators
576 
577  friend Ostream& operator<<(Ostream&, const particle&);
578 
579  friend bool operator==(const particle& pA, const particle& pB);
580 
581  friend bool operator!=(const particle& pA, const particle& pB);
582 };
583 
584 
585 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586 
587 } // End namespace Foam
588 
589 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
590 
591 #include "particleI.H"
592 
593 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
594 
595 #ifdef NoRepository
596  #include "particleTemplates.C"
597 #endif
598 
599 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
600 
601 #endif
602 
603 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
bool onBoundary() const
Is the particle on the boundary/(or outside the domain)?
Definition: particleI.H:810
label & tetPt()
Return current tet face particle is in.
Definition: particleI.H:628
A tetrahedron primitive.
Definition: tetrahedron.H:62
void findTris(const vector &position, DynamicList< label > &faceList, const tetPointRef &tet, const FixedList< vector, 4 > &tetAreas, const FixedList< label, 4 > &tetPlaneBasePtIs, const scalar tol) const
Find the tet tri faces between position and tet centre.
Definition: particleI.H:32
scalar currentTime() const
Return the particle current time.
Definition: particleI.H:858
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
uint8_t direction
Definition: direction.H:46
static const scalar trackingCorrectionTol
Fraction of distance to tet centre to move a particle to.
Definition: particle.H:326
scalar movingTetLambda(const vector &from, const vector &to, const label triI, const vector &tetArea, const label tetPlaneBasePtI, const label celli, const label tetFacei, const label tetPtI, const scalar tol) const
Find the lambda value for a moving tri face.
Definition: particleI.H:141
void prepareForParallelTransfer(const label patchi, TrackData &td)
Convert global addressing to the processor patch.
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:566
const double e
Elementary charge.
Definition: doubleFloat.H:78
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
friend bool operator!=(const particle &pA, const particle &pB)
Macros for adding to particle property lists.
void initCellFacePt()
Check the stored cell value (setting if necessary) and.
Definition: particleI.H:670
scalar tetLambda(const vector &from, const vector &to, const label triI, const vector &tetArea, const label tetPlaneBasePtI, const label celli, const label tetFacei, const label tetPtI, const scalar tol) const
Find the lambda value for the line to-from across the.
Definition: particleI.H:69
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &tetIs)
Overridable function to handle the particle hitting a wallPatch.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label tetPtI_
Index of the point on the face that defines the decomposed.
Definition: particle.H:158
void crossEdgeConnectedFace(const label &celli, label &tetFacei, label &tetPtI, const edge &e)
Cross the from the given face across the given edge of the.
Definition: particleI.H:458
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:634
void hitWallFaces(const CloudType &td, const vector &from, const vector &to, scalar &lambdaMin, tetIndices &closestTetIs)
label faceInterpolation() const
Return the index of the face to be used in the interpolation.
Definition: particleI.H:894
label track(const vector &endPosition, TrackData &td)
Track particle to end of trajectory.
static const scalar lambdaDistanceToleranceCoeff
Fraction of the cell volume to use in determining tolerance values.
Definition: particle.H:330
particle(const polyMesh &mesh, const vector &position, const label celli, const label tetFacei, const label tetPtI)
Construct from components.
Definition: particle.C:48
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:123
vector position_
Position of particle.
Definition: particle.H:140
bool hitPatch(const polyPatch &, TrackData &td, const label patchi, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a.
void hitCyclicAMIPatch(const cyclicAMIPolyPatch &, TrackData &td, const vector &direction)
Overridable function to handle the particle hitting a cyclicAMIPatch.
label origId() const
Return const access to the particle id on originating processor.
Definition: particleI.H:840
label origId_
Local particle id on originating processor.
Definition: particle.H:164
Base particle class.
Definition: particle.H:78
bool internalFace(const label facei) const
Is this global face an internal face?
Definition: particleI.H:866
Neighbour processor patch.
friend Ostream & operator<<(Ostream &, const particle &)
tetPointRef currentTet() const
Return the geometry of the current tet that the.
Definition: particleI.H:640
scalar & stepFraction()
Return the fraction of time-step completed.
Definition: particleI.H:816
CloudType & cloud()
Return a reference to the cloud.
Definition: particle.H:125
label origProc_
Originating processor id.
Definition: particle.H:161
const vector & position() const
Return current particle position.
Definition: particleI.H:586
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
label patch(const label facei) const
Which patch is particle on.
Definition: particleI.H:878
TypeName("particle")
Runtime type information.
static void writeFields(const CloudType &c)
Write the fields associated with the owner cloud.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:368
Cyclic plane patch.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
bool softImpact() const
Return the impact model to be used, soft or hard (default).
Definition: particleI.H:852
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
Wedge front and back plane patch.
vector normal() const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:646
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
Cyclic patch for Arbitrary Mesh Interface (AMI)
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:112
void hitWedgePatch(const wedgePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a wedgePatch.
scalar stepFraction_
Fraction of time-step completed.
Definition: particle.H:149
static const scalar minStepFractionTol
Minimum stepFraction tolerance.
Definition: particle.H:333
void hitFace(TrackData &td)
Overridable function to handle the particle hitting a face.
label celli_
Index of the cell it is in.
Definition: particle.H:143
label tetFacei_
Index of the face that owns the decomposed tet that the.
Definition: particle.H:153
label & face()
Return current face particle is on otherwise -1.
Definition: particleI.H:664
label origProc() const
Return const access to the originating processor id.
Definition: particleI.H:828
void tetNeighbour(label triI)
Modify the tet owner data by crossing triI.
Definition: particleI.H:321
Factory class to read-construct particles used for.
Definition: particle.H:375
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Intrusive doubly-linked list.
Definition: IDLList.H:47
scalar trackToFace(const vector &endPosition, TrackData &td)
Track particle to a given position and returns 1.0 if the.
bool boundaryFace(const label facei) const
Is this global face a boundary face?
Definition: particleI.H:872
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:89
label patchi
const dimensionedScalar c
Speed of light in a vacuum.
label & cell()
Return current cell particle is in.
Definition: particleI.H:604
friend bool operator==(const particle &pA, const particle &pB)
label n
static label particleCount_
Cumulative particle counter - used to provode unique ID.
Definition: particle.H:319
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void correctAfterParallelTransfer(const label patchi, TrackData &td)
Convert processor patch addressing to the global equivalents.
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that.
Definition: particle.C:131
label & tetFace()
Return current tet face particle is in.
Definition: particleI.H:616
volScalarField & p
#define DefinePropertyList(str)
bool switchProcessor
Flag to switch processor.
Definition: particle.H:109
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label patchFace(const label patchi, const label facei) const
Which face of this patch is this particle on.
Definition: particleI.H:885
virtual ~particle()
Destructor.
Definition: particle.H:394
TrackingData(CloudType &cloud)
Definition: particle.H:116
bool operator!=(const particle &, const particle &)
Definition: particle.C:145
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
vector oldNormal() const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:652
void hitCyclicPatch(const cyclicPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a cyclicPatch.
Namespace for OpenFOAM.
void hitSymmetryPatch(const symmetryPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
label facei_
Face index if the particle is on a face otherwise -1.
Definition: particle.H:146