parLagrangianDistributorTemplates.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2018-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
30#include "Time.H"
31#include "IOobjectList.H"
33#include "cloud.H"
34#include "CompactIOField.H"
35#include "DynamicList.H"
36#include "passivePositionParticleCloud.H"
37
38// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39
40template<class Container>
42(
43 const IOobjectList& objects,
44 const wordRes& selectedFields
45)
46{
47 wordList fieldNames =
48 (
49 selectedFields.empty()
50 ? objects.names<Container>()
51 : objects.names<Container>(selectedFields)
52 );
53
54 // Parallel synchronise
55 // - Combine names from all processors
56
57 Pstream::combineGather(fieldNames, ListOps::uniqueEqOp<word>());
58 Pstream::broadcast(fieldNames);
59
60 // Sort for consistent order on all processors
61 Foam::sort(fieldNames);
62
63 return fieldNames;
64}
65
66
67template<class Type>
69(
70 const mapDistributeBase& map,
71 const word& cloudName,
72 const IOobjectList& objects,
73 const wordRes& selectedFields
74) const
75{
76 typedef IOField<Type> fieldType;
77
78 const wordList fieldNames
79 (
80 filterObjects<fieldType>
81 (
82 objects,
83 selectedFields
84 )
85 );
86
87
88 bool reconstruct = false;
89
90 label nFields = 0;
91 for (const word& objectName : fieldNames)
92 {
93 if (!nFields)
94 {
95 // Performing an all-to-one (reconstruct)?
97 (
98 (!map.constructSize() || Pstream::master()),
99 andOp<bool>()
100 );
101 }
102
103 if (verbose_)
104 {
105 if (!nFields)
106 {
107 Info<< " Distributing lagrangian "
108 << fieldType::typeName << "s\n" << nl;
109 }
110 Info<< " " << objectName << nl;
111 }
112 ++nFields;
113
114 // Read if present
115 IOField<Type> field
116 (
117 IOobject
118 (
119 objectName,
120 srcMesh_.time().timeName(),
122 srcMesh_,
125 false
126 ),
127 label(0)
128 );
129
130 map.distribute(field);
131
132
133 const IOobject fieldIO
134 (
135 objectName,
136 tgtMesh_.time().timeName(),
138 tgtMesh_,
141 false
142 );
143
144 if (field.size())
145 {
146 IOField<Type>(fieldIO, std::move(field)).write();
147 }
148 else if (!reconstruct)
149 {
150 // When running with -overwrite it should also delete the old
151 // files. Below works but is not optimal.
152
153 const fileName fldName(fieldIO.objectPath());
154 Foam::rm(fldName);
155 }
156 }
157
158 if (nFields && verbose_) Info<< endl;
159 return nFields;
160}
161
162
163template<class Type>
165(
166 const mapDistributeBase& map,
167 const word& cloudName,
168 const IOobjectList& objects,
169 const wordRes& selectedFields
170) const
171{
172 typedef CompactIOField<Field<Type>, Type> fieldType;
173
174 DynamicList<word> fieldNames;
175
176 fieldNames.append
177 (
178 filterObjects<fieldType>
179 (
180 objects,
181 selectedFields
182 )
183 );
184
185 // Append IOField Field names
186 fieldNames.append
187 (
188 filterObjects<IOField<Field<Type>>>
189 (
190 objects,
191 selectedFields
192 )
193 );
194
195 bool reconstruct = false;
196
197 label nFields = 0;
198 for (const word& objectName : fieldNames)
199 {
200 if (!nFields)
201 {
202 // Performing an all-to-one (reconstruct)?
204 (
205 (!map.constructSize() || Pstream::master()),
206 andOp<bool>()
207 );
208 }
209
210 if (verbose_)
211 {
212 if (!nFields)
213 {
214 Info<< " Distributing lagrangian "
215 << fieldType::typeName << "s\n" << nl;
216 }
217 Info<< " " << objectName << nl;
218 }
219 ++nFields;
220
221 // Read if present
222 CompactIOField<Field<Type>, Type> field
223 (
224 IOobject
225 (
226 objectName,
227 srcMesh_.time().timeName(),
229 srcMesh_,
232 false
233 ),
234 label(0)
235 );
236
237 // Distribute
238 map.distribute(field);
239
240 // Write
241 const IOobject fieldIO
242 (
243 objectName,
244 tgtMesh_.time().timeName(),
246 tgtMesh_,
249 false
250 );
251
252 if (field.size())
253 {
254 CompactIOField<Field<Type>, Type>
255 (
256 fieldIO,
257 std::move(field)
258 ).write();
259 }
260 else if (!reconstruct)
261 {
262 // When running with -overwrite it should also delete the old
263 // files. Below works but is not optimal.
264
265 const fileName fldName(fieldIO.objectPath());
266 Foam::rm(fldName);
267 }
268 }
269
270 if (nFields && verbose_) Info<< endl;
271 return nFields;
272}
273
274
275template<class Container>
277(
278 const passivePositionParticleCloud& cloud,
279 const IOobjectList& objects,
280 const wordRes& selectedFields
281)
282{
283 const word fieldClassName(Container::typeName);
284
285 const wordList fieldNames
286 (
287 filterObjects<Container>
288 (
289 objects,
290 selectedFields
291 )
292 );
293
294 label nFields = 0;
295 for (const word& objectName : fieldNames)
296 {
297 if (verbose_)
298 {
299 if (!nFields)
300 {
301 Info<< " Reading lagrangian "
302 << Container::typeName << "s\n" << nl;
303 }
304 Info<< " " << objectName << nl;
305 }
306 ++nFields;
307
308 // Read if present
309 Container* fieldPtr = new Container
310 (
311 IOobject
312 (
313 objectName,
314 cloud.time().timeName(),
315 cloud,
318 ),
319 label(0)
320 );
321
322 fieldPtr->store();
323 }
324
325 if (nFields && verbose_) Info<< endl;
326 return nFields;
327}
328
329
330template<class Container>
332(
333 const mapDistributeBase& map,
334 passivePositionParticleCloud& cloud
335) const
336{
337 HashTable<Container*> fields
338 (
339 cloud.lookupClass<Container>()
340 );
341
342 bool reconstruct = false;
343
344 label nFields = 0;
345 forAllIters(fields, iter)
346 {
347 Container& field = *(iter.val());
348
349 if (!nFields)
350 {
351 // Performing an all-to-one (reconstruct)?
353 (
354 (!map.constructSize() || Pstream::master()),
355 andOp<bool>()
356 );
357 }
358
359 if (verbose_)
360 {
361 if (!nFields)
362 {
363 Info<< " Distributing lagrangian "
364 << Container::typeName << "s\n" << nl;
365 }
366 Info<< " " << field.name() << nl;
367 }
368 ++nFields;
369
370 map.distribute(field);
371
372 const IOobject fieldIO
373 (
374 field.name(),
375 tgtMesh_.time().timeName(),
376 cloud::prefix/cloud.name(),
377 tgtMesh_,
380 false
381 );
382
383 if (field.size())
384 {
385 Container(fieldIO, std::move(field)).write();
386 }
387 else if (!reconstruct)
388 {
389 // When running with -overwrite it should also delete the old
390 // files. Below works but is not optimal.
391
392 const fileName fldName(fieldIO.objectPath());
393 Foam::rm(fldName);
394 }
395 }
396
397 if (nFields && verbose_) Info<< endl;
398 return nFields;
399}
400
401
402// ************************************************************************* //
static void combineGather(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag, const label comm)
static void broadcast(Type &value, const label comm=UPstream::worldComm)
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
static wordList filterObjects(const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Pick up any fields of a given type.
label distributeStoredFields(const mapDistributeBase &map, passivePositionParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
label distributeFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fieldFields.
label distributeFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fields.
splitCell * master() const
Definition: splitCell.H:113
rDeltaTY field()
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: MSwindows.C:1012
List< word > wordList
A List of words.
Definition: fileName.H:63
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:342
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
runTime write()
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
const word cloudName(propsDict.get< word >("cloud"))