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  for (const polyPatch& pp : polyMesh_.boundaryMesh())
45  {
46  const auto* camipp = isA<cyclicAMIPolyPatch>(pp);
47 
48  if (camipp && camipp->owner() && camipp->AMI().distributed())
49  {
51  << "Particle tracking across AMI patches is only currently "
52  << "supported for cases where the AMI patches reside on a "
53  << "single processor" << abort(FatalError);
54  break;
55  }
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
62 template<class ParticleType>
64 (
65  const polyMesh& pMesh,
66  const Foam::zero,
67  const word& cloudName
68 )
69 :
70  cloud(pMesh, cloudName),
71  polyMesh_(pMesh),
72  geometryType_(cloud::geometryType::COORDINATES)
73 {
74  checkPatches();
75 
76  (void)polyMesh_.tetBasePtIs();
77  (void)polyMesh_.oldCellCentres();
78 }
79 
80 
81 template<class ParticleType>
83 (
84  const polyMesh& pMesh,
85  const word& cloudName,
86  const IDLList<ParticleType>& particles
87 )
88 :
89  Cloud<ParticleType>(pMesh, Foam::zero{}, cloudName)
90 {
91  if (particles.size())
92  {
93  IDLList<ParticleType>::operator=(particles);
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
100 template<class ParticleType>
102 {
103  this->append(pPtr);
104 }
105 
106 
107 template<class ParticleType>
109 {
110  delete(this->remove(&p));
111 }
112 
113 
114 template<class ParticleType>
116 {
117  for (ParticleType& p : *this)
118  {
119  if (p.cell() == -1)
120  {
122  << "deleting lost particle at position " << p.position()
123  << endl;
124 
125  deleteParticle(p);
126  }
127  }
128 }
129 
130 
131 template<class ParticleType>
132 void Foam::Cloud<ParticleType>::cloudReset(const Cloud<ParticleType>& c)
133 {
134  // Reset particle count and particles only
135  // - not changing the cloud object registry or reference to the polyMesh
136  ParticleType::particleCount_ = 0;
138 }
139 
140 
141 template<class ParticleType>
142 template<class TrackCloudType>
144 (
145  TrackCloudType& cloud,
146  typename ParticleType::trackingData& td,
147  const scalar trackTime
148 )
149 {
150  const polyBoundaryMesh& pbm = pMesh().boundaryMesh();
151  const globalMeshData& pData = polyMesh_.globalData();
152 
153  // Which patches are processor patches
154  const labelList& procPatches = pData.processorPatches();
155 
156  // Indexing of equivalent patch on neighbour processor into the
157  // procPatches list on the neighbour
158  const labelList& procPatchNeighbours = pData.processorPatchNeighbours();
159 
160  // Which processors this processor is connected to
161  const labelList& neighbourProcs = pData.topology().procNeighbours();
162 
163  // Initialise the stepFraction moved for the particles
164  for (ParticleType& p : *this)
165  {
166  p.reset();
167  }
168 
169  // Clear the global positions as these are about to change
170  globalPositionsPtr_.clear();
171 
172 
173  // For v2112 and earlier: pre-assembled lists of particles
174  // to be transferred and target patch on a per processor basis.
175  // Apart from memory overhead of assembling the lists this adds
176  // allocations/de-allocation when building linked-lists.
177 
178  // Now stream particle transfer tuples directly into PstreamBuffers.
179  // Use a local cache of UOPstream wrappers for the formatters
180  // (since there are potentially many particles being shifted about).
181 
182 
183  // Allocate transfer buffers,
184  // automatic clearStorage when UIPstream closes is disabled.
185  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
186  pBufs.allowClearRecv(false);
187 
188  // Cache of opened UOPstream wrappers
189  PtrList<UOPstream> UOPstreamPtrs(Pstream::nProcs());
190 
191  // While there are particles to transfer
192  while (true)
193  {
194  // Reset transfer buffers
195  pBufs.clear();
196 
197  // Rewind existing streams
198  forAll(UOPstreamPtrs, proci)
199  {
200  auto* osptr = UOPstreamPtrs.get(proci);
201  if (osptr)
202  {
203  osptr->rewind();
204  }
205  }
206 
207  // Loop over all particles
208  for (ParticleType& p : *this)
209  {
210  // Move the particle
211  const bool keepParticle = p.move(cloud, td, trackTime);
212 
213  // If the particle is to be kept
214  // (i.e. it hasn't passed through an inlet or outlet)
215  if (keepParticle)
216  {
217  if (td.switchProcessor)
218  {
219  #ifdef FULLDEBUG
220  if
221  (
222  !Pstream::parRun()
223  || !p.onBoundaryFace()
224  || procPatchNeighbours[p.patch()] < 0
225  )
226  {
228  << "Switch processor flag is true when no parallel "
229  << "transfer is possible. This is a bug."
230  << exit(FatalError);
231  }
232  #endif
233 
234  const label patchi = p.patch();
235 
236  const label toProci =
237  (
238  refCast<const processorPolyPatch>(pbm[patchi])
239  .neighbProcNo()
240  );
241 
242  // Get/create output stream
243  auto* osptr = UOPstreamPtrs.get(toProci);
244  if (!osptr)
245  {
246  osptr = new UOPstream(toProci, pBufs);
247  UOPstreamPtrs.set(toProci, osptr);
248  }
249 
250  p.prepareForParallelTransfer();
251 
252  // Tuple: (patchi particle)
253  (*osptr) << procPatchNeighbours[patchi] << p;
254 
255  // Can now remove from my list
256  deleteParticle(p);
257  }
258  }
259  else
260  {
261  deleteParticle(p);
262  }
263  }
264 
265  if (!Pstream::parRun())
266  {
267  break;
268  }
269 
270  pBufs.finishedNeighbourSends(neighbourProcs);
271 
272  if (!returnReduceOr(pBufs.hasRecvData()))
273  {
274  // No parcels to transfer
275  break;
276  }
277 
278  // Retrieve from receive buffers
279  for (const label proci : neighbourProcs)
280  {
281  if (pBufs.recvDataCount(proci))
282  {
283  UIPstream is(proci, pBufs);
284 
285  // Read out each (patchi particle) tuple
286  while (!is.eof())
287  {
288  label patchi = pTraits<label>(is);
289  auto* newp = new ParticleType(polyMesh_, is);
290 
291  // The real patch index
292  patchi = procPatches[patchi];
293 
294  (*newp).correctAfterParallelTransfer(patchi, td);
295  addParticle(newp);
296  }
297  }
298  }
299  }
300 }
301 
302 
303 template<class ParticleType>
304 void Foam::Cloud<ParticleType>::autoMap(const mapPolyMesh& mapper)
305 {
306  if (!globalPositionsPtr_)
307  {
309  << "Global positions are not available. "
310  << "Cloud::storeGlobalPositions has not been called."
311  << exit(FatalError);
312  }
313 
314  // Reset stored data that relies on the mesh
315  // polyMesh_.clearCellTree();
316  cellWallFacesPtr_.clear();
317 
318  // Ask for the tetBasePtIs to trigger all processors to build
319  // them, otherwise, if some processors have no particles then
320  // there is a comms mismatch.
321  (void)polyMesh_.tetBasePtIs();
322  (void)polyMesh_.oldCellCentres();
323 
324  const vectorField& positions = globalPositionsPtr_();
325 
326  label i = 0;
327  for (ParticleType& p : *this)
328  {
329  p.autoMap(positions[i], mapper);
330  ++i;
331  }
332 }
333 
334 
335 template<class ParticleType>
337 {
338  OFstream os
339  (
340  this->db().time().path()/this->name() + "_positions.obj"
341  );
342 
343  for (const ParticleType& p : *this)
344  {
345  const point position(p.position());
346  os << "v "
347  << position.x() << ' '
348  << position.y() << ' '
349  << position.z() << nl;
350  }
351 }
352 
353 
354 template<class ParticleType>
356 {
357  // Store the global positions for later use by autoMap. It would be
358  // preferable not to need this. If the mapPolyMesh object passed to autoMap
359  // had a copy of the old mesh then the global positions could be recovered
360  // within autoMap, and this pre-processing would not be necessary.
361 
362  globalPositionsPtr_.reset(new vectorField(this->size()));
363  vectorField& positions = globalPositionsPtr_();
364 
365  label i = 0;
366  for (const ParticleType& p : *this)
367  {
368  positions[i] = p.position();
369  ++i;
370  }
371 }
372 
373 
374 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
375 
376 #include "CloudIO.C"
377 
378 // ************************************************************************* //
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:125
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:598
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:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:297
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:94
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void storeGlobalPositions() const
Call this before a topology change.
Definition: Cloud.C:348
const word cloudName(propsDict.get< word >("cloud"))
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
Cloud(const polyMesh &mesh, const Foam::zero, const word &cloudName)
Construct without particles.
Definition: Cloud.C:57
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:108
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:137
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.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition: Cloud.C:101
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:329
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.