parFvFieldDistributorTemplates.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) 2016-2023 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 
29 #include "parFvFieldDistributor.H"
30 #include "Time.H"
31 #include "PtrList.H"
32 #include "fvPatchFields.H"
33 #include "emptyFvPatch.H"
34 #include "emptyFvPatchField.H"
35 #include "emptyFvsPatchField.H"
36 #include "IOobjectList.H"
37 #include "mapDistributePolyMesh.H"
38 #include "processorFvPatch.H"
39 
40 #include "distributedFieldMapper.H"
42 
43 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
44 
45 template<class Type>
48 (
49  const DimensionedField<Type, volMesh>& fld
50 ) const
51 {
52  // Create internalField by remote mapping
53 
55  (
57  distMap_.cellMap()
58  );
59 
60  auto tfield = tmp<DimensionedField<Type, volMesh>>::New
61  (
62  IOobject
63  (
64  fld.name(),
65  tgtMesh_.time().timeName(),
66  fld.local(),
67  tgtMesh_,
70  ),
71  tgtMesh_,
72  fld.dimensions(),
73  Field<Type>(fld, mapper)
74  );
75 
76  tfield.ref().oriented() = fld.oriented();
77 
78  return tfield;
79 }
80 
81 
82 template<class Type>
85 (
86  const GeometricField<Type, fvPatchField, volMesh>& fld
87 ) const
88 {
89  // Create internalField by remote mapping
91  (
93  distMap_.cellMap()
94  );
95 
96  DimensionedField<Type, volMesh> internalField
97  (
98  IOobject
99  (
100  fld.name(),
101  tgtMesh_.time().timeName(),
102  fld.local(),
103  tgtMesh_,
106  ),
107  tgtMesh_,
108  fld.dimensions(),
109  Field<Type>(fld.internalField(), mapper)
110  );
111 
112  internalField.oriented() = fld.oriented();
113 
114 
115  // Create patchFields by remote mapping
116  // Note: patchFields still on source mesh, not target mesh
117 
118  PtrList<fvPatchField<Type>> oldPatchFields(fld.mesh().boundary().size());
119 
120  const auto& bfld = fld.boundaryField();
121 
122  forAll(bfld, patchi)
123  {
124  if (patchFaceMaps_.set(patchi))
125  {
126  // Clone local patch field
127  oldPatchFields.set(patchi, bfld[patchi].clone());
128 
130  (
132  patchFaceMaps_[patchi]
133  );
134 
135  // Map into local copy
136  oldPatchFields[patchi].autoMap(mapper);
137  }
138  }
139 
140 
141  // Clone the oldPatchFields onto the target patches. This is just to reset
142  // the reference to the patch, size and content stay the same.
143 
144  PtrList<fvPatchField<Type>> newPatchFields(tgtMesh_.boundary().size());
145 
146  forAll(oldPatchFields, patchi)
147  {
148  if (oldPatchFields.set(patchi))
149  {
150  const auto& pfld = oldPatchFields[patchi];
151 
152  labelList dummyMap(identity(pfld.size()));
153  directFvPatchFieldMapper dummyMapper(dummyMap);
154 
155  newPatchFields.set
156  (
157  patchi,
159  (
160  pfld,
161  tgtMesh_.boundary()[patchi],
163  dummyMapper
164  )
165  );
166  }
167  }
168 
169  // Add some empty patches on remaining patches
170  // (... probably processor patches)
171 
172  forAll(newPatchFields, patchi)
173  {
174  if (!newPatchFields.set(patchi))
175  {
176  newPatchFields.set
177  (
178  patchi,
180  (
182  tgtMesh_.boundary()[patchi],
184  )
185  );
186  }
187  }
188 
189  // Return geometric field
190 
191  return tmp<GeometricField<Type, fvPatchField, volMesh>>::New
192  (
193  std::move(internalField),
194  newPatchFields
195  );
196 }
197 
198 
199 template<class Type>
202 (
203  const GeometricField<Type, fvsPatchField, surfaceMesh>& fld
204 ) const
205 {
206  // Create internalField by remote mapping
208  (
210  distMap_.faceMap()
211  );
212 
213 
214  Field<Type> primitiveField;
215  {
216  // Create flat field of internalField + all patch fields
217  Field<Type> flatFld(fld.mesh().nFaces(), Type(Zero));
218  SubList<Type>(flatFld, fld.internalField().size())
219  = fld.internalField();
220 
221  for (const fvsPatchField<Type>& fvp : fld.boundaryField())
222  {
223  SubList<Type>(flatFld, fvp.size(), fvp.patch().start()) = fvp;
224  }
225 
226  // Map all faces
227  primitiveField = Field<Type>(flatFld, mapper, fld.is_oriented());
228 
229  // Trim to internal faces (note: could also have special mapper)
230  primitiveField.resize
231  (
232  min
233  (
234  primitiveField.size(),
235  tgtMesh_.nInternalFaces()
236  )
237  );
238  }
239 
240 
241  DimensionedField<Type, surfaceMesh> internalField
242  (
243  IOobject
244  (
245  fld.name(),
246  tgtMesh_.time().timeName(),
247  fld.local(),
248  tgtMesh_,
251  ),
252  tgtMesh_,
253  fld.dimensions(),
254  std::move(primitiveField)
255  );
256 
257  internalField.oriented() = fld.oriented();
258 
259 
260  // Create patchFields by remote mapping
261  // Note: patchFields still on source mesh, not target mesh
262 
263  PtrList<fvsPatchField<Type>> oldPatchFields(fld.mesh().boundary().size());
264 
265  const auto& bfld = fld.boundaryField();
266 
267  forAll(bfld, patchi)
268  {
269  if (patchFaceMaps_.set(patchi))
270  {
271  // Clone local patch field
272  oldPatchFields.set(patchi, bfld[patchi].clone());
273 
275  (
277  patchFaceMaps_[patchi]
278  );
279 
280  // Map into local copy
281  oldPatchFields[patchi].autoMap(mapper);
282  }
283  }
284 
285 
286  PtrList<fvsPatchField<Type>> newPatchFields(tgtMesh_.boundary().size());
287 
288  // Clone the patchFields onto the base patches. This is just to reset
289  // the reference to the patch, size and content stay the same.
290  forAll(oldPatchFields, patchi)
291  {
292  if (oldPatchFields.set(patchi))
293  {
294  const fvsPatchField<Type>& pfld = oldPatchFields[patchi];
295 
296  labelList dummyMap(identity(pfld.size()));
297  directFvPatchFieldMapper dummyMapper(dummyMap);
298 
299  newPatchFields.set
300  (
301  patchi,
303  (
304  pfld,
305  tgtMesh_.boundary()[patchi],
307  dummyMapper
308  )
309  );
310  }
311  }
312 
313  // Add some empty patches on remaining patches
314  // (... probably processor patches)
315  forAll(newPatchFields, patchi)
316  {
317  if (!newPatchFields.set(patchi))
318  {
319  newPatchFields.set
320  (
321  patchi,
323  (
325  tgtMesh_.boundary()[patchi],
327  )
328  );
329  }
330  }
331 
332 
333  // Return geometric field
334  return tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
335  (
336  std::move(internalField),
337  newPatchFields
338  );
339 }
340 
341 
342 template<class Type>
345 (
346  const IOobject& fieldObject
347 ) const
348 {
349  // Read field
350  DimensionedField<Type, volMesh> fld
351  (
352  fieldObject,
353  srcMesh_
354  );
355 
356  // Distribute
357  return distributeField(fld);
358 }
359 
360 
361 template<class Type>
364 (
365  const IOobject& fieldObject
366 ) const
367 {
368  // Read field
369  GeometricField<Type, fvPatchField, volMesh> fld
370  (
371  fieldObject,
372  srcMesh_
373  );
374 
375  // Distribute
376  return distributeField(fld);
377 }
378 
379 
380 template<class Type>
383 (
384  const IOobject& fieldObject
385 ) const
386 {
387  // Read field
388  GeometricField<Type, fvsPatchField, surfaceMesh> fld
389  (
390  fieldObject,
391  srcMesh_
392  );
393 
394  // Distribute
395  return distributeField(fld);
396 }
397 
398 
399 template<class Type>
401 (
402  const IOobjectList& objects,
403  const wordRes& selectedFields
404 ) const
405 {
406  typedef DimensionedField<Type, volMesh> fieldType;
407 
408  // Available fields, sorted order
409  const wordList fieldNames =
410  (
411  selectedFields.empty()
412  ? objects.sortedNames<fieldType>()
413  : objects.sortedNames<fieldType>(selectedFields)
414  );
415 
416  label nFields = 0;
417  for (const word& fieldName : fieldNames)
418  {
419  if ("cellDist" == fieldName)
420  {
421  // There is an odd chance this is an internal field
422  continue;
423  }
424  if (verbose_)
425  {
426  if (!nFields)
427  {
428  Info<< " Reconstructing "
429  << fieldType::typeName << "s\n" << nl;
430  }
431  Info<< " " << fieldName << nl;
432  }
433  ++nFields;
434 
435  tmp<fieldType> tfld
436  (
437  distributeInternalField<Type>(*(objects[fieldName]))
438  );
439 
440  if (isWriteProc_.good())
441  {
442  if (isWriteProc_)
443  {
444  tfld().write();
445  }
446  }
447  else if (writeHandler_ && writeHandler_->good())
448  {
449  auto oldHandler = fileOperation::fileHandler(writeHandler_);
450  const label oldComm = UPstream::commWorld(fileHandler().comm());
451 
452  tfld().write();
453 
454  writeHandler_ = fileOperation::fileHandler(oldHandler);
455  UPstream::commWorld(oldComm);
456  }
457  }
458 
459  if (nFields && verbose_) Info<< endl;
460  return nFields;
461 }
462 
463 
464 template<class Type>
466 (
467  const IOobjectList& objects,
468  const wordRes& selectedFields
469 ) const
470 {
471  typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
472 
473  // Available fields, sorted order
474  const wordList fieldNames =
475  (
476  selectedFields.empty()
477  ? objects.sortedNames<fieldType>()
478  : objects.sortedNames<fieldType>(selectedFields)
479  );
480 
481  label nFields = 0;
482  for (const word& fieldName : fieldNames)
483  {
484  if ("cellDist" == fieldName)
485  {
486  continue;
487  }
488  if (verbose_)
489  {
490  if (!nFields)
491  {
492  Info<< " Reconstructing "
493  << fieldType::typeName << "s\n" << nl;
494  }
495  Info<< " " << fieldName << nl;
496  }
497  ++nFields;
498 
499  tmp<fieldType> tfld
500  (
501  distributeVolumeField<Type>(*(objects[fieldName]))
502  );
503 
504  if (isWriteProc_.good())
505  {
506  if (isWriteProc_)
507  {
508  tfld().write();
509  }
510  }
511  else if (writeHandler_ && writeHandler_->good())
512  {
513  auto oldHandler = fileOperation::fileHandler(writeHandler_);
514  const label oldComm = UPstream::commWorld(fileHandler().comm());
515 
516  tfld().write();
517 
518  writeHandler_ = fileOperation::fileHandler(oldHandler);
519  UPstream::commWorld(oldComm);
520  }
521  }
522 
523  if (nFields && verbose_) Info<< endl;
524  return nFields;
525 }
526 
527 
528 template<class Type>
530 (
531  const IOobjectList& objects,
532  const wordRes& selectedFields
533 ) const
534 {
535  typedef GeometricField<Type, fvsPatchField, surfaceMesh> fieldType;
536 
537  // Available fields, sorted order
538  const wordList fieldNames =
539  (
540  selectedFields.empty()
541  ? objects.sortedNames<fieldType>()
542  : objects.sortedNames<fieldType>(selectedFields)
543  );
544 
545  label nFields = 0;
546  for (const word& fieldName : fieldNames)
547  {
548  if (verbose_)
549  {
550  if (!nFields)
551  {
552  Info<< " Reconstructing "
553  << fieldType::typeName << "s\n" << nl;
554  }
555  Info<< " " << fieldName << nl;
556  }
557  ++nFields;
558 
559  tmp<fieldType> tfld
560  (
561  distributeSurfaceField<Type>(*(objects[fieldName]))
562  );
563 
564  if (isWriteProc_.good())
565  {
566  if (isWriteProc_)
567  {
568  tfld().write();
569  }
570  }
571  else if (writeHandler_ && writeHandler_->good())
572  {
573  auto oldHandler = fileOperation::fileHandler(writeHandler_);
574  const label oldComm = UPstream::commWorld(fileHandler().comm());
575 
576  tfld().write();
577 
578  writeHandler_ = fileOperation::fileHandler(oldHandler);
579  UPstream::commWorld(oldComm);
580  }
581  }
582 
583  if (nFields && verbose_) Info<< endl;
584  return nFields;
585 }
586 
587 
588 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > distributeSurfaceField(const IOobject &fieldObject) const
Read and distribute surface field.
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
label distributeVolumeFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected volume fields.
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
label distributeInternalFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, redistribute and write all/selected volume internal fields.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
static const char *const typeName
Typename for Field.
Definition: Field.H:86
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
const mapDistribute & cellMap() const noexcept
Cell distribute map.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
DistributedFieldMapper< directFieldMapper > distributedFieldMapper
A directFieldMapper with distributed (with local or remote) quantities.
static const UList< label > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:45
DistributedFieldMapper< directFvPatchFieldMapper > distributedFvPatchFieldMapper
A directFvPatchFieldMapper with direct mapping from local or remote quantities.
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds)
Definition: UPstream.H:431
tmp< DimensionedField< Type, volMesh > > distributeInternalField(const IOobject &fieldObject) const
Read and distribute volume internal field.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
tmp< DimensionedField< Type, volMesh > > distributeField(const DimensionedField< Type, volMesh > &) const
Redistribute volume internal field.
List< word > wordList
List of word.
Definition: fileName.H:58
tmp< GeometricField< Type, fvPatchField, volMesh > > distributeVolumeField(const IOobject &fieldObject) const
Read and distribute volume field.
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
static const fileOperation & fileHandler()
Return the current file handler. Will create the default file handler if necessary.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
label distributeSurfaceFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, reconstruct and write all/selected surface fields.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133