Cloud.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017, 2020 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "Cloud.H"
30 #include "processorPolyPatch.H"
31 #include "globalMeshData.H"
33 #include "mapPolyMesh.H"
34 #include "Time.H"
35 #include "OFstream.H"
36 #include "wallPolyPatch.H"
37 #include "cyclicAMIPolyPatch.H"
38 
39 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40 
41 template<class ParticleType>
43 {
44  const polyBoundaryMesh& pbm = polyMesh_.boundaryMesh();
45  bool ok = true;
46  for (const polyPatch& pp : pbm)
47  {
48  if (isA<cyclicAMIPolyPatch>(pp))
49  {
50  const cyclicAMIPolyPatch& cami =
51  refCast<const cyclicAMIPolyPatch>(pp);
52 
53  if (cami.owner())
54  {
55  ok = ok && (cami.AMI().singlePatchProc() != -1);
56  }
57  }
58  }
59 
60  if (!ok)
61  {
63  << "Particle tracking across AMI patches is only currently "
64  << "supported for cases where the AMI patches reside on a "
65  << "single processor" << abort(FatalError);
66  }
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
72 template<class ParticleType>
74 (
75  const polyMesh& pMesh,
76  const word& cloudName,
77  const IDLList<ParticleType>& particles
78 )
79 :
80  cloud(pMesh, cloudName),
82  polyMesh_(pMesh),
83  labels_(),
84  globalPositionsPtr_(),
85  geometryType_(cloud::geometryType::COORDINATES)
86 {
87  checkPatches();
88  polyMesh_.oldCellCentres();
89 
90  // Ask for the tetBasePtIs to trigger all processors to build
91  // them, otherwise, if some processors have no particles then
92  // there is a comms mismatch.
93  polyMesh_.tetBasePtIs();
94 
95  if (particles.size())
96  {
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
104 template<class ParticleType>
106 {
107  this->append(pPtr);
108 }
109 
110 
111 template<class ParticleType>
113 {
114  delete(this->remove(&p));
115 }
116 
117 
118 template<class ParticleType>
120 {
121  for (ParticleType& p : *this)
122  {
123  if (p.cell() == -1)
124  {
126  << "deleting lost particle at position " << p.position()
127  << endl;
128 
129  deleteParticle(p);
130  }
131  }
132 }
133 
134 
135 template<class ParticleType>
137 {
138  // Reset particle count and particles only
139  // - not changing the cloud object registry or reference to the polyMesh
140  ParticleType::particleCount_ = 0;
142 }
143 
144 
145 template<class ParticleType>
146 template<class TrackCloudType>
148 (
149  TrackCloudType& cloud,
150  typename ParticleType::trackingData& td,
151  const scalar trackTime
152 )
153 {
154  const polyBoundaryMesh& pbm = pMesh().boundaryMesh();
155  const globalMeshData& pData = polyMesh_.globalData();
156 
157  // Which patches are processor patches
158  const labelList& procPatches = pData.processorPatches();
159 
160  // Indexing of equivalent patch on neighbour processor into the
161  // procPatches list on the neighbour
162  const labelList& procPatchNeighbours = pData.processorPatchNeighbours();
163 
164  // Which processors this processor is connected to
165  const labelList& neighbourProcs = pData[Pstream::myProcNo()];
166 
167  // Indexing from the processor number into the neighbourProcs list
168  labelList neighbourProcIndices(Pstream::nProcs(), -1);
169 
170  forAll(neighbourProcs, i)
171  {
172  neighbourProcIndices[neighbourProcs[i]] = i;
173  }
174 
175  // Initialise the stepFraction moved for the particles
176  forAllIters(*this, pIter)
177  {
178  pIter().reset();
179  }
180 
181  // List of lists of particles to be transferred for all of the
182  // neighbour processors
183  List<IDLList<ParticleType>> particleTransferLists
184  (
185  neighbourProcs.size()
186  );
187 
188  // List of destination processorPatches indices for all of the
189  // neighbour processors
190  List<DynamicList<label>> patchIndexTransferLists
191  (
192  neighbourProcs.size()
193  );
194 
195  // Allocate transfer buffers
196  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
197 
198  // Clear the global positions as there are about to change
199  globalPositionsPtr_.clear();
200 
201  // While there are particles to transfer
202  while (true)
203  {
204  particleTransferLists = IDLList<ParticleType>();
205  forAll(patchIndexTransferLists, i)
206  {
207  patchIndexTransferLists[i].clear();
208  }
209 
210  // Loop over all particles
211  for (ParticleType& p : *this)
212  {
213  // Move the particle
214  bool keepParticle = p.move(cloud, td, trackTime);
215 
216  // If the particle is to be kept
217  // (i.e. it hasn't passed through an inlet or outlet)
218  if (keepParticle)
219  {
220  if (td.switchProcessor)
221  {
222  #ifdef FULLDEBUG
223  if
224  (
225  !Pstream::parRun()
226  || !p.onBoundaryFace()
227  || procPatchNeighbours[p.patch()] < 0
228  )
229  {
231  << "Switch processor flag is true when no parallel "
232  << "transfer is possible. This is a bug."
233  << exit(FatalError);
234  }
235  #endif
236 
237  const label patchi = p.patch();
238 
239  const label n = neighbourProcIndices
240  [
241  refCast<const processorPolyPatch>
242  (
243  pbm[patchi]
244  ).neighbProcNo()
245  ];
246 
247  p.prepareForParallelTransfer();
248 
249  particleTransferLists[n].append(this->remove(&p));
250 
251  patchIndexTransferLists[n].append
252  (
253  procPatchNeighbours[patchi]
254  );
255  }
256  }
257  else
258  {
259  deleteParticle(p);
260  }
261  }
262 
263  if (!Pstream::parRun())
264  {
265  break;
266  }
267 
268 
269  // Clear transfer buffers
270  pBufs.clear();
271 
272  // Stream into send buffers
273  forAll(particleTransferLists, i)
274  {
275  if (particleTransferLists[i].size())
276  {
277  UOPstream particleStream
278  (
279  neighbourProcs[i],
280  pBufs
281  );
282 
283  particleStream
284  << patchIndexTransferLists[i]
285  << particleTransferLists[i];
286  }
287  }
288 
289 
290  // Start sending. Sets number of bytes transferred
291  labelList allNTrans(Pstream::nProcs());
292  pBufs.finishedSends(allNTrans);
293 
294 
295  bool transferred = false;
296 
297  for (const label n : allNTrans)
298  {
299  if (n)
300  {
301  transferred = true;
302  break;
303  }
304  }
305  reduce(transferred, orOp<bool>());
306 
307  if (!transferred)
308  {
309  break;
310  }
311 
312  // Retrieve from receive buffers
313  for (const label neighbProci : neighbourProcs)
314  {
315  label nRec = allNTrans[neighbProci];
316 
317  if (nRec)
318  {
319  UIPstream particleStream(neighbProci, pBufs);
320 
321  labelList receivePatchIndex(particleStream);
322 
323  IDLList<ParticleType> newParticles
324  (
325  particleStream,
326  typename ParticleType::iNew(polyMesh_)
327  );
328 
329  label pI = 0;
330 
331  for (ParticleType& newp : newParticles)
332  {
333  label patchi = procPatches[receivePatchIndex[pI++]];
334 
335  newp.correctAfterParallelTransfer(patchi, td);
336 
337  addParticle(newParticles.remove(&newp));
338  }
339  }
340  }
341  }
342 }
343 
344 
345 template<class ParticleType>
347 {
348  if (!globalPositionsPtr_)
349  {
351  << "Global positions are not available. "
352  << "Cloud::storeGlobalPositions has not been called."
353  << exit(FatalError);
354  }
355 
356  // Reset stored data that relies on the mesh
357  // polyMesh_.clearCellTree();
358  cellWallFacesPtr_.clear();
359 
360  // Ask for the tetBasePtIs to trigger all processors to build
361  // them, otherwise, if some processors have no particles then
362  // there is a comms mismatch.
363  polyMesh_.tetBasePtIs();
364  polyMesh_.oldCellCentres();
365 
366  const vectorField& positions = globalPositionsPtr_();
367 
368  label i = 0;
369  forAllIters(*this, iter)
370  {
371  iter().autoMap(positions[i], mapper);
372  ++i;
373  }
374 }
375 
376 
377 template<class ParticleType>
379 {
380  OFstream pObj
381  (
382  this->db().time().path()/this->name() + "_positions.obj"
383  );
384 
385  forAllConstIters(*this, pIter)
386  {
387  const ParticleType& p = pIter();
388  const point position(p.position());
389  pObj<< "v " << position.x() << " " << position.y() << " "
390  << position.z() << nl;
391  }
392 
393  pObj.flush();
394 }
395 
396 
397 template<class ParticleType>
399 {
400  // Store the global positions for later use by autoMap. It would be
401  // preferable not to need this. If the mapPolyMesh object passed to autoMap
402  // had a copy of the old mesh then the global positions could be recovered
403  // within autoMap, and this pre-processing would not be necessary.
404 
405  globalPositionsPtr_.reset(new vectorField(this->size()));
406 
407  vectorField& positions = globalPositionsPtr_();
408 
409  label i = 0;
410  forAllConstIters(*this, iter)
411  {
412  positions[i] = iter().position();
413  ++i;
414  }
415 }
416 
417 
418 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
419 
420 #include "CloudIO.C"
421 
422 // ************************************************************************* //
Foam::Cloud::addParticle
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:105
p
volScalarField & p
Definition: createFieldRefs.H:8
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::UOPstream
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:57
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const
Return list of processor patch labels.
Definition: globalMeshData.H:381
globalMeshData.H
Foam::Cloud::writePositions
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:378
mapPolyMesh.H
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:87
wallPolyPatch.H
Cloud.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
CloudIO.C
append
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
PstreamCombineReduceOps.H
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::Cloud::storeGlobalPositions
void storeGlobalPositions() const
Call this before a topology change.
Definition: Cloud.C:398
Foam::Cloud::deleteParticle
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition: Cloud.C:112
Foam::Cloud::autoMap
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:346
Foam::Cloud::deleteLostParticles
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:119
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::ILList
Template class for intrusive linked lists.
Definition: ILList.H:52
Foam::Cloud::move
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:148
Foam::FatalError
error FatalError
processorPolyPatch.H
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:107
reduce
reduce(hasMovingMesh, orOp< bool >())
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::globalMeshData::processorPatchNeighbours
const labelList & processorPatchNeighbours() const
Return processorPatchIndices of the neighbours.
Definition: globalMeshData.H:396
cyclicAMIPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::Cloud::cloudReset
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition: Cloud.C:136
Time.H
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:57
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::Vector< scalar >
Foam::List< label >
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:54
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::Cloud::Cloud
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
Definition: Cloud.C:74
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::UIPstream
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:56
Foam::orOp
Definition: ops.H:234
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::OSstream::flush
virtual void flush()
Flush stream.
Definition: OSstream.C:269