parFvFieldReconstructorFields.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-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 "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 
43 
44 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
45 
46 template<class Type>
49 (
50  const DimensionedField<Type, volMesh>& fld
51 ) const
52 {
53  distributedUnallocatedDirectFieldMapper mapper
54  (
56  distMap_.cellMap()
57  );
58 
59  Field<Type> internalField(fld, mapper);
60 
61  // Construct a volField
62  IOobject baseIO
63  (
64  fld.name(),
65  baseMesh_.time().timeName(),
66  fld.local(),
67  baseMesh_,
70  );
71 
72  auto tfield = tmp<DimensionedField<Type, volMesh>>::New
73  (
74  baseIO,
75  baseMesh_,
76  fld.dimensions(),
77  internalField
78  );
79 
80  tfield.ref().oriented() = fld.oriented();
81 
82  return tfield;
83 }
84 
85 
86 template<class Type>
89 (
90  const IOobject& fieldIoObject
91 ) const
92 {
93  // Read the field
94  DimensionedField<Type, volMesh> fld
95  (
96  fieldIoObject,
97  procMesh_
98  );
99 
100  // Distribute onto baseMesh
101  return reconstructFvVolumeInternalField(fld);
102 }
103 
104 
105 // Reconstruct a field onto the baseMesh
106 template<class Type>
109 (
110  const GeometricField<Type, fvPatchField, volMesh>& fld
111 ) const
112 {
113  // Create the internalField by remote mapping
114  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
115 
116  distributedUnallocatedDirectFieldMapper mapper
117  (
119  distMap_.cellMap()
120  );
121 
122  Field<Type> internalField(fld.internalField(), mapper);
123 
124 
125 
126  // Create the patchFields by remote mapping
127  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
128  // Note: patchFields still on mesh, not baseMesh
129 
130  PtrList<fvPatchField<Type>> patchFields(fld.mesh().boundary().size());
131 
132  const typename GeometricField<Type, fvPatchField, volMesh>::Boundary&
133  bfld = fld.boundaryField();
134 
135  forAll(bfld, patchI)
136  {
137  if (patchFaceMaps_.set(patchI))
138  {
139  // Clone local patch field
140  patchFields.set(patchI, bfld[patchI].clone());
141 
142  distributedUnallocatedDirectFvPatchFieldMapper mapper
143  (
145  patchFaceMaps_[patchI]
146  );
147 
148  // Map into local copy
149  patchFields[patchI].autoMap(mapper);
150  }
151  }
152 
153 
154  PtrList<fvPatchField<Type>> basePatchFields
155  (
156  baseMesh_.boundary().size()
157  );
158 
159  // Clone the patchFields onto the base patches. This is just to reset
160  // the reference to the patch, size and content stay the same.
161  forAll(patchFields, patchI)
162  {
163  if (patchFields.set(patchI))
164  {
165  const fvPatch& basePatch = baseMesh_.boundary()[patchI];
166 
167  const fvPatchField<Type>& pfld = patchFields[patchI];
168 
169  labelList dummyMap(identity(pfld.size()));
170  directFvPatchFieldMapper dummyMapper(dummyMap);
171 
172  basePatchFields.set
173  (
174  patchI,
176  (
177  pfld,
178  basePatch,
180  dummyMapper
181  )
182  );
183  }
184  }
185 
186  // Add some empty patches on remaining patches (tbd.probably processor
187  // patches)
188  forAll(basePatchFields, patchI)
189  {
190  if (patchI >= patchFields.size() || !patchFields.set(patchI))
191  {
192  basePatchFields.set
193  (
194  patchI,
196  (
197  emptyFvPatchField<Type>::typeName,
198  baseMesh_.boundary()[patchI],
200  )
201  );
202  }
203  }
204 
205  // Construct a volField
206  IOobject baseIO
207  (
208  fld.name(),
209  baseMesh_.time().timeName(),
210  fld.local(),
211  baseMesh_,
214  );
215 
216  auto tfield = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
217  (
218  baseIO,
219  baseMesh_,
220  fld.dimensions(),
221  internalField,
222  basePatchFields
223  );
224 
225  tfield.ref().oriented()= fld.oriented();
226 
227  return tfield;
228 }
229 
230 
231 template<class Type>
234 (
235  const IOobject& fieldIoObject
236 ) const
237 {
238  // Read the field
239  GeometricField<Type, fvPatchField, volMesh> fld
240  (
241  fieldIoObject,
242  procMesh_
243  );
244 
245  // Distribute onto baseMesh
246  return reconstructFvVolumeField(fld);
247 }
248 
249 
250 template<class Type>
253 (
254  const GeometricField<Type, fvsPatchField, surfaceMesh>& fld
255 ) const
256 {
257  // Create the internalField by remote mapping
258  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
259 
260  distributedUnallocatedDirectFieldMapper mapper
261  (
263  distMap_.faceMap()
264  );
265 
266  // Create flat field of internalField + all patch fields
267  Field<Type> flatFld(fld.mesh().nFaces(), Type(Zero));
268  SubList<Type>(flatFld, fld.internalField().size()) = fld.internalField();
269  forAll(fld.boundaryField(), patchI)
270  {
271  const fvsPatchField<Type>& fvp = fld.boundaryField()[patchI];
272 
273  SubList<Type>(flatFld, fvp.size(), fvp.patch().start()) = fvp;
274  }
275 
276  // Map all faces
277  Field<Type> internalField(flatFld, mapper, fld.oriented()());
278 
279  // Trim to internal faces (note: could also have special mapper)
280  internalField.setSize
281  (
282  min
283  (
284  internalField.size(),
285  baseMesh_.nInternalFaces()
286  )
287  );
288 
289 
290  // Create the patchFields by remote mapping
291  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
292  // Note: patchFields still on mesh, not baseMesh
293 
294  PtrList<fvsPatchField<Type>> patchFields(fld.mesh().boundary().size());
295 
296  const typename GeometricField<Type, fvsPatchField, surfaceMesh>::Boundary&
297  bfld = fld.boundaryField();
298 
299  forAll(bfld, patchI)
300  {
301  if (patchFaceMaps_.set(patchI))
302  {
303  // Clone local patch field
304  patchFields.set(patchI, bfld[patchI].clone());
305 
306  distributedUnallocatedDirectFvPatchFieldMapper mapper
307  (
309  patchFaceMaps_[patchI]
310  );
311 
312  // Map into local copy
313  patchFields[patchI].autoMap(mapper);
314  }
315  }
316 
317 
318  PtrList<fvsPatchField<Type>> basePatchFields
319  (
320  baseMesh_.boundary().size()
321  );
322 
323  // Clone the patchFields onto the base patches. This is just to reset
324  // the reference to the patch, size and content stay the same.
325  forAll(patchFields, patchI)
326  {
327  if (patchFields.set(patchI))
328  {
329  const fvPatch& basePatch = baseMesh_.boundary()[patchI];
330 
331  const fvsPatchField<Type>& pfld = patchFields[patchI];
332 
333  labelList dummyMap(identity(pfld.size()));
334  directFvPatchFieldMapper dummyMapper(dummyMap);
335 
336  basePatchFields.set
337  (
338  patchI,
340  (
341  pfld,
342  basePatch,
344  dummyMapper
345  )
346  );
347  }
348  }
349 
350  // Add some empty patches on remaining patches (tbd.probably processor
351  // patches)
352  forAll(basePatchFields, patchI)
353  {
354  if (patchI >= patchFields.size() || !patchFields.set(patchI))
355  {
356  basePatchFields.set
357  (
358  patchI,
360  (
361  emptyFvsPatchField<Type>::typeName,
362  baseMesh_.boundary()[patchI],
364  )
365  );
366  }
367  }
368 
369  // Construct a volField
370  IOobject baseIO
371  (
372  fld.name(),
373  baseMesh_.time().timeName(),
374  fld.local(),
375  baseMesh_,
378  );
379 
380  auto tfield = tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
381  (
382  baseIO,
383  baseMesh_,
384  fld.dimensions(),
385  internalField,
386  basePatchFields
387  );
388 
389  tfield.ref().oriented() = fld.oriented();
390 
391  return tfield;
392 }
393 
394 
395 template<class Type>
398 (
399  const IOobject& fieldIoObject
400 ) const
401 {
402  // Read the field
403  GeometricField<Type, fvsPatchField, surfaceMesh> fld
404  (
405  fieldIoObject,
406  procMesh_
407  );
408 
409  return reconstructFvSurfaceField(fld);
410 }
411 
412 
413 template<class Type>
415 (
416  const IOobjectList& objects,
417  const wordRes& selectedFields
418 ) const
419 {
420  typedef DimensionedField<Type, volMesh> fieldType;
421 
422  // Available fields, sorted order
423  const wordList fieldNames =
424  (
425  selectedFields.empty()
426  ? objects.sortedNames<fieldType>()
427  : objects.sortedNames<fieldType>(selectedFields)
428  );
429 
430  label nFields = 0;
431  for (const word& fieldName : fieldNames)
432  {
433  if (!nFields++)
434  {
435  Info<< " Reconstructing "
436  << fieldType::typeName << "s\n" << nl;
437  }
438 
439  Info<< " " << fieldName << nl;
440 
441  tmp<fieldType> tfld
442  (
443  reconstructFvVolumeInternalField<Type>(*(objects[fieldName]))
444  );
445  if (isWriteProc_)
446  {
447  tfld().write();
448  }
449  }
450 
451  if (nFields) Info<< endl;
452  return nFields;
453 }
454 
455 
456 template<class Type>
458 (
459  const IOobjectList& objects,
460  const wordRes& selectedFields
461 ) const
462 {
463  typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
464 
465  // Available fields, sorted order
466  const wordList fieldNames =
467  (
468  selectedFields.empty()
469  ? objects.sortedNames<fieldType>()
470  : objects.sortedNames<fieldType>(selectedFields)
471  );
472 
473  label nFields = 0;
474  for (const word& fieldName : fieldNames)
475  {
476  if ("cellDist" == fieldName)
477  {
478  continue;
479  }
480  if (!nFields++)
481  {
482  Info<< " Reconstructing "
483  << fieldType::typeName << "s\n" << nl;
484  }
485  Info<< " " << fieldName << nl;
486 
487  tmp<fieldType> tfld
488  (
489  reconstructFvVolumeField<Type>(*(objects[fieldName]))
490  );
491  if (isWriteProc_)
492  {
493  tfld().write();
494  }
495  }
496 
497  if (nFields) Info<< endl;
498  return nFields;
499 }
500 
501 
502 template<class Type>
504 (
505  const IOobjectList& objects,
506  const wordRes& selectedFields
507 ) const
508 {
509  typedef GeometricField<Type, fvsPatchField, surfaceMesh> fieldType;
510 
511  // Available fields, sorted order
512  const wordList fieldNames =
513  (
514  selectedFields.empty()
515  ? objects.sortedNames<fieldType>()
516  : objects.sortedNames<fieldType>(selectedFields)
517  );
518 
519  label nFields = 0;
520  for (const word& fieldName : fieldNames)
521  {
522  if (!nFields++)
523  {
524  Info<< " Reconstructing "
525  << fieldType::typeName << "s\n" << nl;
526  }
527  Info<< " " << fieldName << nl;
528 
529  tmp<fieldType> tfld
530  (
531  reconstructFvSurfaceField<Type>(*(objects[fieldName]))
532  );
533  if (isWriteProc_)
534  {
535  tfld().write();
536  }
537  }
538 
539  if (nFields) Info<< endl;
540  return nFields;
541 }
542 
543 
544 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalFields
label reconstructFvVolumeInternalFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, reconstruct and write all/selected volume internal fields.
distributedUnallocatedDirectFvPatchFieldMapper.H
processorFvPatch.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalField
tmp< DimensionedField< Type, volMesh > > reconstructFvVolumeInternalField(const DimensionedField< Type, volMesh > &) const
Reconstruct volume internal field.
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::DimensionedField::null
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
Definition: DimensionedFieldI.H:33
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
emptyFvsPatchField.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::mapDistributePolyMesh::cellMap
const mapDistribute & cellMap() const
Cell distribute map.
Definition: mapDistributePolyMesh.H:223
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::parFvFieldReconstructor::reconstructFvSurfaceFields
label reconstructFvSurfaceFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, reconstruct and write all/selected surface fields.
mapDistributePolyMesh.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
directFvPatchFieldMapper.H
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::fvsPatchField::New
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
Definition: fvsPatchFieldNew.C:74
Foam::parFvFieldReconstructor::reconstructFvSurfaceField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructFvSurfaceField(const GeometricField< Type, fvsPatchField, surfaceMesh > &) const
Reconstruct surface field.
fieldNames
const wordRes fieldNames(propsDict.getOrDefault< wordRes >("fields", wordRes()))
emptyFvPatch.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvPatchField::New
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
Definition: fvPatchFieldNew.C:88
Time.H
Foam::parFvFieldReconstructor::reconstructFvVolumeFields
label reconstructFvVolumeFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes()) const
Read, reconstruct and write all/selected volume fields.
emptyFvPatchField.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
fvPatchFields.H
Foam::parFvFieldReconstructor::reconstructFvVolumeField
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructFvVolumeField(const GeometricField< Type, fvPatchField, volMesh > &fld) const
Reconstruct volume field.
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
PtrList.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
distributedUnallocatedDirectFieldMapper.H
Foam::UList::null
static const UList< T > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:53
Foam::IOobject::NO_READ
Definition: IOobject.H:188
parFvFieldReconstructor.H