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_;
105 static const scalar negativeSpaceDisplacementFactor;
124 template <
class TrackCloudType>
156 scalar stepFraction_;
172 inline void stationaryTetGeometry
194 void stationaryTetReverseTransform
204 inline void movingTetGeometry
206 const scalar endStepFraction,
219 const scalar endStepFraction
228 void movingTetReverseTransform
230 const scalar endStepFraction,
255 void changeTet(
const label tetTriI);
258 void changeFace(
const label tetTriI);
268 void changeToMasterPatch();
279 const bool boundaryFail,
280 const string boundaryMsg
290 template<
class TrackCloudType>
294 template<
class TrackCloudType>
299 template<
class TrackCloudType>
303 template<
class TrackCloudType>
307 template<
class TrackCloudType>
311 template<
class TrackCloudType>
316 template<
class TrackCloudType>
321 template<
class TrackCloudType>
330 template<
class TrackCloudType>
334 template<
class TrackCloudType>
348 "(coordinatesa coordinatesb coordinatesc coordinatesd) " 349 "celli tetFacei tetPti facei stepFraction origProc origId" 364 const label tetFacei,
486 inline bool onFace()
const;
512 const vector& displacement,
513 const scalar fraction
519 const vector& displacement,
520 const scalar fraction
526 const vector& displacement,
527 const scalar fraction
535 const vector& displacement,
536 const scalar fraction,
543 const vector& displacement,
544 const scalar fraction,
551 const vector& displacement,
552 const scalar fraction,
559 template<
class TrackCloudType>
563 TrackCloudType&
cloud,
568 template<
class TrackCloudType>
572 const scalar fraction,
573 TrackCloudType& cloud,
631 const label procCell,
632 const label procTetFace
645 template<
class TrackCloudType>
649 template<
class TrackCloudType>
Template class for intrusive linked lists.
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
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.
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.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
scalar trackToCell(const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
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.
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.
label cell() const
Return current cell particle is in.
An ordered pair of two objects of type <T> with first() and second() elements.
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
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.
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 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.
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.
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.
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.
const barycentric & coordinates() const
Return current particle coordinates.
static label particleCount_
Cumulative particle counter - used to provode 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?
Mesh consisting of general polyhedral cells.
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.
#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.
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.
dimensionSet transform(const dimensionSet &)
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.