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 -------------------------------------------------------------------------------
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 
30 #include "Time.H"
31 #include "IOobjectList.H"
32 #include "mapDistributePolyMesh.H"
33 #include "cloud.H"
34 #include "CompactIOField.H"
35 #include "DynamicList.H"
36 #include "passivePositionParticleCloud.H"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 template<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 - combine names from all processors
55  Pstream::combineReduce(fieldNames, ListOps::uniqueEqOp<word>());
56  Foam::sort(fieldNames); // Consistent order
57 
58  return fieldNames;
59 }
60 
61 
62 template<class Type>
64 (
65  const mapDistributeBase& map,
66  const word& cloudName,
67  const IOobjectList& objects,
68  const wordRes& selectedFields
69 ) const
70 {
71  typedef IOField<Type> fieldType;
72 
73  const wordList fieldNames
74  (
75  filterObjects<fieldType>
76  (
77  objects,
78  selectedFields
79  )
80  );
81 
82 
83  bool reconstruct = false;
84 
85  label nFields = 0;
86  for (const word& objectName : fieldNames)
87  {
88  if (!nFields)
89  {
90  // Performing an all-to-one (reconstruct)?
91  reconstruct =
92  returnReduceAnd(!map.constructSize() || Pstream::master());
93  }
94 
95  if (verbose_)
96  {
97  if (!nFields)
98  {
99  Info<< " Distributing lagrangian "
100  << fieldType::typeName << "s\n" << nl;
101  }
102  Info<< " " << objectName << nl;
103  }
104  ++nFields;
105 
106  // Read if present
107  IOField<Type> field
108  (
109  IOobject
110  (
111  objectName,
112  srcMesh_.time().timeName(),
114  srcMesh_,
118  ),
119  label(0)
120  );
121 
122  map.distribute(field);
123 
124 
125  const IOobject fieldIO
126  (
127  objectName,
128  tgtMesh_.time().timeName(),
130  tgtMesh_,
134  );
135 
136  if (field.size())
137  {
138  IOField<Type>(fieldIO, std::move(field)).write();
139  }
140  else if (!reconstruct)
141  {
142  // When running with -overwrite it should also delete the old
143  // files. Below works but is not optimal.
144 
145  const fileName fldName(fieldIO.objectPath());
146  Foam::rm(fldName);
147  }
148  }
149 
150  if (nFields && verbose_) Info<< endl;
151  return nFields;
152 }
153 
154 
155 template<class Type>
157 (
158  const mapDistributeBase& map,
159  const word& cloudName,
160  const IOobjectList& objects,
161  const wordRes& selectedFields
162 ) const
163 {
164  typedef CompactIOField<Field<Type>, Type> fieldType;
165 
166  DynamicList<word> fieldNames;
167 
168  fieldNames.append
169  (
170  filterObjects<fieldType>
171  (
172  objects,
173  selectedFields
174  )
175  );
176 
177  // Append IOField Field names
178  fieldNames.append
179  (
180  filterObjects<IOField<Field<Type>>>
181  (
182  objects,
183  selectedFields
184  )
185  );
186 
187  bool reconstruct = false;
188 
189  label nFields = 0;
190  for (const word& objectName : fieldNames)
191  {
192  if (!nFields)
193  {
194  // Performing an all-to-one (reconstruct)?
195  reconstruct =
196  returnReduceAnd(!map.constructSize() || Pstream::master());
197  }
198 
199  if (verbose_)
200  {
201  if (!nFields)
202  {
203  Info<< " Distributing lagrangian "
204  << fieldType::typeName << "s\n" << nl;
205  }
206  Info<< " " << objectName << nl;
207  }
208  ++nFields;
209 
210  // Read if present
211  CompactIOField<Field<Type>, Type> field
212  (
213  IOobject
214  (
215  objectName,
216  srcMesh_.time().timeName(),
218  srcMesh_,
222  ),
223  label(0)
224  );
225 
226  // Distribute
227  map.distribute(field);
228 
229  // Write
230  const IOobject fieldIO
231  (
232  objectName,
233  tgtMesh_.time().timeName(),
235  tgtMesh_,
239  );
240 
241  if (field.size())
242  {
243  CompactIOField<Field<Type>, Type>
244  (
245  fieldIO,
246  std::move(field)
247  ).write();
248  }
249  else if (!reconstruct)
250  {
251  // When running with -overwrite it should also delete the old
252  // files. Below works but is not optimal.
253 
254  const fileName fldName(fieldIO.objectPath());
255  Foam::rm(fldName);
256  }
257  }
258 
259  if (nFields && verbose_) Info<< endl;
260  return nFields;
261 }
262 
263 
264 template<class Container>
266 (
267  const passivePositionParticleCloud& cloud,
268  const IOobjectList& objects,
269  const wordRes& selectedFields
270 )
271 {
272  const word fieldClassName(Container::typeName);
273 
274  const wordList fieldNames
275  (
276  filterObjects<Container>
277  (
278  objects,
279  selectedFields
280  )
281  );
282 
283  label nFields = 0;
284  for (const word& objectName : fieldNames)
285  {
286  if (verbose_)
287  {
288  if (!nFields)
289  {
290  Info<< " Reading lagrangian "
291  << Container::typeName << "s\n" << nl;
292  }
293  Info<< " " << objectName << nl;
294  }
295  ++nFields;
296 
297  // Read if present
298  Container* fieldPtr = new Container
299  (
300  IOobject
301  (
302  objectName,
303  cloud.time().timeName(),
304  cloud,
308  ),
309  label(0)
310  );
311 
312  fieldPtr->store();
313  }
314 
315  if (nFields && verbose_) Info<< endl;
316  return nFields;
317 }
318 
319 
320 template<class Container>
322 (
323  const mapDistributeBase& map,
324  passivePositionParticleCloud& cloud
325 ) const
326 {
327  HashTable<Container*> fields
328  (
329  cloud.lookupClass<Container>()
330  );
331 
332  bool reconstruct = false;
333 
334  label nFields = 0;
335  forAllIters(fields, iter)
336  {
337  Container& field = *(iter.val());
338 
339  if (!nFields)
340  {
341  // Performing an all-to-one (reconstruct)?
342  reconstruct =
343  returnReduceAnd(!map.constructSize() || Pstream::master());
344  }
345 
346  if (verbose_)
347  {
348  if (!nFields)
349  {
350  Info<< " Distributing lagrangian "
351  << Container::typeName << "s\n" << nl;
352  }
353  Info<< " " << field.name() << nl;
354  }
355  ++nFields;
356 
357  map.distribute(field);
358 
359  const IOobject fieldIO
360  (
361  field.name(),
362  tgtMesh_.time().timeName(),
363  cloud::prefix/cloud.name(),
364  tgtMesh_,
368  );
369 
370  if (field.size())
371  {
372  Container(fieldIO, std::move(field)).write();
373  }
374  else if (!reconstruct)
375  {
376  // When running with -overwrite it should also delete the old
377  // files. Below works but is not optimal.
378 
379  const fileName fldName(fieldIO.objectPath());
380  Foam::rm(fldName);
381  }
382  }
383 
384  if (nFields && verbose_) Info<< endl;
385  return nFields;
386 }
387 
388 
389 // ************************************************************************* //
rDeltaTY field()
label distributeFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fields.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
Ignore writing from objectRegistry::writeObject()
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const word cloudName(propsDict.get< word >("cloud"))
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:348
Reading is optional [identical to LAZY_READ].
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:329
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
label distributeFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fieldFields.
static wordList filterObjects(const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Pick up any fields of a given type.
List< word > wordList
List of word.
Definition: fileName.H:58
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1037
static label readFields(const passivePositionParticleCloud &cloud, const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Read and store all fields of a cloud.
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Reading is optional [identical to READ_IF_PRESENT].
Request registration (bool: true)
Do not request registration (bool: false)
label distributeStoredFields(const mapDistributeBase &map, passivePositionParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: POSIX.C:1404
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:93