38 template<
class ParticleType>
41 const polyBoundaryMesh& pbm = polyMesh_.boundaryMesh();
45 if (isA<cyclicAMIPolyPatch>(pbm[patchi]))
47 const cyclicAMIPolyPatch& cami =
48 refCast<const cyclicAMIPolyPatch>(pbm[
patchi]);
50 ok = ok && cami.singlePatchProc() != -1;
57 <<
"Particle tracking across AMI patches is only currently " 58 <<
"supported for cases where the AMI patches reside on a " 66 template<
class ParticleType>
70 const word& cloudName,
74 cloud(pMesh, cloudName),
84 polyMesh_.tetBasePtIs();
95 template<
class ParticleType>
102 template<
class ParticleType>
105 delete(this->
remove(&
p));
109 template<
class ParticleType>
114 ParticleType&
p = pIter();
119 <<
"deleting lost particle at position " << p.position()
128 template<
class ParticleType>
133 ParticleType::particleCount_ = 0;
138 template<
class ParticleType>
139 template<
class TrackCloudType>
142 TrackCloudType& cloud,
143 typename ParticleType::trackingData& td,
144 const scalar trackTime
158 const labelList& neighbourProcs = pData[Pstream::myProcNo()];
161 labelList neighbourProcIndices(Pstream::nProcs(), -1);
165 neighbourProcIndices[neighbourProcs[i]] = i;
171 pIter().stepFraction() = 0;
178 neighbourProcs.
size()
185 neighbourProcs.
size()
192 globalPositionsPtr_.clear();
198 forAll(patchIndexTransferLists, i)
200 patchIndexTransferLists[i].clear();
206 ParticleType&
p = pIter();
209 bool keepParticle = p.move(cloud, td, trackTime);
215 if (td.switchProcessor)
221 || !p.onBoundaryFace()
222 || procPatchNeighbours[p.patch()] < 0
226 <<
"Switch processor flag is true when no parallel " 227 <<
"transfer is possible. This is a bug." 232 const label patchi = p.patch();
234 const label n = neighbourProcIndices
236 refCast<const processorPolyPatch>
242 p.prepareForParallelTransfer();
244 particleTransferLists[
n].
append(this->
remove(&p));
246 patchIndexTransferLists[
n].append
248 procPatchNeighbours[patchi]
258 if (!Pstream::parRun())
268 forAll(particleTransferLists, i)
270 if (particleTransferLists[i].size())
279 << patchIndexTransferLists[i]
280 << particleTransferLists[i];
287 pBufs.finishedSends(allNTrans);
290 bool transferred =
false;
310 label neighbProci = neighbourProcs[i];
312 label nRec = allNTrans[neighbProci];
316 UIPstream particleStream(neighbProci, pBufs);
318 labelList receivePatchIndex(particleStream);
323 typename ParticleType::iNew(polyMesh_)
330 ParticleType& newp = newpIter();
332 label patchi = procPatches[receivePatchIndex[pI++]];
334 newp.correctAfterParallelTransfer(patchi, td);
336 addParticle(newParticles.
remove(&newp));
344 template<
class ParticleType>
347 if (!globalPositionsPtr_.valid())
350 <<
"Global positions are not available. " 351 <<
"Cloud::storeGlobalPositions has not been called." 358 polyMesh_.tetBasePtIs();
360 const vectorField& positions = globalPositionsPtr_();
365 iter().autoMap(positions[i], mapper);
371 template<
class ParticleType>
376 this->db().time().
path()/this->
name() +
"_positions.obj" 381 const ParticleType&
p = pIter();
382 pObj<<
"v " << p.position().x() <<
" " << p.position().y() <<
" " 383 << p.position().z() <<
nl;
390 template<
class ParticleType>
398 globalPositionsPtr_.reset(
new vectorField(this->size()));
405 positions[i] = iter().position();
Template class for intrusive linked lists.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
#define forAll(list, i)
Loop across all elements in list.
const labelList & processorPatches() const
Return list of processor patch labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
void size(const label)
Override size to be inconsistent with allocated storage.
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.
volVectorField vectorField(fieldObject, mesh)
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Input inter-processor communications stream operating on external buffer.
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
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.
void append(const T &)
Append an element at the end of the list.
virtual void flush()
Flush stream.
T * remove(T *p)
Remove and return element.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
errorManip< error > abort(error &err)
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
Base cloud calls templated on particle type.
Output inter-processor communications stream operating on external buffer.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
const labelList & processorPatchNeighbours() const
Return processorPatchIndices of the neighbours.
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 move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Mesh consisting of general polyhedral cells.
void writePositions() const
Write positions to <cloudName>_positions.obj file.