parLagrangianRedistributorFields.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 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"
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 template<class Container>
41 (
42  const IOobjectList& objects,
43  const wordRes& selectedFields
44 )
45 {
47  (
48  selectedFields.empty()
49  ? objects.names<Container>()
50  : objects.names<Container>(selectedFields)
51  );
52 
53  // Parallel synchronise
54  // - Combine names from all processors
55 
56  Pstream::combineGather(fieldNames, ListOps::uniqueEqOp<word>());
58 
59  // Sort for consistent order on all processors
61 
62  return fieldNames;
63 }
64 
65 
66 template<class Type>
68 (
69  const mapDistributeBase& map,
70  const word& cloudName,
71  const IOobjectList& objects,
72  const wordRes& selectedFields
73 ) const
74 {
75  typedef IOField<Type> fieldType;
76 
77  const wordList fieldNames
78  (
79  filterObjects<IOField<Type>>
80  (
81  objects,
82  selectedFields
83  )
84  );
85 
86  label nFields = 0;
87  for (const word& objectName : fieldNames)
88  {
89  if (!nFields)
90  {
91  Info<< " Redistributing lagrangian "
92  << fieldType::typeName << "s\n" << endl;
93  }
94  Info<< " " << objectName << endl;
95 
96  // Read if present
97  IOField<Type> field
98  (
99  IOobject
100  (
101  objectName,
102  srcMesh_.time().timeName(),
104  srcMesh_,
107  false
108  ),
109  label(0)
110  );
111 
112  map.distribute(field);
113 
114 
115  const IOobject fieldIO
116  (
117  objectName,
118  tgtMesh_.time().timeName(),
120  tgtMesh_,
123  false
124  );
125 
126  if (field.size())
127  {
128  IOField<Type>
129  (
130  fieldIO,
131  std::move(field)
132  ).write();
133  }
134  else
135  {
136  // When running with -overwrite it should also delete the old
137  // files. Below works but is not optimal.
138 
139  const fileName fldName(fieldIO.objectPath());
140  Foam::rm(fldName);
141  }
142  }
143 
144  return nFields;
145 }
146 
147 
148 template<class Type>
150 (
151  const mapDistributeBase& map,
152  const word& cloudName,
153  const IOobjectList& objects,
154  const wordRes& selectedFields
155 ) const
156 {
157  typedef CompactIOField<Field<Type>, Type> fieldType;
158 
160  (
161  filterObjects<fieldType>
162  (
163  objects,
164  selectedFields
165  )
166  );
167 
168  // Append IOField Field names
169  {
170  wordList ioFieldNames
171  (
172  filterObjects<IOField<Field<Type>>>
173  (
174  objects,
175  selectedFields
176  )
177  );
178  fieldNames.append(ioFieldNames);
179  }
180 
181  label nFields = 0;
182  for (const word& objectName : fieldNames)
183  {
184  if (!nFields++)
185  {
186  Info<< " Redistributing lagrangian "
187  << fieldType::typeName << "s\n" << nl;
188  }
189  Info<< " " << objectName << nl;
190 
191  // Read if present
192  CompactIOField<Field<Type>, Type> field
193  (
194  IOobject
195  (
196  objectName,
197  srcMesh_.time().timeName(),
199  srcMesh_,
202  false
203  ),
204  label(0)
205  );
206 
207  // Distribute
208  map.distribute(field);
209 
210  // Write
211  const IOobject fieldIO
212  (
213  objectName,
214  tgtMesh_.time().timeName(),
216  tgtMesh_,
219  false
220  );
221 
222  if (field.size())
223  {
224  CompactIOField<Field<Type>, Type>
225  (
226  fieldIO,
227  std::move(field)
228  ).write();
229  }
230  else
231  {
232  // When running with -overwrite it should also delete the old
233  // files. Below works but is not optimal.
234 
235  const fileName fldName(fieldIO.objectPath());
236  Foam::rm(fldName);
237  }
238  }
239 
240  if (nFields) Info<< endl;
241  return nFields;
242 }
243 
244 
245 template<class Container>
247 (
248  const passivePositionParticleCloud& cloud,
249  const IOobjectList& objects,
250  const wordRes& selectedFields
251 )
252 {
253  const word fieldClassName(Container::typeName);
254 
255  const wordList fieldNames
256  (
257  filterObjects<Container>
258  (
259  objects,
260  selectedFields
261  )
262  );
263 
264  label nFields = 0;
265  for (const word& objectName : fieldNames)
266  {
267  if (!nFields++)
268  {
269  Info<< " Reading lagrangian "
270  << Container::typeName << "s\n" << nl;
271  }
272  Info<< " " << objectName << nl;
273 
274  // Read if present
275  Container* fieldPtr = new Container
276  (
277  IOobject
278  (
279  objectName,
280  cloud.time().timeName(),
281  cloud,
284  ),
285  label(0)
286  );
287 
288  fieldPtr->store();
289  }
290 
291  return nFields;
292 }
293 
294 
295 template<class Container>
297 (
298  const mapDistributeBase& map,
299  passivePositionParticleCloud& cloud
300 ) const
301 {
302  HashTable<Container*> fields
303  (
304  cloud.lookupClass<Container>()
305  );
306 
307  label nFields = 0;
308  forAllIters(fields, iter)
309  {
310  Container& field = *(iter.val());
311 
312  if (!nFields++)
313  {
314  Info<< " Redistributing lagrangian "
315  << Container::typeName << "s\n" << endl;
316  }
317  Info<< " " << field.name() << endl;
318 
319  map.distribute(field);
320 
321  const IOobject fieldIO
322  (
323  field.name(),
324  tgtMesh_.time().timeName(),
325  cloud::prefix/cloud.name(),
326  tgtMesh_,
329  false
330  );
331 
332  if (field.size())
333  {
334  Container
335  (
336  fieldIO,
337  std::move(field)
338  ).write();
339  }
340  else
341  {
342  // When running with -overwrite it should also delete the old
343  // files. Below works but is not optimal.
344 
345  const fileName fldName(fieldIO.objectPath());
346  Foam::rm(fldName);
347  }
348  }
349 
350  return nFields;
351 }
352 
353 
354 // ************************************************************************* //
Foam::parLagrangianRedistributor::redistributeFields
label redistributeFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fields.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
cloudName
const word cloudName(propsDict.get< word >("cloud"))
Foam::parLagrangianRedistributor::redistributeStoredFields
label redistributeStoredFields(const mapDistributeBase &map, passivePositionParticleCloud &cloud) const
Redistribute and write stored lagrangian fields.
cloud.H
Foam::rm
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: MSwindows.C:1004
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::parLagrangianRedistributor::redistributeFieldFields
label redistributeFieldFields(const mapDistributeBase &map, const word &cloudName, const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected lagrangian fieldFields.
parLagrangianRedistributor.H
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
mapDistributePolyMesh.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::parLagrangianRedistributor::filterObjects
static wordList filterObjects(const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Pick up any fields of a given type.
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
field
rDeltaTY field()
passivePositionParticleCloud.H
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::Pstream::combineGather
static void combineGather(const List< commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:48
Time.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::parLagrangianRedistributor::readFields
static label readFields(const passivePositionParticleCloud &cloud, const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Read and store all fields of a cloud.
Foam::List< word >
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::IOobject::NO_READ
Definition: IOobject.H:188
CompactIOField.H
Foam::Pstream::combineScatter
static void combineScatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:183
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97