Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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>
42void Foam::Cloud<ParticleType>::checkPatches() const
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
62template<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();
76 (void)polyMesh_.tetBasePtIs();
77 (void)polyMesh_.oldCellCentres();
78}
79
80
81template<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
100template<class ParticleType>
102{
103 this->append(pPtr);
104}
105
106
107template<class ParticleType>
109{
110 delete(this->remove(&p));
111}
112
113
114template<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
126 }
127 }
128}
129
130
131template<class ParticleType>
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
141template<class ParticleType>
142template<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;
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 (
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 (
239 .neighbProcNo()
240 );
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();
252 // Tuple: (patchi particle)
253 (*osptr) << procPatchNeighbours[patchi] << p;
254
255 // Can now remove from my list
257 }
258 }
259 else
260 {
261 deleteParticle(p);
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 }
298 }
299 }
300}
301
302
303template<class ParticleType>
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
335template<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 }
352
353
354template<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// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Base cloud calls templated on particle type.
Definition Cloud.H:64
void writePositions() const
Write positions to <cloudName>_positions.obj file.
Definition Cloud.C:329
void deleteParticle(ParticleType &p)
Remove particle from cloud and delete.
Definition Cloud.C:101
void deleteLostParticles()
Remove lost particles from cloud and delete.
Definition Cloud.C:108
Cloud(const polyMesh &mesh, const Foam::zero, const word &cloudName)
Construct without particles.
Definition Cloud.C:57
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition Cloud.C:137
void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
Definition Cloud.C:297
void cloudReset(const Cloud< ParticleType > &c)
Reset the particles.
Definition Cloud.C:125
const polyMesh & pMesh() const noexcept
Return the polyMesh reference.
Definition Cloud.H:193
cloud::geometryType geometryType_
Geometry type.
Definition Cloud.H:121
void storeGlobalPositions() const
Call this before a topology change.
Definition Cloud.C:348
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition Cloud.C:94
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
void operator=(const ILList< DLListBase, T > &lst)
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
bool eof() const noexcept
True if end of input seen.
Definition IOstream.H:289
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
Buffers for inter-processor communications streams (UOPstream, UIPstream).
bool allowClearRecv() const noexcept
Is clearStorage of individual receive buffer by external hooks allowed? (default: true).
bool hasRecvData() const
True if any (local) recv buffers have unconsumed data. Must call finishedSends() or other finished....
label recvDataCount(const label proci) const
Number of unconsumed receive bytes for the specified processor. Must call finishedSends() or other fi...
void clear()
Clear all send/recv buffers and reset states.
void finishedNeighbourSends(const labelUList &neighProcs, const bool wait=true)
Mark the send phase as being finished, with communication being limited to a known subset of send/rec...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
bool get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition UList.H:868
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
geometryType
Cloud geometry type (internal or IO representations).
Definition cloud.H:63
cloud(const cloud &)=delete
No copy construct.
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
Return processorPatchIndices of the neighbours processor patches. -1 if not running parallel.
const processorTopology & topology() const noexcept
The processor to processor topology.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const Time & time() const noexcept
Return time registry.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const labelList & procNeighbours() const
The neighbour processor connections (ascending order) associated with the local rank.
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#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))))
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
ILList< DLListBase, T > IDLList
Definition IDLList.H:39
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299