39 template<
class ParticleType>
42 const polyBoundaryMesh& pbm = polyMesh_.boundaryMesh();
46 if (isA<cyclicAMIPolyPatch>(pbm[patchi]))
48 const cyclicAMIPolyPatch& cami =
49 refCast<const cyclicAMIPolyPatch>(pbm[
patchi]);
51 ok = ok && cami.singlePatchProc() != -1;
58 <<
"Particle tracking across AMI patches is only currently " 59 <<
"supported for cases where the AMI patches reside on a " 65 template<
class ParticleType>
71 const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
75 if (Pstream::parRun())
79 if (isA<processorPolyPatch>(pbm[patchi]))
81 const processorPolyPatch& ppp =
82 refCast<const processorPolyPatch>(pbm[
patchi]);
84 result[
patchi] = ppp.neighbProcNo();
93 template<
class ParticleType>
99 const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
103 if (Pstream::parRun())
105 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
109 if (isA<processorPolyPatch>(pbm[patchi]))
111 const processorPolyPatch& ppp =
112 refCast<const processorPolyPatch>(pbm[
patchi]);
114 UOPstream(ppp.neighbProcNo(), pBufs)()
119 pBufs.finishedSends();
123 if (isA<processorPolyPatch>(pbm[patchi]))
125 const processorPolyPatch& ppp =
126 refCast<const processorPolyPatch>(pbm[
patchi]);
128 UIPstream(ppp.neighbProcNo(), pBufs)()
138 template<
class ParticleType>
141 const polyMesh& pMesh
144 const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
150 if (isA<nonConformalCyclicPolyPatch>(pbm[patchi]))
152 const nonConformalCyclicPolyPatch& nccPp =
153 refCast<const nonConformalCyclicPolyPatch>(pbm[
patchi]);
155 result[nccPp.origPatchID()].append(patchi);
163 template<
class ParticleType>
166 const polyBoundaryMesh& pbm = polyMesh_.boundaryMesh();
168 forAll(patchNonConformalCyclicPatches_, patchi)
170 forAll(patchNonConformalCyclicPatches_[patchi], i)
172 const label nccPatchi =
173 patchNonConformalCyclicPatches_[
patchi][i];
175 const nonConformalCyclicPolyPatch& nccPp =
176 refCast<const nonConformalCyclicPolyPatch>(pbm[nccPatchi]);
178 if (nccPp.owner()) nccPp.rays();
186 template<
class ParticleType>
190 const word& cloudName,
194 cloud(pMesh, cloudName),
197 patchNbrProc_(patchNbrProc(pMesh)),
198 patchNbrProcPatch_(patchNbrProcPatch(pMesh)),
199 patchNonConformalCyclicPatches_(patchNonConformalCyclicPatches(pMesh)),
200 globalPositionsPtr_()
207 polyMesh_.tetBasePtIs();
208 polyMesh_.oldCellCentres();
210 if (particles.size())
219 template<
class ParticleType>
226 template<
class ParticleType>
229 delete(this->
remove(&
p));
233 template<
class ParticleType>
238 ParticleType&
p = pIter();
243 <<
"deleting lost particle at position " << p.position()
252 template<
class ParticleType>
257 ParticleType::particleCount_ = 0;
262 template<
class ParticleType>
263 template<
class TrackCloudType>
266 TrackCloudType& cloud,
267 typename ParticleType::trackingData& td,
268 const scalar trackTime
272 globalPositionsPtr_.clear();
294 forAll(sendParticles, proci)
296 sendParticles[proci].
clear();
297 sendPatchIndices[proci].
clear();
303 ParticleType&
p = pIter();
306 bool keepParticle = p.move(cloud, td, trackTime);
311 if (td.sendToProc != -1)
314 if (!Pstream::parRun() || !p.onBoundaryFace())
317 <<
"Switch processor flag is true when no parallel " 318 <<
"transfer is possible. This is a bug." 323 p.prepareForParallelTransfer(td);
325 sendParticles[td.sendToProc].
append(this->
remove(&p));
327 sendPatchIndices[td.sendToProc].
append(td.sendToPatch);
337 if (!Pstream::parRun())
346 forAll(sendParticles, proci)
348 if (sendParticles[proci].size())
353 << sendPatchIndices[proci]
354 << sendParticles[proci];
359 labelList receiveSizes(Pstream::nProcs());
363 bool transferred =
false;
364 forAll(receiveSizes, proci)
366 if (receiveSizes[proci])
379 forAll(receiveSizes, proci)
381 if (receiveSizes[proci])
385 const labelList receivePatchIndices(particleStream);
390 typename ParticleType::iNew(polyMesh_)
397 const label patchi = receivePatchIndices[i ++];
399 ParticleType&
p = iter();
403 p.correctAfterParallelTransfer(td);
405 addParticle(newParticles.
remove(&p));
413 template<
class ParticleType>
416 if (!globalPositionsPtr_.valid())
419 <<
"Global positions are not available. " 420 <<
"Cloud::storeGlobalPositions has not been called." 427 polyMesh_.tetBasePtIs();
428 polyMesh_.oldCellCentres();
430 const vectorField& positions = globalPositionsPtr_();
435 iter().autoMap(positions[i], mapper);
441 template<
class ParticleType>
446 this->db().time().path()/this->
name() +
"_positions.obj" 451 const ParticleType&
p = pIter();
452 pObj<<
"v " << p.position().x() <<
" " << p.position().y() <<
" " 453 << p.position().z() <<
nl;
460 template<
class ParticleType>
468 globalPositionsPtr_.reset(
new vectorField(this->size()));
475 positions[i] = iter().position();
Template class for intrusive linked lists.
List< labelList > labelListList
A List of labelList.
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
#define forAll(list, i)
Loop across all elements in list.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
const labelList & patchNbrProcPatch() const
Map from patch index to the corresponding patch index on the.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
const labelList & patchNbrProc() const
Map from patch index to the neighbouring processor index.
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void autoMap(const polyTopoChangeMap &)
Remap the cells of particles corresponding to the.
volVectorField vectorField(fieldObject, mesh)
const dimensionedScalar c
Speed of light in a vacuum.
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
Input inter-processor communications stream operating on external buffer.
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
const labelListList & patchNonConformalCyclicPatches() const
Return map from patch index to connected non-conformal cyclics.
void storeGlobalPositions() const
Call this before a topology change. Stores the particles global.
void clear()
Clear the list, i.e. set size to zero.
A class for handling words, derived from string.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void append(const T &)
Append an element at the end of the list.
virtual void flush()
Flush stream.
List< label > labelList
A List of labels.
T * remove(T *p)
Remove and return element.
errorManip< error > abort(error &err)
Base cloud calls templated on particle type.
Output inter-processor communications stream operating on external buffer.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
void deleteLostParticles()
Remove lost particles from cloud and delete.
void clear()
Clear storage and reset.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
#define WarningInFunction
Report a warning using Foam::Warning.
Mesh consisting of general polyhedral cells.
void writePositions() const
Write positions to <cloudName>_positions.obj file.