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-2023 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"
32 #include "PstreamBuffers.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  bool ok = true;
45  for (const polyPatch& pp : polyMesh_.boundaryMesh())
46  {
47  const auto* camipp = isA<cyclicAMIPolyPatch>(pp);
48 
49  if (camipp && camipp->owner())
50  {
51  ok = (camipp->AMI().singlePatchProc() != -1);
52  if (!ok)
53  {
54  break;
55  }
56  }
57  }
58 
59  if (!ok)
60  {
62  << "Particle tracking across AMI patches is only currently "
63  << "supported for cases where the AMI patches reside on a "
64  << "single processor" << abort(FatalError);
65  }
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class ParticleType>
73 (
74  const polyMesh& pMesh,
75  const word& cloudName,
76  const IDLList<ParticleType>& particles
77 )
78 :
79  cloud(pMesh, cloudName),
80  IDLList<ParticleType>(),
81  polyMesh_(pMesh),
82  labels_(),
83  globalPositionsPtr_(),
84  geometryType_(cloud::geometryType::COORDINATES)
85 {
86  checkPatches();
87  (void)polyMesh_.oldCellCentres();
88 
89  // Ask for the tetBasePtIs to trigger all processors to build
90  // them, otherwise, if some processors have no particles then
91  // there is a comms mismatch.
92  (void)polyMesh_.tetBasePtIs();
93 
94  if (particles.size())
95  {
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 template<class ParticleType>
105 {
106  this->append(pPtr);
107 }
108 
109 
110 template<class ParticleType>
112 {
113  delete(this->remove(&p));
114 }
115 
116 
117 template<class ParticleType>
119 {
120  for (ParticleType& p : *this)
121  {
122  if (p.cell() == -1)
123  {
125  << "deleting lost particle at position " << p.position()
126  << endl;
127 
128  deleteParticle(p);
129  }
130  }
131 }
132 
133 
134 template<class ParticleType>
135 void Foam::Cloud<ParticleType>::cloudReset(const Cloud<ParticleType>& c)
136 {
137  // Reset particle count and particles only
138  // - not changing the cloud object registry or reference to the polyMesh
139  ParticleType::particleCount_ = 0;
141 }
142 
143 
144 template<class ParticleType>
145 template<class TrackCloudType>
147 (
148  TrackCloudType& cloud,
149  typename ParticleType::trackingData& td,
150  const scalar trackTime
151 )
152 {
153  const polyBoundaryMesh& pbm = pMesh().boundaryMesh();
154  const globalMeshData& pData = polyMesh_.globalData();
155 
156  // Which patches are processor patches
157  const labelList& procPatches = pData.processorPatches();
158 
159  // Indexing of equivalent patch on neighbour processor into the
160  // procPatches list on the neighbour
161  const labelList& procPatchNeighbours = pData.processorPatchNeighbours();
162 
163  // Which processors this processor is connected to
164  const labelList& neighbourProcs = pData.topology().procNeighbours();
165 
166  // Initialise the stepFraction moved for the particles
167  for (ParticleType& p : *this)
168  {
169  p.reset();
170  }
171 
172  // Clear the global positions as these are about to change
173  globalPositionsPtr_.clear();
174 
175 
176  // For v2112 and earlier: pre-assembled lists of particles
177  // to be transferred and target patch on a per processor basis.
178  // Apart from memory overhead of assembling the lists this adds
179  // allocations/de-allocation when building linked-lists.
180 
181  // Now stream particle transfer tuples directly into PstreamBuffers.
182  // Use a local cache of UOPstream wrappers for the formatters
183  // (since there are potentially many particles being shifted about).
184 
185 
186  // Allocate transfer buffers,
187  // automatic clearStorage when UIPstream closes is disabled.
188  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
189  pBufs.allowClearRecv(false);
190 
191  // Cache of opened UOPstream wrappers
192  PtrList<UOPstream> UOPstreamPtrs(Pstream::nProcs());
193 
194  // While there are particles to transfer
195  while (true)
196  {
197  // Reset transfer buffers
198  pBufs.clear();
199 
200  // Rewind existing streams
201  forAll(UOPstreamPtrs, proci)
202  {
203  auto* osptr = UOPstreamPtrs.get(proci);
204  if (osptr)
205  {
206  osptr->rewind();
207  }
208  }
209 
210  // Loop over all particles
211  for (ParticleType& p : *this)
212  {
213  // Move the particle
214  const 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 toProci =
240  (
241  refCast<const processorPolyPatch>(pbm[patchi])
242  .neighbProcNo()
243  );
244 
245  // Get/create output stream
246  auto* osptr = UOPstreamPtrs.get(toProci);
247  if (!osptr)
248  {
249  osptr = new UOPstream(toProci, pBufs);
250  UOPstreamPtrs.set(toProci, osptr);
251  }
252 
253  p.prepareForParallelTransfer();
254 
255  // Tuple: (patchi particle)
256  (*osptr) << procPatchNeighbours[patchi] << p;
257 
258  // Can now remove from my list
259  deleteParticle(p);
260  }
261  }
262  else
263  {
264  deleteParticle(p);
265  }
266  }
267 
268  if (!Pstream::parRun())
269  {
270  break;
271  }
272 
273  pBufs.finishedNeighbourSends(neighbourProcs);
274 
275  if (!returnReduceOr(pBufs.hasRecvData()))
276  {
277  // No parcels to transfer
278  break;
279  }
280 
281  // Retrieve from receive buffers
282  for (const label proci : neighbourProcs)
283  {
284  if (pBufs.recvDataCount(proci))
285  {
286  UIPstream is(proci, pBufs);
287 
288  // Read out each (patchi particle) tuple
289  while (!is.eof())
290  {
291  label patchi = pTraits<label>(is);
292  auto* newp = new ParticleType(polyMesh_, is);
293 
294  // The real patch index
295  patchi = procPatches[patchi];
296 
297  (*newp).correctAfterParallelTransfer(patchi, td);
298  addParticle(newp);
299  }
300  }
301  }
302  }
303 }
304 
305 
306 template<class ParticleType>
307 void Foam::Cloud<ParticleType>::autoMap(const mapPolyMesh& mapper)
308 {
309  if (!globalPositionsPtr_)
310  {
312  << "Global positions are not available. "
313  << "Cloud::storeGlobalPositions has not been called."
314  << exit(FatalError);
315  }
316 
317  // Reset stored data that relies on the mesh
318  // polyMesh_.clearCellTree();
319  cellWallFacesPtr_.clear();
320 
321  // Ask for the tetBasePtIs to trigger all processors to build
322  // them, otherwise, if some processors have no particles then
323  // there is a comms mismatch.
324  (void)polyMesh_.tetBasePtIs();
325  (void)polyMesh_.oldCellCentres();
326 
327  const vectorField& positions = globalPositionsPtr_();
328 
329  label i = 0;
330  for (ParticleType& p : *this)
331  {
332  p.autoMap(positions[i], mapper);
333  ++i;
334  }
335 }
336 
337 
338 template<class ParticleType>
340 {
341  OFstream os
342  (
343  this->db().time().path()/this->name() + "_positions.obj"
344  );
345 
346  for (const ParticleType& p : *this)
347  {
348  const point position(p.position());
349  os << "v "
350  << position.x() << ' '
351  << position.y() << ' '
352  << position.z() << nl;
353  }
354 }
355 
356 
357 template<class ParticleType>
359 {
360  // Store the global positions for later use by autoMap. It would be
361  // preferable not to need this. If the mapPolyMesh object passed to autoMap
362  // had a copy of the old mesh then the global positions could be recovered
363  // within autoMap, and this pre-processing would not be necessary.
364 
365  globalPositionsPtr_.reset(new vectorField(this->size()));
366  vectorField& positions = globalPositionsPtr_();
367 
368  label i = 0;
369  for (const ParticleType& p : *this)
370  {
371  positions[i] = p.position();
372  ++i;
373  }
374 }
375 
376 
377 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
378 
379 #include "CloudIO.C"
380 
381 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:45
const polyBoundaryMesh & pbm
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition: Cloud.C:128
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Cloud(const polyMesh &mesh, const word &cloudName, const IDLList< ParticleType > &particles)
Construct from mesh and a list of particles.
Definition: Cloud.C:66
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:300
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:97
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
void storeGlobalPositions() const
Call this before a topology change.
Definition: Cloud.C:351
const word cloudName(propsDict.get< word >("cloud"))
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:584
A class for handling words, derived from Foam::string.
Definition: word.H:63
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Base cloud calls templated on particle type.
Definition: Cloud.H:51
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
OBJstream os(runTime.globalPath()/outputName)
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:111
const labelList & processorPatchNeighbours() const noexcept
Return processorPatchIndices of the neighbours processor patches. -1 if not running parallel...
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
const labelList & procNeighbours() const
The neighbour processor connections (ascending order) associated with the local rank.
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:140
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
const processorTopology & topology() const noexcept
The processor to processor topology.
const dimensionedScalar c
Speed of light in a vacuum.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition: Cloud.C:104
volScalarField & p
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:332
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())