59 class cyclicPolyPatch;
60 class cyclicAMIPolyPatch;
61 class cyclicACMIPolyPatch;
62 class cyclicRepeatAMIPolyPatch;
63 class processorPolyPatch;
64 class symmetryPlanePolyPatch;
65 class symmetryPolyPatch;
77 bool operator==(
const particle&,
const particle&);
79 bool operator!=(
const particle&,
const particle&);
92 static const std::size_t sizeofPosition_;
95 static const std::size_t sizeofFields_;
99 static const label maxNBehind_;
118 template <
class TrackCloudType>
150 scalar stepFraction_;
178 inline void stationaryTetGeometry
200 void stationaryTetReverseTransform
210 inline void movingTetGeometry
212 const scalar endStepFraction,
225 const scalar endStepFraction
234 void movingTetReverseTransform
236 const scalar endStepFraction,
261 void changeTet(
const label tetTriI);
264 void changeFace(
const label tetTriI);
274 void changeToMasterPatch();
284 const bool boundaryFail,
285 const string boundaryMsg
295 template<
class TrackCloudType>
299 template<
class TrackCloudType>
304 template<
class TrackCloudType>
308 template<
class TrackCloudType>
312 template<
class TrackCloudType>
316 template<
class TrackCloudType>
319 const vector& displacement,
320 const scalar fraction,
321 TrackCloudType&
cloud,
327 template<
class TrackCloudType>
330 const vector& displacement,
331 const scalar fraction,
332 TrackCloudType& cloud,
338 template<
class TrackCloudType>
341 const vector& displacement,
342 const scalar fraction,
343 TrackCloudType& cloud,
348 template<
class TrackCloudType>
352 template<
class TrackCloudType>
366 "(coordinatesa coordinatesb coordinatesc coordinatesd) " 367 "celli tetFacei tetPti facei stepFraction " 368 "behind nBehind origProc origId" 383 const label tetFacei,
505 inline bool onFace()
const;
535 const vector& displacement,
536 const scalar fraction
542 const vector& displacement,
543 const scalar fraction
549 const vector& displacement,
550 const scalar fraction
558 const vector& displacement,
559 const scalar fraction,
566 const vector& displacement,
567 const scalar fraction,
574 const vector& displacement,
575 const scalar fraction,
582 template<
class TrackCloudType>
585 const vector& displacement,
586 const scalar fraction,
587 TrackCloudType&
cloud,
593 template<
class TrackCloudType>
596 const vector& displacement,
597 const scalar fraction,
598 TrackCloudType& cloud,
603 template<
class TrackCloudType>
606 const vector& displacement,
607 const scalar fraction,
608 TrackCloudType& cloud,
662 const label procCell,
663 const label procTetFace
676 template<
class TrackCloudType>
680 template<
class TrackCloudType>
Template class for intrusive linked lists.
void hitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
label origProc() const
Return the originating processor ID.
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
const polyMesh & mesh() const
Return the mesh database.
trackingData(const TrackCloudType &cloud)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
virtual autoPtr< particle > clone() const
Construct a clone.
A 1D vector of objects of type <T> with a fixed size <Size>.
bool onBoundaryFace() const
Is the particle on a boundary face?
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.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but stops when a tet triangle is hit. On.
scalar trackToAndHitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
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.
const dimensionedScalar & c
Speed of light in a vacuum.
scalar trackToCell(const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label tetPt() const
Return current tet face particle is in.
friend Ostream & operator<<(Ostream &, const particle &)
vector normal() const
Return the normal of the tri on tetFacei_ for the.
void patchData(vector &normal, vector &displacement) const
Get the normal and displacement of the current patch location.
label cell() const
Return current cell particle is in.
void hitFaceNoChangeToMasterPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
As above, but does not change the master patch. Needed in order for.
void prepareForInteractionListReferral(const transformer &transform)
Break the topology and store the particle position so that the.
An ordered pair of two objects of type <T> with first() and second() elements.
TypeName("particle")
Runtime type information.
label face() const
Return current face particle is on otherwise -1.
void reset()
Set step fraction and behind data to zero in preparation for a new.
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.
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
A cloud is a collection of lagrangian particles.
void hitCyclicACMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
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...
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
label origId() const
Return the particle ID on the originating processor.
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
bool onFace() const
Is the particle on a face?
label tetFace() const
Return current tet face particle is in.
label getNewParticleID() const
Get unique particle creation id.
Factory class to read-construct particles used for.
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.
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void writePosition(Ostream &) const
Write the particle position and cell.
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.
const barycentric & coordinates() const
Return current particle coordinates.
static label particleCount_
Cumulative particle counter - used to provide unique ID.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool onInternalFace() const
Is the particle on an internal face?
void hitCyclicAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a cyclicAMIPatch.
Mesh consisting of general polyhedral cells.
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
#define DefinePropertyList(str)
virtual ~particle()
Destructor.
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
bool operator!=(const particle &, const particle &)
vector position() const
Return current particle position.
label patch() const
Return the index of patch that the particle is on.
void hitCyclicRepeatAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
dimensionSet transform(const dimensionSet &)
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.