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-2015 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
106 
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
316  DefinePropertyList("(Px Py Pz) cellI tetFaceI tetPtI origProc origId");
317 
318  //- Cumulative particle counter - used to provode unique ID
319  static label particleCount_;
320 
321  //- Fraction of distance to tet centre to move a particle to
322  // 'rescue' it from a tracking problem
323  static const scalar trackingCorrectionTol;
324 
325  //- Fraction of the cell volume to use in determining tolerance values
326  // for the denominator and numerator of lambda
327  static const scalar lambdaDistanceToleranceCoeff;
328 
329  //- Minimum stepFraction tolerance
330  static const scalar minStepFractionTol;
331 
332 
333  // Constructors
334 
335  //- Construct from components
336  particle
337  (
338  const polyMesh& mesh,
339  const vector& position,
340  const label cellI,
341  const label tetFaceI,
342  const label tetPtI
343  );
344 
345  //- Construct from components, tetFaceI_ and tetPtI_ are not
346  // supplied so they will be deduced by a search
347  particle
348  (
349  const polyMesh& mesh,
350  const vector& position,
351  const label cellI,
352  bool doCellFacePt = true
353  );
354 
355  //- Construct from Istream
356  particle(const polyMesh& mesh, Istream&, bool readFields = true);
357 
358  //- Construct as a copy
359  particle(const particle& p);
360 
361  //- Construct as a copy with refernce to a new mesh
362  particle(const particle& p, const polyMesh& mesh);
363 
364  //- Construct a clone
365  virtual autoPtr<particle> clone() const
366  {
367  return autoPtr<particle>(new particle(*this));
368  }
369 
370  //- Factory class to read-construct particles used for
371  // parallel transfer
372  class iNew
373  {
374  const polyMesh& mesh_;
375 
376  public:
377 
378  iNew(const polyMesh& mesh)
379  :
380  mesh_(mesh)
381  {}
382 
384  {
385  return autoPtr<particle>(new particle(mesh_, is, true));
386  }
387  };
388 
389 
390  //- Destructor
391  virtual ~particle()
392  {}
393 
394 
395  // Member Functions
396 
397  // Access
398 
399  //- Get unique particle creation id
400  inline label getNewParticleID() const;
401 
402  //- Return the mesh database
403  inline const polyMesh& mesh() const;
404 
405  //- Return current particle position
406  inline const vector& position() const;
407 
408  //- Return current particle position
409  inline vector& position();
410 
411  //- Return current cell particle is in
412  inline label& cell();
413 
414  //- Return current cell particle is in
415  inline label cell() const;
416 
417  //- Return current tet face particle is in
418  inline label& tetFace();
419 
420  //- Return current tet face particle is in
421  inline label tetFace() const;
422 
423  //- Return current tet face particle is in
424  inline label& tetPt();
425 
426  //- Return current tet face particle is in
427  inline label tetPt() const;
428 
429  //- Return the indices of the current tet that the
430  // particle occupies.
431  inline tetIndices currentTetIndices() const;
432 
433  //- Return the geometry of the current tet that the
434  // particle occupies.
435  inline tetPointRef currentTet() const;
436 
437  //- Return the normal of the tri on tetFaceI_ for the
438  // current tet.
439  inline vector normal() const;
440 
441  //- Return the normal of the tri on tetFaceI_ for the
442  // current tet at the start of the timestep, i.e. based
443  // on oldPoints
444  inline vector oldNormal() const;
445 
446  //- Return current face particle is on otherwise -1
447  inline label& face();
448 
449  //- Return current face particle is on otherwise -1
450  inline label face() const;
451 
452  //- Return the impact model to be used, soft or hard (default).
453  inline bool softImpact() const;
454 
455  //- Return the particle current time
456  inline scalar currentTime() const;
457 
458 
459  // Check
460 
461  //- Check the stored cell value (setting if necessary) and
462  // initialise the tetFace and tetPt values
463  inline void initCellFacePt();
464 
465  //- Is the particle on the boundary/(or outside the domain)?
466  inline bool onBoundary() const;
467 
468  //- Is this global face an internal face?
469  inline bool internalFace(const label faceI) const;
470 
471  //- Is this global face a boundary face?
472  inline bool boundaryFace(const label faceI) const;
473 
474  //- Which patch is particle on
475  inline label patch(const label faceI) const;
476 
477  //- Which face of this patch is this particle on
478  inline label patchFace
479  (
480  const label patchI,
481  const label faceI
482  ) const;
483 
484  //- Return the fraction of time-step completed
485  inline scalar& stepFraction();
486 
487  //- Return the fraction of time-step completed
488  inline scalar stepFraction() const;
489 
490  //- Return const access to the originating processor id
491  inline label origProc() const;
492 
493  //- Return the originating processor id for manipulation
494  inline label& origProc();
495 
496  //- Return const access to the particle id on originating processor
497  inline label origId() const;
498 
499  //- Return the particle id on originating processor for manipulation
500  inline label& origId();
501 
502 
503  // Track
504 
505  //- Track particle to end of trajectory
506  // or until it hits the boundary.
507  // On entry 'stepFraction()' should be set to the fraction of the
508  // time-step at which the tracking starts and on exit it contains
509  // the fraction of the time-step completed.
510  // Returns the boundary face index if the track stops at the
511  // boundary, -1 otherwise.
512  template<class TrackData>
513  label track(const vector& endPosition, TrackData& td);
514 
515  //- Track particle to a given position and returns 1.0 if the
516  // trajectory is completed without hitting a face otherwise
517  // stops at the face and returns the fraction of the trajectory
518  // completed.
519  // on entry 'stepFraction()' should be set to the fraction of the
520  // time-step at which the tracking starts.
521  template<class TrackData>
522  scalar trackToFace(const vector& endPosition, TrackData& td);
523 
524  //- Return the index of the face to be used in the interpolation
525  // routine
526  inline label faceInterpolation() const;
527 
528 
529  // Transformations
530 
531  //- Transform the physical properties of the particle
532  // according to the given transformation tensor
533  virtual void transformProperties(const tensor& T);
534 
535  //- Transform the physical properties of the particle
536  // according to the given separation vector
537  virtual void transformProperties(const vector& separation);
538 
539  //- The nearest distance to a wall that
540  // the particle can be in the n direction
541  virtual scalar wallImpactDistance(const vector& n) const;
542 
543 
544  // Parallel transfer
545 
546  //- Convert global addressing to the processor patch
547  // local equivalents
548  template<class TrackData>
549  void prepareForParallelTransfer(const label patchI, TrackData& td);
550 
551  //- Convert processor patch addressing to the global equivalents
552  // and set the cellI to the face-neighbour
553  template<class TrackData>
554  void correctAfterParallelTransfer(const label patchI, TrackData& td);
555 
556 
557  // I-O
558 
559  //- Read the fields associated with the owner cloud
560  template<class CloudType>
561  static void readFields(CloudType& c);
562 
563  //- Write the fields associated with the owner cloud
564  template<class CloudType>
565  static void writeFields(const CloudType& c);
566 
567  //- Write the particle position and cell
568  void writePosition(Ostream&) const;
569 
570 
571  // Friend Operators
572 
573  friend Ostream& operator<<(Ostream&, const particle&);
574 
575  friend bool operator==(const particle& pA, const particle& pB);
576 
577  friend bool operator!=(const particle& pA, const particle& pB);
578 };
579 
580 
581 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
582 
583 } // End namespace Foam
584 
585 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586 
587 #include "particleI.H"
588 
589 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
590 
591 #ifdef NoRepository
592 # include "particleTemplates.C"
593 #endif
594 
595 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
596 
597 #endif
598 
599 // ************************************************************************* //
Foam::particle::boundaryFace
bool boundaryFace(const label faceI) const
Is this global face a boundary face?
Definition: particleI.H:872
Foam::Tensor
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
Definition: complexI.H:224
Foam::particle::hitPatch
bool hitPatch(const polyPatch &, TrackData &td, const label patchI, const scalar trackFraction, const tetIndices &tetIs)
Overridable function to handle the particle hitting a.
Definition: particleTemplates.C:932
Foam::particle::stepFraction_
scalar stepFraction_
Fraction of time-step completed.
Definition: particle.H:149
Foam::particle::onBoundary
bool onBoundary() const
Is the particle on the boundary/(or outside the domain)?
Definition: particleI.H:810
Foam::particle::face
label & face()
Return current face particle is on otherwise -1.
Definition: particleI.H:658
p
p
Definition: pEqn.H:62
Foam::particle::origId
label origId() const
Return const access to the particle id on originating processor.
Definition: particleI.H:840
Foam::particle::hitWedgePatch
void hitWedgePatch(const wedgePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a wedgePatch.
Definition: particleTemplates.C:946
Foam::wedgePolyPatch
Wedge front and back plane patch.
Definition: wedgePolyPatch.H:48
Foam::particle::minStepFractionTol
static const scalar minStepFractionTol
Minimum stepFraction tolerance.
Definition: particle.H:329
Foam::particle::iNew
Factory class to read-construct particles used for.
Definition: particle.H:371
Foam::particle::hitFace
void hitFace(TrackData &td)
Overridable function to handle the particle hitting a face.
Definition: particleTemplates.C:926
Foam::particle::wallImpactDistance
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that.
Definition: particle.C:131
Foam::cyclicPolyPatch
Cyclic plane patch.
Definition: cyclicPolyPatch.H:63
Foam::particle::trackToFace
scalar trackToFace(const vector &endPosition, TrackData &td)
Track particle to a given position and returns 1.0 if the.
Foam::particle::crossEdgeConnectedFace
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
Foam::particle::stepFraction
scalar & stepFraction()
Return the fraction of time-step completed.
Definition: particleI.H:816
Foam::DynamicList< label >
Foam::particle::TrackingData::keepParticle
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:112
particleTemplates.C
Foam::particle::TrackingData::cloud
CloudType & cloud()
Return a reference to the cloud.
Definition: particle.H:125
Foam::particle::getNewParticleID
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:566
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::particle::softImpact
bool softImpact() const
Return the impact model to be used, soft or hard (default).
Definition: particleI.H:852
Cloud.H
Foam::particle::tetNeighbour
void tetNeighbour(label triI)
Modify the tet owner data by crossing triI.
Definition: particleI.H:321
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::particle::tetLambda
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
Foam::particle::currentTetIndices
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:634
Foam::particle::internalFace
bool internalFace(const label faceI) const
Is this global face an internal face?
Definition: particleI.H:866
Foam::particle::faceInterpolation
label faceInterpolation() const
Return the index of the face to be used in the interpolation.
Definition: particleI.H:894
Foam::particle::operator==
friend bool operator==(const particle &pA, const particle &pB)
Foam::particle::faceI_
label faceI_
Face index if the particle is on a face otherwise -1.
Definition: particle.H:146
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
faceList.H
Foam::particle::particleCount_
static label particleCount_
Cumulative particle counter - used to provode unique ID.
Definition: particle.H:315
Foam::particle::writePosition
void writePosition(Ostream &) const
Write the particle position and cell.
Definition: particleIO.C:89
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
polyMeshTetDecomposition.H
OFstream.H
Foam::particle::~particle
virtual ~particle()
Destructor.
Definition: particle.H:390
Foam::particle::cell
label & cell()
Return current cell particle is in.
Definition: particleI.H:598
n
label n
Definition: TABSMDCalcMethod2.H:31
particleMacros.H
Macros for adding to particle property lists.
Foam::particle::tetFace
label & tetFace()
Return current tet face particle is in.
Definition: particleI.H:610
Foam::particle::hitSymmetryPatch
void hitSymmetryPatch(const symmetryPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Definition: particleTemplates.C:978
Foam::label
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
Foam::IDLList
Intrusive doubly-linked list.
Definition: IDLList.H:47
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::particle::tetPt
label & tetPt()
Return current tet face particle is in.
Definition: particleI.H:622
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::particle::hitSymmetryPlanePatch
void hitSymmetryPlanePatch(const symmetryPlanePolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Definition: particleTemplates.C:964
Foam::symmetryPlanePolyPatch
Symmetry-plane patch.
Definition: symmetryPlanePolyPatch.H:48
Foam::particle::hitCyclicPatch
void hitCyclicPatch(const cyclicPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a cyclicPatch.
Definition: particleTemplates.C:992
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::symmetryPolyPatch
Symmetry patch for non-planar or multi-plane patches.
Definition: symmetryPolyPatch.H:51
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
Foam::particle::mesh
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:580
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
Foam::particle::clone
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:364
Foam::particle::iNew::operator()
autoPtr< particle > operator()(Istream &is) const
Definition: particle.H:382
Foam::operator!=
bool operator!=(const particle &, const particle &)
Definition: particle.C:147
Foam::particle::movingTetLambda
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
Foam::wallPolyPatch
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Foam::particle::trackingCorrectionTol
static const scalar trackingCorrectionTol
Fraction of distance to tet centre to move a particle to.
Definition: particle.H:322
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::particle::correctAfterParallelTransfer
void correctAfterParallelTransfer(const label patchI, TrackData &td)
Convert processor patch addressing to the global equivalents.
Definition: particleTemplates.C:53
Foam::particle::particle
particle(const polyMesh &mesh, const vector &position, const label cellI, const label tetFaceI, const label tetPtI)
Construct from components.
Definition: particle.C:48
Foam::particle::tetPtI_
label tetPtI_
Index of the point on the face that defines the decomposed.
Definition: particle.H:158
Foam::particle::findTris
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
particleI.H
Foam::particle::operator!=
friend bool operator!=(const particle &pA, const particle &pB)
Foam::particle::cellI_
label cellI_
Index of the cell it is in.
Definition: particle.H:143
pointField.H
Foam::particle::lambdaDistanceToleranceCoeff
static const scalar lambdaDistanceToleranceCoeff
Fraction of the cell volume to use in determining tolerance values.
Definition: particle.H:326
Foam::particle::TrackingData::switchProcessor
bool switchProcessor
Flag to switch processor.
Definition: particle.H:109
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:73
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::cloud
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::particle::hitWallPatch
void hitWallPatch(const wallPolyPatch &, TrackData &td, const tetIndices &tetIs)
Overridable function to handle the particle hitting a wallPatch.
Definition: particleTemplates.C:1105
Foam::particle::initCellFacePt
void initCellFacePt()
Check the stored cell value (setting if necessary) and.
Definition: particleI.H:670
DefinePropertyList
#define DefinePropertyList(str)
Definition: particleMacros.H:42
Foam::particle::hitCyclicAMIPatch
void hitCyclicAMIPatch(const cyclicAMIPolyPatch &, TrackData &td, const vector &direction)
Overridable function to handle the particle hitting a cyclicAMIPatch.
Definition: particleTemplates.C:1040
Foam::particle::TrackingData::cloud_
CloudType & cloud_
Reference to the cloud containing (this) particle.
Definition: particle.H:99
Foam::particle::transformProperties
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:123
Foam::Vector< scalar >
Foam::particle::position_
vector position_
Position of particle.
Definition: particle.H:140
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::particle::TrackingData::cloudType
CloudType cloudType
Definition: particle.H:106
Foam::particle::origProc
label origProc() const
Return const access to the originating processor id.
Definition: particleI.H:828
Foam::particle::hitWallFaces
void hitWallFaces(const CloudType &td, const vector &from, const vector &to, scalar &lambdaMin, tetIndices &closestTetIs)
Definition: particleTemplates.C:717
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::particle
Base particle class.
Definition: particle.H:78
Foam::particle::track
label track(const vector &endPosition, TrackData &td)
Track particle to end of trajectory.
Foam::particle::operator<<
friend Ostream & operator<<(Ostream &, const particle &)
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::particle::position
const vector & position() const
Return current particle position.
Definition: particleI.H:586
vector.H
Foam::particle::oldNormal
vector oldNormal() const
Return the normal of the tri on tetFaceI_ for the.
Definition: particleI.H:652
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::particle::readFields
static void readFields(CloudType &c)
Read the fields associated with the owner cloud.
Definition: particleTemplates.C:129
Foam::particle::hitProcessorPatch
void hitProcessorPatch(const processorPolyPatch &, TrackData &td)
Overridable function to handle the particle hitting a.
Definition: particleTemplates.C:1099
Foam::particle::tetFaceI_
label tetFaceI_
Index of the face that owns the decomposed tet that the.
Definition: particle.H:153
Foam::particle::normal
vector normal() const
Return the normal of the tri on tetFaceI_ for the.
Definition: particleI.H:646
Foam::particle::currentTime
scalar currentTime() const
Return the particle current time.
Definition: particleI.H:858
Foam::particle::sizeofFields_
static const std::size_t sizeofFields_
Size in bytes of the fields.
Definition: particle.H:88
Foam::particle::iNew::mesh_
const polyMesh & mesh_
Definition: particle.H:373
Foam::particle::origProc_
label origProc_
Originating processor id.
Definition: particle.H:161
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
FixedList.H
Foam::particle::TypeName
TypeName("particle")
Runtime type information.
Foam::particle::iNew::iNew
iNew(const polyMesh &mesh)
Definition: particle.H:377
Foam::particle::prepareForParallelTransfer
void prepareForParallelTransfer(const label patchI, TrackData &td)
Convert global addressing to the processor patch.
Definition: particleTemplates.C:41
Foam::particle::writeFields
static void writeFields(const CloudType &c)
Write the fields associated with the owner cloud.
Definition: particleTemplates.C:159
Foam::particle::patchFace
label patchFace(const label patchI, const label faceI) const
Which face of this patch is this particle on.
Definition: particleI.H:885
Foam::particle::origId_
label origId_
Local particle id on originating processor.
Definition: particle.H:164
Foam::particle::patch
label patch(const label faceI) const
Which patch is particle on.
Definition: particleI.H:878
Foam::particle::TrackingData
Definition: particle.H:94
IDLList.H
Foam::particle::TrackingData::TrackingData
TrackingData(CloudType &cloud)
Definition: particle.H:116
Foam::particle::mesh_
const polyMesh & mesh_
Reference to the polyMesh database.
Definition: particle.H:137
tetrahedron.H
Foam::tetrahedron
A tetrahedron primitive.
Definition: tetrahedron.H:62
Foam::particle::currentTet
tetPointRef currentTet() const
Return the geometry of the current tet that the.
Definition: particleI.H:640
Foam::particle::sizeofPosition_
static const std::size_t sizeofPosition_
Size in bytes of the position data.
Definition: particle.H:85
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:51