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-2023 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 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 template<class Type>
46 {
47  void operator()(IDLList<Type>& x, const IDLList<Type>& y) const
48  {
49  if (y.size())
50  {
51  forAllConstIter(typename IDLList<Type>, y, iter)
52  {
53  x.append(new Type(iter()));
54  }
55  }
56  }
57 };
58 
59 }
60 
61 
62 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
63 
64 template<class ParticleType>
66 (
67  const polyMesh& pMesh
68 )
69 {
70  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
71 
72  labelList result(pbm.size(), -1);
73 
74  if (Pstream::parRun())
75  {
76  forAll(pbm, patchi)
77  {
78  if (isA<processorPolyPatch>(pbm[patchi]))
79  {
80  const processorPolyPatch& ppp =
81  refCast<const processorPolyPatch>(pbm[patchi]);
82 
83  result[patchi] = ppp.neighbProcNo();
84  }
85  }
86  }
87 
88  return result;
89 }
90 
91 
92 template<class ParticleType>
94 (
95  const polyMesh& pMesh
96 )
97 {
98  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
99 
100  labelList result(pbm.size(), -1);
101 
102  if (Pstream::parRun())
103  {
104  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
105 
106  forAll(pbm, patchi)
107  {
108  if (isA<processorPolyPatch>(pbm[patchi]))
109  {
110  const processorPolyPatch& ppp =
111  refCast<const processorPolyPatch>(pbm[patchi]);
112 
113  UOPstream(ppp.neighbProcNo(), pBufs)()
114  << ppp.index();
115  }
116  }
117 
118  pBufs.finishedSends();
119 
120  forAll(pbm, patchi)
121  {
122  if (isA<processorPolyPatch>(pbm[patchi]))
123  {
124  const processorPolyPatch& ppp =
125  refCast<const processorPolyPatch>(pbm[patchi]);
126 
127  UIPstream(ppp.neighbProcNo(), pBufs)()
128  >> result[patchi];
129  }
130  }
131  }
132 
133  return result;
134 }
135 
136 
137 template<class ParticleType>
139 (
140  const polyMesh& pMesh
141 )
142 {
143  const polyBoundaryMesh& pbm = pMesh.boundaryMesh();
144 
145  labelListList result(pbm.size(), labelList());
146 
147  forAll(pbm, patchi)
148  {
149  if (isA<nonConformalCyclicPolyPatch>(pbm[patchi]))
150  {
151  const nonConformalCyclicPolyPatch& nccPp =
152  refCast<const nonConformalCyclicPolyPatch>(pbm[patchi]);
153 
154  result[nccPp.origPatchID()].append(patchi);
155  }
156  }
157 
158  return result;
159 }
160 
161 
162 template<class ParticleType>
164 {
165  const polyBoundaryMesh& pbm = pMesh_.boundaryMesh();
166 
167  forAll(patchNonConformalCyclicPatches_, patchi)
168  {
169  forAll(patchNonConformalCyclicPatches_[patchi], i)
170  {
171  const label nccPatchi =
172  patchNonConformalCyclicPatches_[patchi][i];
173 
174  const nonConformalCyclicPolyPatch& nccPp =
175  refCast<const nonConformalCyclicPolyPatch>(pbm[nccPatchi]);
176 
177  if (nccPp.owner()) nccPp.rays();
178  }
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
184 
185 template<class ParticleType>
187 (
188  const polyMesh& pMesh,
189  const word& cloudName,
190  const IDLList<ParticleType>& particles
191 )
192 :
193  cloud(pMesh, cloudName),
194  IDLList<ParticleType>(),
195  pMesh_(pMesh),
196  patchNbrProc_(patchNbrProc(pMesh)),
197  patchNbrProcPatch_(patchNbrProcPatch(pMesh)),
198  patchNonConformalCyclicPatches_(patchNonConformalCyclicPatches(pMesh)),
199  globalPositionsPtr_(),
200  timeIndex_(-1)
201 {
202  // Ask for the tetBasePtIs and oldCellCentres to trigger all processors to
203  // build them, otherwise, if some processors have no particles then there
204  // is a comms mismatch.
205  pMesh_.tetBasePtIs();
206  pMesh_.oldCellCentres();
207 
208  if (particles.size())
209  {
211  }
212 }
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 template<class ParticleType>
219 {
220  this->append(pPtr);
221 }
222 
223 
224 template<class ParticleType>
226 {
227  delete(this->remove(&p));
228 }
229 
230 
231 template<class ParticleType>
233 {
234  forAllIter(typename Cloud<ParticleType>, *this, pIter)
235  {
236  ParticleType& p = pIter();
237 
238  if (p.cell() == -1)
239  {
241  << "deleting lost particle at position "
242  << p.position(pMesh_) << endl;
243 
244  deleteParticle(p);
245  }
246  }
247 }
248 
249 
250 template<class ParticleType>
252 {
253  // Reset particle count and particles only
254  // - not changing the cloud object registry or reference to the polyMesh
255  ParticleType::particleCount_ = 0;
257 }
258 
259 
260 template<class ParticleType>
262 {
263  forAllIter(typename Cloud<ParticleType>, *this, pIter)
264  {
265  pIter().reset(0);
266  }
267 
268  timeIndex_ = pMesh_.time().timeIndex();
269 }
270 
271 
272 template<class ParticleType>
273 template<class TrackCloudType>
275 (
276  TrackCloudType& cloud,
277  typename ParticleType::trackingData& td
278 )
279 {
280  // If the time has changed, modify the particles accordingly
281  if (timeIndex_ != pMesh_.time().timeIndex())
282  {
283  changeTimeStep();
284  }
285 
286  // Clear the global positions as these are about to change
287  globalPositionsPtr_.clear();
288 
289  // Ensure rays are available for non conformal transfers
290  storeRays();
291 
292  // Create transfer buffers
294 
295  // Create lists of particles and patch indices to transfer
297  List<DynamicList<label>> sendPatchIndices(Pstream::nProcs());
298 
299  // While there are particles to transfer
300  while (true)
301  {
302  // Clear the transfer lists
303  forAll(sendParticles, proci)
304  {
305  sendParticles[proci].clear();
306  sendPatchIndices[proci].clear();
307  }
308 
309  // Loop over all particles
310  forAllIter(typename Cloud<ParticleType>, *this, pIter)
311  {
312  ParticleType& p = pIter();
313 
314  // Move the particle
315  const bool keepParticle = p.move(cloud, td);
316 
317  // If the particle is to be kept
318  if (keepParticle)
319  {
320  if (td.sendToProc != -1)
321  {
322  #ifdef FULLDEBUG
323  if (!Pstream::parRun() || !p.onBoundaryFace(pMesh_))
324  {
326  << "Switch processor flag is true when no parallel "
327  << "transfer is possible. This is a bug."
328  << exit(FatalError);
329  }
330  #endif
331 
332  p.prepareForParallelTransfer(cloud, td);
333 
334  sendParticles[td.sendToProc].append(this->remove(&p));
335 
336  sendPatchIndices[td.sendToProc].append(td.sendToPatch);
337  }
338  }
339  else
340  {
341  deleteParticle(p);
342  }
343  }
344 
345  // If running in serial then everything has been moved, so finish
346  if (!Pstream::parRun())
347  {
348  break;
349  }
350 
351  // Clear transfer buffers
352  pBufs.clear();
353 
354  // Stream into send buffers
355  forAll(sendParticles, proci)
356  {
357  if (sendParticles[proci].size())
358  {
359  UOPstream particleStream(proci, pBufs);
360 
361  particleStream
362  << sendPatchIndices[proci]
363  << sendParticles[proci];
364  }
365  }
366 
367  // Start sending. Sets number of bytes transferred.
368  labelList receiveSizes(Pstream::nProcs());
369  pBufs.finishedSends(receiveSizes);
370 
371  // Determine if any particles were transferred. If not, then finish.
372  bool transferred = false;
373  forAll(receiveSizes, proci)
374  {
375  if (receiveSizes[proci])
376  {
377  transferred = true;
378  break;
379  }
380  }
381  reduce(transferred, orOp<bool>());
382  if (!transferred)
383  {
384  break;
385  }
386 
387  // Retrieve from receive buffers and add into the cloud
388  forAll(receiveSizes, proci)
389  {
390  if (receiveSizes[proci])
391  {
392  UIPstream particleStream(proci, pBufs);
393 
394  const labelList receivePatchIndices(particleStream);
395 
396  IDLList<ParticleType> newParticles(particleStream);
397 
398  label i = 0;
399 
400  forAllIter(typename Cloud<ParticleType>, newParticles, iter)
401  {
402  const label patchi = receivePatchIndices[i ++];
403 
404  ParticleType& p = iter();
405 
406  td.sendToPatch = patchi;
407 
408  p.correctAfterParallelTransfer(cloud, td);
409 
410  addParticle(newParticles.remove(&p));
411  }
412  }
413  }
414  }
415 }
416 
417 
418 template<class ParticleType>
420 {
421  if (map.reverseCellMap().empty()) return;
422 
423  // Ask for the tetBasePtIs to trigger all processors to build
424  // them, otherwise, if some processors have no particles then
425  // there is a comms mismatch.
426  pMesh_.tetBasePtIs();
427  pMesh_.oldCellCentres();
428 
429  if (!globalPositionsPtr_.valid())
430  {
432  << "Global positions are not available. "
433  << "Cloud::storeGlobalPositions has not been called."
434  << exit(FatalError);
435  }
436 
437  const vectorField& positions = globalPositionsPtr_();
438 
439  label particlei = 0;
440  forAllIter(typename Cloud<ParticleType>, *this, iter)
441  {
442  const label celli = map.reverseCellMap()[iter().cell()];
443  iter().map(pMesh_, positions[particlei++], celli);
444  }
445 }
446 
447 
448 template<class ParticleType>
450 {
451  // Ask for the tetBasePtIs to trigger all processors to build
452  // them, otherwise, if some processors have no particles then
453  // there is a comms mismatch.
454  pMesh_.tetBasePtIs();
455  pMesh_.oldCellCentres();
456 
457  // Update cached mesh indexing
458  patchNbrProc_ = patchNbrProc(pMesh_);
459  patchNbrProcPatch_ = patchNbrProcPatch(pMesh_);
460  patchNonConformalCyclicPatches_ = patchNonConformalCyclicPatches(pMesh_);
461 
462  if (!globalPositionsPtr_.valid())
463  {
465  << "Global positions are not available. "
466  << "Cloud::storeGlobalPositions has not been called."
467  << exit(FatalError);
468  }
469 
470  const vectorField& positions = globalPositionsPtr_();
471 
472  // Loop the particles. Map those that remain on this processor, and
473  // transfer others into send arrays.
474  List<DynamicList<label>> sendCellIndices(Pstream::nProcs());
475  List<DynamicList<point>> sendPositions(Pstream::nProcs());
477  {
478  label particlei = 0;
479  forAllIter(typename Cloud<ParticleType>, *this, iter)
480  {
481  const point& pos = positions[particlei ++];
482 
483  const remote tgtProcCell =
484  map.mapper().srcToTgtPoint(iter().cell(), pos);
485 
486  if (tgtProcCell == remote())
487  {
489  << "Particle at " << pos << " mapped to a location outside "
490  << "of the new mesh. This particle will be removed." << nl;
491  }
492  else if (tgtProcCell.proci == Pstream::myProcNo())
493  {
494  iter().map(pMesh_, pos, tgtProcCell.elementi);
495  }
496  else
497  {
498  sendCellIndices[tgtProcCell.proci].append(tgtProcCell.elementi);
499  sendPositions[tgtProcCell.proci].append(pos);
500  sendParticles[tgtProcCell.proci].append(this->remove(iter));
501  }
502  }
503  }
504 
505  // If serial then there is nothing more to do
506  if (!Pstream::parRun())
507  {
508  return;
509  }
510 
511  // Create transfer buffers
513 
514  // Stream into send buffers
515  forAll(sendParticles, proci)
516  {
517  if (sendParticles[proci].size())
518  {
519  UOPstream particleStream(proci, pBufs);
520 
521  particleStream
522  << sendCellIndices[proci]
523  << sendPositions[proci]
524  << sendParticles[proci];
525  }
526  }
527 
528  // Finish sending
529  labelList receiveSizes(Pstream::nProcs());
530  pBufs.finishedSends(receiveSizes);
531 
532  // Retrieve from receive buffers and map into the new mesh
533  forAll(sendParticles, proci)
534  {
535  if (receiveSizes[proci])
536  {
537  UIPstream particleStream(proci, pBufs);
538 
539  const labelList receiveCellIndices(particleStream);
540  const List<point> receivePositions(particleStream);
541  IDLList<ParticleType> receiveParticles(particleStream);
542 
543  label particlei = 0;
544  forAllIter(typename Cloud<ParticleType>, receiveParticles, iter)
545  {
546  const label celli = receiveCellIndices[particlei];
547  const vector& pos = receivePositions[particlei ++];
548 
549  iter().map(pMesh_, pos, celli);
550  this->append(receiveParticles.remove(iter));
551  }
552  }
553  }
554 }
555 
556 
557 template<class ParticleType>
559 {
560  // Ask for the tetBasePtIs to trigger all processors to build
561  // them, otherwise, if some processors have no particles then
562  // there is a comms mismatch.
563  pMesh_.tetBasePtIs();
564  pMesh_.oldCellCentres();
565 
566  // Update cached mesh indexing
567  patchNbrProc_ = patchNbrProc(pMesh_);
568  patchNbrProcPatch_ = patchNbrProcPatch(pMesh_);
569  patchNonConformalCyclicPatches_ = patchNonConformalCyclicPatches(pMesh_);
570 
571  if (!globalPositionsPtr_.valid())
572  {
574  << "Global positions are not available. "
575  << "Cloud::storeGlobalPositions has not been called."
576  << exit(FatalError);
577  }
578 
579  const vectorField& positions = globalPositionsPtr_();
580 
581  // Distribute the global positions
582  List<List<point>> cellParticlePositions(map.nOldCells());
583  {
584  labelList cellParticleis(map.nOldCells(), 0);
585  forAllIter(typename Cloud<ParticleType>, *this, iter)
586  {
587  cellParticleis[iter().cell()] ++;
588  }
589  forAll(cellParticlePositions, celli)
590  {
591  cellParticlePositions[celli].resize(cellParticleis[celli]);
592  }
593 
594  label particlei = 0;
595  cellParticleis = 0;
596  forAllIter(typename Cloud<ParticleType>, *this, iter)
597  {
598  const label celli = iter().cell();
599  label& cellParticlei = cellParticleis[celli];
600 
601  cellParticlePositions[celli][cellParticlei ++] =
602  positions[particlei ++];
603  }
604  }
606  (
608  List<labelPair>(),
609  pMesh_.nCells(),
610  map.cellMap().subMap(),
611  false,
612  map.cellMap().constructMap(),
613  false,
614  cellParticlePositions,
616  flipOp(),
617  List<point>()
618  );
619 
620  // Distribute the particles
621  List<IDLList<ParticleType>> cellParticles(map.nOldCells());
622  forAllIter(typename Cloud<ParticleType>, *this, iter)
623  {
624  cellParticles[iter().cell()].append(this->remove(iter));
625  }
627  (
629  List<labelPair>(),
630  pMesh_.nCells(),
631  map.cellMap().subMap(),
632  false,
633  map.cellMap().constructMap(),
634  false,
635  cellParticles,
637  flipOp(),
639  );
640 
641  // Locate the particles within the new mesh
642  forAll(cellParticles, celli)
643  {
644  label cellParticlei = 0;
645  forAllIter(typename IDLList<ParticleType>, cellParticles[celli], iter)
646  {
647  const point& pos = cellParticlePositions[celli][cellParticlei++];
648 
649  iter().map(pMesh_, pos, celli);
650 
651  this->append(cellParticles[celli].remove(iter));
652  }
653  }
654 }
655 
656 
657 template<class ParticleType>
659 {
660  OFstream pObj
661  (
662  this->db().time().path()/this->name() + "_positions.obj"
663  );
664 
665  forAllConstIter(typename Cloud<ParticleType>, *this, pIter)
666  {
667  const point pos = pIter().position(pMesh_);
668 
669  pObj<< "v " << pos.x() << " " << pos.y() << " " << pos.z() << nl;
670  }
671 
672  pObj.flush();
673 }
674 
675 
676 template<class ParticleType>
678 {
679  // Store the global positions for later use by mapping functions. It would
680  // be preferable not to need this. If the objects passed to had a copy of
681  // the old mesh then the global positions could be recovered within the
682  // mapping functions, and this pre-processing would not be necessary.
683 
684  globalPositionsPtr_.reset(new vectorField(this->size()));
685 
686  vectorField& positions = globalPositionsPtr_();
687 
688  label particlei = 0;
689  forAllConstIter(typename Cloud<ParticleType>, *this, iter)
690  {
691  positions[particlei++] = iter().position(pMesh_);
692  }
693 }
694 
695 
696 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
697 
698 #include "CloudIO.C"
699 
700 // ************************************************************************* //
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:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Base cloud calls templated on particle type.
Definition: Cloud.H:74
void changeTimeStep()
Change the particles' state from the end of the previous time.
Definition: Cloud.C:261
void deleteParticle(ParticleType &)
Remove particle from cloud and delete.
Definition: Cloud.C:225
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:658
const labelListList & patchNonConformalCyclicPatches() const
Return map from patch index to connected non-conformal cyclics.
Definition: Cloud.H:184
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: Cloud.C:419
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: Cloud.C:558
const labelList & patchNbrProcPatch() const
Map from patch index to the corresponding patch index on the.
Definition: Cloud.H:178
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:232
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: Cloud.C:449
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition: Cloud.C:251
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
Definition: Cloud.C:187
const labelList & patchNbrProc() const
Map from patch index to the neighbouring processor index.
Definition: Cloud.H:171
void storeGlobalPositions() const
Call this before a topology change. Stores the particles global.
Definition: Cloud.C:677
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:218
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td)
Move the particles.
Definition: Cloud.C:275
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:86
virtual void flush()
Flush stream.
Definition: OSstream.C:207
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.
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 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
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
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
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...
label origPatchID() const
Original patchID.
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:1077
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1416
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:66
label proci
Processor index.
Definition: remote.H:63
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
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:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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.
error FatalError
static const char nl
Definition: Ostream.H:260
volScalarField & p
void operator()(IDLList< Type > &x, const IDLList< Type > &y) const
Definition: Cloud.C:47
List operator to append one list onto another.
Definition: ListOps.H:307
const word cloudName(propsDict.lookup("cloudName"))