Go to the documentation of this file.
40 #ifndef wallBoundedParticle_H
41 #define wallBoundedParticle_H
54 class wallBoundedParticle;
56 Ostream&
operator<<(Ostream&,
const wallBoundedParticle&);
57 Ostream&
operator<<(Ostream&,
const InfoProxy<wallBoundedParticle>&);
82 template <
class TrackCloudType>
85 const TrackCloudType&
cloud,
163 const label tetFacei,
175 bool newFormat =
true
233 template<
class TrackCloudType>
236 TrackCloudType&
cloud,
245 template<
class TrackCloudType>
248 TrackCloudType&
cloud,
250 const scalar trackFraction
255 template<
class TrackCloudType>
258 TrackCloudType&
cloud,
263 template<
class TrackCloudType>
266 TrackCloudType&
cloud,
284 template<
class CloudType>
288 template<
class CloudType>
const bitSet & isWallPatch_
static const std::size_t sizeofFields_
Size in bytes of the fields.
A helper class for outputting values to Ostream.
void hitProcessorPatch(TrackCloudType &cloud, trackingData &td)
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
void patchInteraction(TrackCloudType &cloud, trackingData &td, const scalar trackFraction)
Do all patch interaction.
wallBoundedParticle(const polyMesh &c, const point &position, const label celli, const label tetFacei, const label tetPti, const label meshEdgeStart, const label diagEdge)
Construct from components.
static void readFields(CloudType &)
Read.
void hitWallPatch(TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
static void writeFields(const CloudType &)
Write.
Factory class to read-construct particles used for.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
label meshEdgeStart_
Particle is on mesh edge:
void crossDiagonalEdge()
Cross diagonal edge into different triangle on same face,cell.
InfoProxy< wallBoundedParticle > info() const
Return info proxy.
Mesh consisting of general polyhedral cells.
scalar trackToEdge(TrackCloudType &cloud, trackingData &td, const vector &endPosition)
Equivalent of trackToFace.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
label meshEdgeStart() const
-1 or label of mesh edge
void crossEdgeConnectedFace(const label &celli, label &tetFacei, label &tetPti, const edge &e)
Replacement for particle::crossEdgeConnectedFace that avoids bombing.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Particle class that tracks on triangles of boundary faces. Use trackToEdge similar to trackToFace on ...
vector position() const
Return current particle position.
const polyMesh & mesh() const
Return the mesh database.
DSMCCloud< dsmcParcel > CloudType
iNew(const polyMesh &mesh)
trackingData(const TrackCloudType &cloud, const bitSet &isWallPatch)
Templated base class for dsmc cloud.
bool isTriAlongTrack(const vector &n, const point &endPosition) const
Is current triangle in the track direction.
autoPtr< particle > clone() const
Construct and return a clone.
Class used to pass tracking data to the trackToFace function.
point localPosition_
Particle position is updated locally as opposed to via track.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A cloud is a registry collection of lagrangian particles.
const dimensionedScalar e
Elementary charge.
const dimensionedScalar c
Speed of light in a vacuum.
label diagEdge() const
-1 or diagonal edge
edge currentEdge() const
Construct current edge.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
scalar trackFaceTri(const vector &n, const vector &endPosition, label &)
Track through single triangle.
label diagEdge_
Particle is on diagonal edge: