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-2021 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  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),
81  polyMesh_(pMesh),
82  labels_(),
83  globalPositionsPtr_(),
84  geometryType_(cloud::geometryType::COORDINATES)
85 {
86  checkPatches();
87  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  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>
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[Pstream::myProcNo()];
165 
166  // Indexing from the processor number into the neighbourProcs list
167  labelList neighbourProcIndices(Pstream::nProcs(), -1);
168 
169  forAll(neighbourProcs, i)
170  {
171  neighbourProcIndices[neighbourProcs[i]] = i;
172  }
173 
174  // Initialise the stepFraction moved for the particles
175  forAllIters(*this, pIter)
176  {
177  pIter().reset();
178  }
179 
180  // List of lists of particles to be transferred for all of the
181  // neighbour processors
182  List<IDLList<ParticleType>> particleTransferLists
183  (
184  neighbourProcs.size()
185  );
186 
187  // List of destination processorPatches indices for all of the
188  // neighbour processors
189  List<DynamicList<label>> patchIndexTransferLists
190  (
191  neighbourProcs.size()
192  );
193 
194  // Allocate transfer buffers
195  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
196 
197  // Clear the global positions as there are about to change
198  globalPositionsPtr_.clear();
199 
200  // While there are particles to transfer
201  while (true)
202  {
203  particleTransferLists = IDLList<ParticleType>();
204  forAll(patchIndexTransferLists, i)
205  {
206  patchIndexTransferLists[i].clear();
207  }
208 
209  // Loop over all particles
210  for (ParticleType& p : *this)
211  {
212  // Move the particle
213  bool keepParticle = p.move(cloud, td, trackTime);
214 
215  // If the particle is to be kept
216  // (i.e. it hasn't passed through an inlet or outlet)
217  if (keepParticle)
218  {
219  if (td.switchProcessor)
220  {
221  #ifdef FULLDEBUG
222  if
223  (
224  !Pstream::parRun()
225  || !p.onBoundaryFace()
226  || procPatchNeighbours[p.patch()] < 0
227  )
228  {
230  << "Switch processor flag is true when no parallel "
231  << "transfer is possible. This is a bug."
232  << exit(FatalError);
233  }
234  #endif
235 
236  const label patchi = p.patch();
237 
238  const label n = neighbourProcIndices
239  [
240  refCast<const processorPolyPatch>
241  (
242  pbm[patchi]
243  ).neighbProcNo()
244  ];
245 
246  p.prepareForParallelTransfer();
247 
248  particleTransferLists[n].append(this->remove(&p));
249 
250  patchIndexTransferLists[n].append
251  (
252  procPatchNeighbours[patchi]
253  );
254  }
255  }
256  else
257  {
258  deleteParticle(p);
259  }
260  }
261 
262  if (!Pstream::parRun())
263  {
264  break;
265  }
266 
267 
268  // Clear transfer buffers
269  pBufs.clear();
270 
271  // Stream into send buffers
272  forAll(particleTransferLists, i)
273  {
274  if (particleTransferLists[i].size())
275  {
276  UOPstream particleStream
277  (
278  neighbourProcs[i],
279  pBufs
280  );
281 
282  particleStream
283  << patchIndexTransferLists[i]
284  << particleTransferLists[i];
285  }
286  }
287 
288 
289  // Start sending. Sets number of bytes transferred
290  labelList allNTrans(Pstream::nProcs());
291  pBufs.finishedSends(allNTrans);
292 
293 
294  bool transferred = false;
295 
296  for (const label n : allNTrans)
297  {
298  if (n)
299  {
300  transferred = true;
301  break;
302  }
303  }
304  reduce(transferred, orOp<bool>());
305 
306  if (!transferred)
307  {
308  break;
309  }
310 
311  // Retrieve from receive buffers
312  for (const label neighbProci : neighbourProcs)
313  {
314  label nRec = allNTrans[neighbProci];
315 
316  if (nRec)
317  {
318  UIPstream particleStream(neighbProci, pBufs);
319 
320  labelList receivePatchIndex(particleStream);
321 
322  IDLList<ParticleType> newParticles
323  (
324  particleStream,
325  typename ParticleType::iNew(polyMesh_)
326  );
327 
328  label pI = 0;
329 
330  for (ParticleType& newp : newParticles)
331  {
332  label patchi = procPatches[receivePatchIndex[pI++]];
333 
334  newp.correctAfterParallelTransfer(patchi, td);
335 
336  addParticle(newParticles.remove(&newp));
337  }
338  }
339  }
340  }
341 }
342 
343 
344 template<class ParticleType>
346 {
347  if (!globalPositionsPtr_)
348  {
350  << "Global positions are not available. "
351  << "Cloud::storeGlobalPositions has not been called."
352  << exit(FatalError);
353  }
354 
355  // Reset stored data that relies on the mesh
356  // polyMesh_.clearCellTree();
357  cellWallFacesPtr_.clear();
358 
359  // Ask for the tetBasePtIs to trigger all processors to build
360  // them, otherwise, if some processors have no particles then
361  // there is a comms mismatch.
362  polyMesh_.tetBasePtIs();
363  polyMesh_.oldCellCentres();
364 
365  const vectorField& positions = globalPositionsPtr_();
366 
367  label i = 0;
368  forAllIters(*this, iter)
369  {
370  iter().autoMap(positions[i], mapper);
371  ++i;
372  }
373 }
374 
375 
376 template<class ParticleType>
378 {
379  OFstream pObj
380  (
381  this->db().time().path()/this->name() + "_positions.obj"
382  );
383 
384  forAllConstIters(*this, pIter)
385  {
386  const ParticleType& p = pIter();
387  const point position(p.position());
388  pObj<< "v " << position.x() << " " << position.y() << " "
389  << position.z() << nl;
390  }
391 
392  pObj.flush();
393 }
394 
395 
396 template<class ParticleType>
398 {
399  // Store the global positions for later use by autoMap. It would be
400  // preferable not to need this. If the mapPolyMesh object passed to autoMap
401  // had a copy of the old mesh then the global positions could be recovered
402  // within autoMap, and this pre-processing would not be necessary.
403 
404  globalPositionsPtr_.reset(new vectorField(this->size()));
405 
406  vectorField& positions = globalPositionsPtr_();
407 
408  label i = 0;
409  forAllConstIters(*this, iter)
410  {
411  positions[i] = iter().position();
412  ++i;
413  }
414 }
415 
416 
417 // * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
418 
419 #include "CloudIO.C"
420 
421 // ************************************************************************* //
Foam::Cloud::addParticle
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:104
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:65
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:63
globalMeshData.H
Foam::Cloud::writePositions
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:377
mapPolyMesh.H
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
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:369
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
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
Definition: globalMeshData.H:381
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::Cloud::storeGlobalPositions
void storeGlobalPositions() const
Call this before a topology change.
Definition: Cloud.C:397
Foam::Cloud::deleteParticle
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition: Cloud.C:111
Foam::Cloud::autoMap
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:345
Foam::Cloud::deleteLostParticles
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:118
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:147
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
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:135
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:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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:55
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
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:73
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::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
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:328
Foam::globalMeshData::processorPatchNeighbours
const labelList & processorPatchNeighbours() const noexcept
Return processorPatchIndices of the neighbours.
Definition: globalMeshData.H:396
Foam::OSstream::flush
virtual void flush()
Flush stream.
Definition: OSstream.C:275