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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
41template<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
71template<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 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
103template<class ParticleType>
105{
106 this->append(pPtr);
107}
108
109
110template<class ParticleType>
112{
113 delete(this->remove(&p));
114}
115
116
117template<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
134template<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
144template<class ParticleType>
145template<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 =
166
167 // Initialise the stepFraction moved for the particles
168 for (ParticleType& p : *this)
169 {
170 p.reset();
171 }
172
173 // Clear the global positions as these are about to change
174 globalPositionsPtr_.clear();
175
176
177 // For v2112 and earlier: pre-assembled lists of particles
178 // to be transferred and target patch on a per processor basis.
179 // Apart from memory overhead of assembling the lists this adds
180 // allocations/de-allocation when building linked-lists.
181
182 // Now stream particle transfer tuples directly into PstreamBuffers.
183 // Use a local cache of UOPstream wrappers for the formatters
184 // (since there are potentially many particles being shifted about).
185
186
187 // Allocate transfer buffers,
188 // automatic clearStorage when UIPstream closes is disabled.
190 pBufs.allowClearRecv(false);
191
192 // Cache of opened UOPstream wrappers
194
195 // While there are particles to transfer
196 while (true)
197 {
198 // Reset transfer buffers
199 pBufs.clear();
200
201 // Rewind existing streams
202 forAll(UOPstreamPtrs, proci)
203 {
204 auto* osptr = UOPstreamPtrs.get(proci);
205 if (osptr)
207 osptr->rewind();
208 }
209 }
210
211 // Loop over all particles
212 for (ParticleType& p : *this)
213 {
214 // Move the particle
215 const bool keepParticle = p.move(cloud, td, trackTime);
216
217 // If the particle is to be kept
218 // (i.e. it hasn't passed through an inlet or outlet)
219 if (keepParticle)
220 {
221 if (td.switchProcessor)
222 {
223 #ifdef FULLDEBUG
224 if
225 (
227 || !p.onBoundaryFace()
228 || procPatchNeighbours[p.patch()] < 0
229 )
230 {
232 << "Switch processor flag is true when no parallel "
233 << "transfer is possible. This is a bug."
234 << exit(FatalError);
235 }
236 #endif
237
238 const label patchi = p.patch();
239
240 const label toProci =
241 (
242 refCast<const processorPolyPatch>(pbm[patchi])
243 .neighbProcNo()
244 );
245
246 // Get/create output stream
247 auto* osptr = UOPstreamPtrs.get(toProci);
248 if (!osptr)
249 {
250 osptr = new UOPstream(toProci, pBufs);
251 UOPstreamPtrs.set(toProci, osptr);
252 }
253
254 p.prepareForParallelTransfer();
255
256 // Tuple: (patchi particle)
257 (*osptr) << procPatchNeighbours[patchi] << p;
258
259 // Can now remove from my list
260 deleteParticle(p);
261 }
262 }
263 else
264 {
265 deleteParticle(p);
266 }
267 }
268
269 if (!Pstream::parRun())
270 {
271 break;
272 }
273
274 pBufs.finishedNeighbourSends(neighbourProcs);
276 if (!returnReduce(pBufs.hasRecvData(), orOp<bool>()))
277 {
278 // No parcels to transfer
279 break;
281
282 // Retrieve from receive buffers
283 for (const label proci : neighbourProcs)
284 {
285 if (pBufs.recvDataCount(proci))
286 {
287 UIPstream is(proci, pBufs);
288
289 // Read out each (patchi particle) tuple
290 while (!is.eof())
291 {
292 label patchi = pTraits<label>(is);
293 auto* newp = new ParticleType(polyMesh_, is);
294
295 // The real patch index
296 patchi = procPatches[patchi];
297
298 (*newp).correctAfterParallelTransfer(patchi, td);
299 addParticle(newp);
300 }
301 }
302 }
303 }
304}
305
306
307template<class ParticleType>
309{
310 if (!globalPositionsPtr_)
311 {
313 << "Global positions are not available. "
314 << "Cloud::storeGlobalPositions has not been called."
315 << exit(FatalError);
316 }
317
318 // Reset stored data that relies on the mesh
319 // polyMesh_.clearCellTree();
320 cellWallFacesPtr_.clear();
321
322 // Ask for the tetBasePtIs to trigger all processors to build
323 // them, otherwise, if some processors have no particles then
324 // there is a comms mismatch.
325 polyMesh_.tetBasePtIs();
326 polyMesh_.oldCellCentres();
327
328 const vectorField& positions = globalPositionsPtr_();
329
330 label i = 0;
331 for (ParticleType& p : *this)
332 {
333 p.autoMap(positions[i], mapper);
334 ++i;
335 }
336}
337
338
339template<class ParticleType>
341{
343 (
344 this->db().time().path()/this->name() + "_positions.obj"
345 );
346
347 for (const ParticleType& p : *this)
348 {
349 const point position(p.position());
350 os << "v "
351 << position.x() << ' '
352 << position.y() << ' '
353 << position.z() << nl;
354 }
355}
356
357
358template<class ParticleType>
360{
361 // Store the global positions for later use by autoMap. It would be
362 // preferable not to need this. If the mapPolyMesh object passed to autoMap
363 // had a copy of the old mesh then the global positions could be recovered
364 // within autoMap, and this pre-processing would not be necessary.
365
366 globalPositionsPtr_.reset(new vectorField(this->size()));
367 vectorField& positions = globalPositionsPtr_();
368
369 label i = 0;
370 for (const ParticleType& p : *this)
371 {
372 positions[i] = p.position();
373 ++i;
374 }
375}
376
377
378// * * * * * * * * * * * * * * * * IOStream operators * * * * * * * * * * * //
379
380#include "CloudIO.C"
381
382// ************************************************************************* //
Base cloud calls templated on particle type.
Definition: Cloud.H:68
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition: Cloud.C:340
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition: Cloud.C:111
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition: Cloud.C:118
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition: Cloud.C:308
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition: Cloud.C:135
void storeGlobalPositions() const
Call this before a topology change.
Definition: Cloud.C:359
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:104
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:403
Template class for intrusive linked lists.
Definition: ILList.H:69
void operator=(const ILList< LListBase, T > &lst)
Copy assignment using the 'clone()' method for each element.
Definition: ILList.C:125
bool eof() const noexcept
True if end of input seen.
Definition: IOstream.H:239
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Output to file stream, using an OSstream.
Definition: OFstream.H:57
Buffers for inter-processor communications streams (UOPstream, UIPstream).
bool allowClearRecv() const noexcept
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
bool hasRecvData() const
label recvDataCount(const label proci) const
void clear()
Clear individual buffers and reset states.
void finishedNeighbourSends(const labelUList &neighProcs, const bool wait=true)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get(const label i) const
Definition: UList.H:528
@ nonBlocking
"nonBlocking"
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
geometryType
Cloud geometry type (internal or IO representations)
Definition: cloud.H:75
virtual void move()=0
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.
const labelList & processorPatchNeighbours() const noexcept
const processorTopology & topology() const noexcept
The processor to processor topology.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:906
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion)
Definition: polyMesh.C:1155
int myProcNo() const noexcept
Return processor number.
const labelListList & procNeighbours() const noexcept
A class for handling words, derived from Foam::string.
Definition: word.H:68
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const word cloudName(propsDict.get< word >("cloud"))