Cloud.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Cloud.H"
27 #include "processorPolyPatch.H"
28 #include "globalMeshData.H"
29 #include "meshToMesh.H"
31 #include "polyTopoChangeMap.H"
32 #include "polyMeshMap.H"
33 #include "polyDistributionMap.H"
34 #include "Time.H"
35 #include "OFstream.H"
36 #include "wallPolyPatch.H"
38 #include "cpuLoad.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 template<class Type>
47 {
48  void operator()(IDLList<Type>& x, const IDLList<Type>& y) const
49  {
50  if (y.size())
51  {
52  forAllConstIter(typename IDLList<Type>, y, iter)
53  {
54  x.append(new Type(iter()));
55  }
56  }
57  }
58 };
59 
60 }
61 
62 
63 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
64 
65 template<class ParticleType>
67 (
68  const polyMesh& pMesh
69 )
70 {
71  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
72 
73  labelList result(pbm.size(), -1);
74 
75  if (Pstream::parRun())
76  {
77  forAll(pbm, patchi)
78  {
79  if (isA<processorPolyPatch>(pbm[patchi]))
80  {
81  const processorPolyPatch& ppp =
82  refCast<const processorPolyPatch>(pbm[patchi]);
83 
84  result[patchi] = ppp.neighbProcNo();
85  }
86  }
87  }
88 
89  return result;
90 }
91 
92 
93 template<class ParticleType>
95 (
96  const polyMesh& pMesh
97 )
98 {
99  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
100 
101  labelList result(pbm.size(), -1);
102 
103  if (Pstream::parRun())
104  {
105  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
106 
107  forAll(pbm, patchi)
108  {
109  if (isA<processorPolyPatch>(pbm[patchi]))
110  {
111  const processorPolyPatch& ppp =
112  refCast<const processorPolyPatch>(pbm[patchi]);
113 
114  UOPstream(ppp.neighbProcNo(), pBufs)()
115  << ppp.index();
116  }
117  }
118 
119  pBufs.finishedSends();
120 
121  forAll(pbm, patchi)
122  {
123  if (isA<processorPolyPatch>(pbm[patchi]))
124  {
125  const processorPolyPatch& ppp =
126  refCast<const processorPolyPatch>(pbm[patchi]);
127 
128  UIPstream(ppp.neighbProcNo(), pBufs)()
129  >> result[patchi];
130  }
131  }
132  }
133 
134  return result;
135 }
136 
137 
138 template<class ParticleType>
141 (
142  const polyMesh& pMesh
143 )
144 {
145  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
146 
147  labelListList result(pbm.size(), labelList());
148 
149  forAll(pbm, patchi)
150  {
151  if (isA<nonConformalCyclicPolyPatch>(pbm[patchi]))
152  {
153  const nonConformalCyclicPolyPatch& nccPp =
154  refCast<const nonConformalCyclicPolyPatch>(pbm[patchi]);
155 
156  result[nccPp.origPatchIndex()].append(patchi);
157  }
158  }
159 
160  return result;
161 }
162 
163 
164 template<class ParticleType>
166 {
167  const polyBoundaryMesh& pbm = pMesh_.boundaryMesh();
168 
169  forAll(patchNonConformalCyclicPatches_, patchi)
170  {
171  forAll(patchNonConformalCyclicPatches_[patchi], i)
172  {
173  const label nccPatchi =
174  patchNonConformalCyclicPatches_[patchi][i];
175 
176  const nonConformalCyclicPolyPatch& nccPp =
177  refCast<const nonConformalCyclicPolyPatch>(pbm[nccPatchi]);
178 
179  if (nccPp.owner()) nccPp.rays();
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
186 
187 template<class ParticleType>
189 (
190  const polyMesh& pMesh,
191  const word& cloudName,
192  const IDLList<ParticleType>& particles
193 )
194 :
195  cloud(pMesh, cloudName),
196  IDLList<ParticleType>(),
197  pMesh_(pMesh),
198  patchNbrProc_(patchNbrProc(pMesh)),
199  patchNbrProcPatch_(patchNbrProcPatch(pMesh)),
200  patchNonConformalCyclicPatches_(patchNonConformalCyclicPatches(pMesh)),
201  globalPositionsPtr_(),
202  timeIndex_(-1)
203 {
204  // Request the tet base points so that they are built on all processors.
205  // Constructing tet base points requires communication, so we can't leave
206  // it until the tracking requests them as those calls are not synchronised.
207  // Some processors might not be doing any tracking at all.
208  pMesh_.tetBasePtIs();
209 
210  // Mark the need to store the old-time cell-centres if the mesh is moving
211  if (!ParticleType::instantaneous)
212  {
213  pMesh_.oldCellCentres();
214  }
215 
216  if (particles.size())
217  {
219  }
220 }
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
225 template<class ParticleType>
227 {
228  this->append(pPtr);
229 }
230 
231 
232 template<class ParticleType>
234 {
235  delete(this->remove(&p));
236 }
237 
238 
239 template<class ParticleType>
241 {
242  label lostCount = 0;
243 
244  forAllIter(typename Cloud<ParticleType>, *this, pIter)
245  {
246  if (pIter().cell() == -1)
247  {
248  deleteParticle(pIter());
249  lostCount ++;
250  }
251  }
252 
253  reduce(lostCount, sumOp<label>());
254  if (lostCount != 0)
255  {
257  << "Cloud " << this->name()
258  << " deleted " << lostCount << " lost particles" << endl;
259  }
260 }
261 
262 
263 template<class ParticleType>
265 (
266  const Cloud<ParticleType>& c
267 )
268 {
269  // Reset particle count and particles only
270  // - not changing the cloud object registry or reference to the polyMesh
271  ParticleType::particleCount = 0;
273 }
274 
275 
276 template<class ParticleType>
278 {
279  forAllIter(typename Cloud<ParticleType>, *this, pIter)
280  {
281  pIter().reset(0);
282  }
283 
284  timeIndex_ = pMesh_.time().timeIndex();
285 }
286 
287 
288 template<class ParticleType>
289 template<class TrackCloudType>
291 (
292  TrackCloudType& cloud,
293  typename ParticleType::trackingData& td
294 )
295 {
296  // If the time has changed, modify the particles accordingly
297  if (!ParticleType::instantaneous && timeIndex_ != pMesh_.time().timeIndex())
298  {
299  changeTimeStep();
300  }
301 
302  // Clear the global positions as these are about to change
303  globalPositionsPtr_.clear();
304 
305  // Ensure rays are available for non conformal transfers
306  storeRays();
307 
308  // Create transfer buffers
310 
311  // Create lists of particles and patch indices to transfer
313  List<DynamicList<label>> sendPatchIndices(Pstream::nProcs());
314 
315  optionalCpuLoad& cloudCpuTime
316  (
317  optionalCpuLoad::New(name() + ":cpuLoad", pMesh_, cloud.cpuLoad())
318  );
319 
320  // While there are particles to transfer
321  while (true)
322  {
323  // Clear the transfer lists
324  forAll(sendParticles, proci)
325  {
326  sendParticles[proci].clear();
327  sendPatchIndices[proci].clear();
328  }
329 
330  if (cloud.cpuLoad())
331  {
332  cloudCpuTime.resetCpuTime();
333  }
334 
335  // Loop over all particles
336  forAllIter(typename Cloud<ParticleType>, *this, pIter)
337  {
338  ParticleType& p = pIter();
339 
340  // Move the particle
341  const bool keepParticle = p.move(cloud, td);
342 
343  if (cloud.cpuLoad())
344  {
345  cloudCpuTime.cpuTimeIncrement(p.cell());
346  }
347 
348  // If the particle is to be kept
349  if (keepParticle)
350  {
351  if (td.sendToProc != -1)
352  {
353  #ifdef FULLDEBUG
354  if (!Pstream::parRun() || !p.onBoundaryFace(pMesh_))
355  {
357  << "Switch processor flag is true when no parallel "
358  << "transfer is possible. This is a bug."
359  << exit(FatalError);
360  }
361  #endif
362 
363  p.prepareForParallelTransfer(cloud, td);
364 
365  sendParticles[td.sendToProc].append(this->remove(&p));
366 
367  sendPatchIndices[td.sendToProc].append(td.sendToPatch);
368  }
369  }
370  else
371  {
372  deleteParticle(p);
373  }
374  }
375 
376  // If running in serial then everything has been moved, so finish
377  if (!Pstream::parRun())
378  {
379  break;
380  }
381 
382  // Clear transfer buffers
383  pBufs.clear();
384 
385  // Stream into send buffers
386  forAll(sendParticles, proci)
387  {
388  if (sendParticles[proci].size())
389  {
390  UOPstream particleStream(proci, pBufs);
391 
392  particleStream
393  << sendPatchIndices[proci]
394  << sendParticles[proci];
395  }
396  }
397 
398  // Start sending. Sets number of bytes transferred.
399  labelList receiveSizes(Pstream::nProcs());
400  pBufs.finishedSends(receiveSizes);
401 
402  // Determine if any particles were transferred. If not, then finish.
403  bool transferred = false;
404  forAll(receiveSizes, proci)
405  {
406  if (receiveSizes[proci])
407  {
408  transferred = true;
409  break;
410  }
411  }
412  reduce(transferred, orOp<bool>());
413  if (!transferred)
414  {
415  break;
416  }
417 
418  // Retrieve from receive buffers and add into the cloud
419  forAll(receiveSizes, proci)
420  {
421  if (receiveSizes[proci])
422  {
423  UIPstream particleStream(proci, pBufs);
424 
425  const labelList receivePatchIndices(particleStream);
426 
427  IDLList<ParticleType> newParticles(particleStream);
428 
429  label i = 0;
430 
431  forAllIter(typename Cloud<ParticleType>, newParticles, iter)
432  {
433  const label patchi = receivePatchIndices[i ++];
434 
435  ParticleType& p = iter();
436 
437  td.sendToPatch = patchi;
438 
439  p.correctAfterParallelTransfer(cloud, td);
440 
441  addParticle(newParticles.remove(&p));
442  }
443  }
444  }
445  }
446 
447  // Warn about any approximate locates
448  Pstream::listCombineGather(td.patchNLocateBoundaryHits, plusEqOp<label>());
449  if (Pstream::master())
450  {
451  forAll(td.patchNLocateBoundaryHits, patchi)
452  {
453  if (td.patchNLocateBoundaryHits[patchi] != 0)
454  {
456  << "Cloud " << name() << " did not accurately locate "
457  << td.patchNLocateBoundaryHits[patchi]
458  << " particles that transferred to patch "
459  << pMesh_.boundaryMesh()[patchi].name() << nl;
460  }
461  }
462  }
463 }
464 
465 
466 template<class ParticleType>
468 (
469  const polyTopoChangeMap& map
470 )
471 {
472  if (map.reverseCellMap().empty()) return;
473 
474  // See comments in the constructor
475  pMesh_.tetBasePtIs();
476  pMesh_.oldCellCentres();
477 
478  if (!globalPositionsPtr_.valid())
479  {
481  << "Global positions are not available. "
482  << "Cloud::storeGlobalPositions has not been called."
483  << exit(FatalError);
484  }
485 
486  const vectorField& positions = globalPositionsPtr_();
487 
488  label lostCount = 0;
489 
490  label particlei = 0;
491  forAllIter(typename Cloud<ParticleType>, *this, iter)
492  {
493  const point& pos = positions[particlei ++];
494 
495  const label celli = map.reverseCellMap()[iter().cell()];
496 
497  if (!iter().locate(pMesh_, pos, celli))
498  {
499  this->remove(iter);
500  lostCount ++;
501  }
502  }
503 
504  reduce(lostCount, sumOp<label>());
505  if (lostCount != 0)
506  {
508  << "Topology change of cloud " << this->name()
509  << " lost " << lostCount << " particles" << endl;
510  }
511 }
512 
513 
514 template<class ParticleType>
516 {
517  // See comments in the constructor
518  pMesh_.tetBasePtIs();
519  pMesh_.oldCellCentres();
520 
521  // Update cached mesh indexing
522  patchNbrProc_ = patchNbrProc(pMesh_);
523  patchNbrProcPatch_ = patchNbrProcPatch(pMesh_);
524  patchNonConformalCyclicPatches_ = patchNonConformalCyclicPatches(pMesh_);
525 
526  if (!globalPositionsPtr_.valid())
527  {
529  << "Global positions are not available. "
530  << "Cloud::storeGlobalPositions has not been called."
531  << exit(FatalError);
532  }
533 
534  const vectorField& positions = globalPositionsPtr_();
535 
536  label lostCount = 0;
537 
538  // Loop the particles. Map those that remain on this processor, and
539  // transfer others into send arrays.
540  List<DynamicList<label>> sendCellIndices(Pstream::nProcs());
541  List<DynamicList<point>> sendPositions(Pstream::nProcs());
543  {
544  label particlei = 0;
545  forAllIter(typename Cloud<ParticleType>, *this, iter)
546  {
547  const point& pos = positions[particlei ++];
548 
549  const remote tgtProcCell =
550  map.mapper().srcToTgtPoint(iter().cell(), pos);
551  const label proci = tgtProcCell.proci;
552  const label celli = tgtProcCell.elementi;
553 
554  if (tgtProcCell == remote())
555  {
556  this->remove(iter);
557  lostCount ++;
558  }
559  else if (proci == Pstream::myProcNo())
560  {
561  if (!iter().locate(pMesh_, pos, celli))
562  {
563  this->remove(iter);
564  lostCount ++;
565  }
566  }
567  else
568  {
569  sendCellIndices[proci].append(celli);
570  sendPositions[proci].append(pos);
571  sendParticles[proci].append(this->remove(iter));
572  }
573  }
574  }
575 
576  // If parallel then send and receive particles that move processes and map
577  // those sent to this process
578  if (Pstream::parRun())
579  {
580  // Create transfer buffers
582 
583  // Stream into send buffers
584  forAll(sendParticles, proci)
585  {
586  if (sendParticles[proci].size())
587  {
588  UOPstream particleStream(proci, pBufs);
589 
590  particleStream
591  << sendCellIndices[proci]
592  << sendPositions[proci]
593  << sendParticles[proci];
594  }
595  }
596 
597  // Finish sending
598  labelList receiveSizes(Pstream::nProcs());
599  pBufs.finishedSends(receiveSizes);
600 
601  // Retrieve from receive buffers and map into the new mesh
602  forAll(sendParticles, proci)
603  {
604  if (receiveSizes[proci])
605  {
606  UIPstream particleStream(proci, pBufs);
607 
608  const labelList receiveCellIndices(particleStream);
609  const List<point> receivePositions(particleStream);
610  IDLList<ParticleType> receiveParticles(particleStream);
611 
612  label particlei = 0;
613  forAllIter(typename Cloud<ParticleType>, receiveParticles, iter)
614  {
615  const label celli = receiveCellIndices[particlei];
616  const vector& pos = receivePositions[particlei ++];
617 
618  if (iter().locate(pMesh_, pos, celli))
619  {
620  this->append(receiveParticles.remove(iter));
621  }
622  else
623  {
624  receiveParticles.remove(iter);
625  lostCount ++;
626  }
627  }
628  }
629  }
630  }
631 
632  reduce(lostCount, sumOp<label>());
633  if (lostCount != 0)
634  {
636  << "Mesh-to-mesh mapping of cloud " << this->name()
637  << " lost " << lostCount << " particles" << endl;
638  }
639 }
640 
641 
642 template<class ParticleType>
644 (
645  const polyDistributionMap& map
646 )
647 {
648  // See comments in the constructor
649  pMesh_.tetBasePtIs();
650  pMesh_.oldCellCentres();
651 
652  // Update cached mesh indexing
653  patchNbrProc_ = patchNbrProc(pMesh_);
654  patchNbrProcPatch_ = patchNbrProcPatch(pMesh_);
655  patchNonConformalCyclicPatches_ = patchNonConformalCyclicPatches(pMesh_);
656 
657  if (!globalPositionsPtr_.valid())
658  {
660  << "Global positions are not available. "
661  << "Cloud::storeGlobalPositions has not been called."
662  << exit(FatalError);
663  }
664 
665  const vectorField& positions = globalPositionsPtr_();
666 
667  // Distribute the global positions
668  List<List<point>> cellParticlePositions(map.nOldCells());
669  {
670  labelList cellParticleis(map.nOldCells(), 0);
671  forAllIter(typename Cloud<ParticleType>, *this, iter)
672  {
673  cellParticleis[iter().cell()] ++;
674  }
675  forAll(cellParticlePositions, celli)
676  {
677  cellParticlePositions[celli].resize(cellParticleis[celli]);
678  }
679 
680  label particlei = 0;
681  cellParticleis = 0;
682  forAllIter(typename Cloud<ParticleType>, *this, iter)
683  {
684  const label celli = iter().cell();
685  label& cellParticlei = cellParticleis[celli];
686 
687  cellParticlePositions[celli][cellParticlei ++] =
688  positions[particlei ++];
689  }
690  }
692  (
694  List<labelPair>(),
695  pMesh_.nCells(),
696  map.cellMap().subMap(),
697  false,
698  map.cellMap().constructMap(),
699  false,
700  cellParticlePositions,
702  flipOp(),
703  List<point>()
704  );
705 
706  // Distribute the particles
707  List<IDLList<ParticleType>> cellParticles(map.nOldCells());
708  forAllIter(typename Cloud<ParticleType>, *this, iter)
709  {
710  cellParticles[iter().cell()].append(this->remove(iter));
711  }
713  (
715  List<labelPair>(),
716  pMesh_.nCells(),
717  map.cellMap().subMap(),
718  false,
719  map.cellMap().constructMap(),
720  false,
721  cellParticles,
723  flipOp(),
725  );
726 
727  label lostCount = 0;
728 
729  // Locate the particles within the new mesh
730  forAll(cellParticles, celli)
731  {
732  label cellParticlei = 0;
733  forAllIter(typename IDLList<ParticleType>, cellParticles[celli], iter)
734  {
735  const point& pos = cellParticlePositions[celli][cellParticlei++];
736 
737  if (iter().locate(pMesh_, pos, celli))
738  {
739  this->append(cellParticles[celli].remove(iter));
740  }
741  else
742  {
743  cellParticles[celli].remove(iter);
744  lostCount ++;
745  }
746  }
747  }
748 
749  reduce(lostCount, sumOp<label>());
750  if (lostCount != 0)
751  {
753  << "Distribution of cloud " << this->name()
754  << " lost " << lostCount << " particles" << endl;
755  }
756 }
757 
758 
759 template<class ParticleType>
761 {
762  OFstream pObj
763  (
764  this->db().time().path()/this->name() + "_positions.obj"
765  );
766 
767  forAllConstIter(typename Cloud<ParticleType>, *this, pIter)
768  {
769  const point pos = pIter().position(pMesh_);
770 
771  pObj<< "v " << pos.x() << " " << pos.y() << " " << pos.z() << nl;
772  }
773 
774  pObj.flush();
775 }
776 
777 
778 template<class ParticleType>
780 {
781  // Store the global positions for later use by mapping functions. It would
782  // be preferable not to need this. If the objects passed to had a copy of
783  // the old mesh then the global positions could be recovered within the
784  // mapping functions, and this pre-processing would not be necessary.
785 
786  globalPositionsPtr_.reset(new vectorField(this->size()));
787 
788  vectorField& positions = globalPositionsPtr_();
789 
790  label particlei = 0;
791  forAllConstIter(typename Cloud<ParticleType>, *this, iter)
792  {
793  positions[particlei++] = iter().position(pMesh_);
794  }
795 }
796 
797 
798 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
799 
800 #include "CloudIO.C"
801 
802 // ************************************************************************* //
scalar y
Combination-Reduction operation for a parallel run. The information from all nodes is collected on th...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Base cloud calls templated on particle type.
Template class for intrusive linked lists.
Definition: ILList.H:67
void operator=(const ILList< LListBase, T > &)
Assignment operator.
Definition: ILList.C:144
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
Output to file stream.
Definition: OFstream.H:87
virtual void flush()
Flush stream.
Definition: OSstream.C:223
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
void clear()
Clear storage and reset.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
T * remove(T *p)
Remove and return element.
Definition: UILList.H:142
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:57
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:58
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
const labelListList & subMap() const
From subsetted data back to original data.
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &, const negateOp &negOp, const int tag=UPstream::msgType())
Distribute data. Note:schedule only used for.
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:51
void changeTimeStep()
Change the particles' state from the end of the previous time.
Definition: Cloud.C:277
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
Definition: Cloud.C:233
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:760
const labelListList & patchNonConformalCyclicPatches() const
Return map from patch index to connected non-conformal cyclics.
Definition: Cloud.H:185
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: Cloud.C:468
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: Cloud.C:644
const labelList & patchNbrProcPatch() const
Map from patch index to the corresponding patch index on the.
Definition: Cloud.H:179
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:240
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: Cloud.C:515
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition: Cloud.C:265
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
Definition: Cloud.C:189
const labelList & patchNbrProc() const
Map from patch index to the neighbouring processor index.
Definition: Cloud.H:172
void storeGlobalPositions() const
Call this before a topology change. Stores the particles global.
Definition: Cloud.C:779
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:226
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:291
remote srcToTgtPoint(const label srcCelli, const point &p) const
Find the target processor and cell associated with a point in a.
Definition: meshToMesh.C:249
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
virtual bool owner() const
Inherit the cyclic owner method.
const patchToPatches::rays & rays() const
Access the rays engine.
label origPatchIndex() const
Original patchID.
virtual void resetCpuTime()
Reset the CPU time (dummy)
Definition: cpuLoad.H:101
virtual void cpuTimeIncrement(const label celli)
Cache the CPU time increment for celli (dummy)
Definition: cpuLoad.H:105
static optionalCpuLoad & New(const word &name, const polyMesh &mesh, const bool loadBalancing)
Construct from polyMesh if loadBalancing is true.
Definition: cpuLoad.C:63
Foam::polyBoundaryMesh.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
label nOldCells() const
Number of cells in mesh before distribution.
const distributionMap & cellMap() const
Cell distribute map.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
const meshToMesh & mapper() const
Return meshToMesh mapper.
Definition: polyMeshMap.H:81
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1046
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1387
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const
Reverse cell map.
Struct for keeping processor, element (cell, face, point) index.
Definition: remote.H:57
label elementi
Element index.
Definition: remote.H:70
label proci
Processor index.
Definition: remote.H:67
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
bool locate(const polyMesh &mesh, const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction, const string &debugPrefix=NullObjectRef< string >())
Initialise the location at the given position. Returns whether or not a.
Definition: tracking.C:1593
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar pos(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
volScalarField & p
void operator()(IDLList< Type > &x, const IDLList< Type > &y) const
Definition: Cloud.C:48
List operator to append one list onto another.
Definition: ListOps.H:304
const word cloudName(propsDict.lookup("cloudName"))