parFvFieldReconstructorReconstructFields.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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "Time.H"
28 #include "PtrList.H"
29 #include "fvPatchFields.H"
30 #include "emptyFvPatch.H"
31 #include "emptyFvPatchField.H"
32 #include "emptyFvsPatchField.H"
33 #include "IOobjectList.H"
34 #include "mapDistributePolyMesh.H"
35 #include "processorFvPatch.H"
36 
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
43 template<class Type>
46 (
48 ) const
49 {
51  (
52  labelUList::null(),
53  distMap_.cellMap()
54  );
55 
56  Field<Type> internalField(fld, mapper);
57 
58  // Construct a volField
59  IOobject baseIO
60  (
61  fld.name(),
62  baseMesh_.time().timeName(),
63  fld.local(),
64  baseMesh_,
65  IOobject::NO_READ,
66  IOobject::NO_WRITE
67  );
68 
70  (
72  (
73  baseIO,
74  baseMesh_,
75  fld.dimensions(),
77  )
78  );
79 }
80 
81 
82 template<class Type>
85 (
86  const IOobject& fieldIoObject
87 ) const
88 {
89  // Read the field
91  (
92  fieldIoObject,
93  procMesh_
94  );
95 
96  // Distribute onto baseMesh
97  return reconstructFvVolumeInternalField(fld);
98 }
99 
100 
101 // Reconstruct a field onto the baseMesh
102 template<class Type>
105 (
107 ) const
108 {
109  // Create the internalField by remote mapping
110  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
111 
113  (
114  labelUList::null(),
115  distMap_.cellMap()
116  );
117 
118  Field<Type> internalField(fld.internalField(), mapper);
119 
120 
121 
122  // Create the patchFields by remote mapping
123  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
124  // Note: patchFields still on mesh, not baseMesh
125 
126  PtrList<fvPatchField<Type> > patchFields(fld.mesh().boundary().size());
127 
128  const typename GeometricField
129  <
130  Type,
131  fvPatchField,
132  volMesh
133  >::GeometricBoundaryField& 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 
143  (
144  labelUList::null(),
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  (
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_,
212  IOobject::NO_READ,
213  IOobject::NO_WRITE
214  );
215 
217  (
219  (
220  baseIO,
221  baseMesh_,
222  fld.dimensions(),
224  basePatchFields
225  )
226  );
227 }
228 
229 
230 template<class Type>
233 (
234  const IOobject& fieldIoObject
235 ) const
236 {
237  // Read the field
239  (
240  fieldIoObject,
241  procMesh_
242  );
243 
244  // Distribute onto baseMesh
245  return reconstructFvVolumeField(fld);
246 }
247 
248 
249 template<class Type>
252 (
254 ) const
255 {
256  // Create the internalField by remote mapping
257  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
258 
260  (
261  labelUList::null(),
262  distMap_.faceMap()
263  );
264 
265  // Create flat field of internalField + all patch fields
266  Field<Type> flatFld(fld.mesh().nFaces(), pTraits<Type>::zero);
267  SubList<Type>(flatFld, fld.internalField().size()).assign
268  (
269  fld.internalField()
270  );
271  forAll(fld.boundaryField(), patchI)
272  {
273  const fvsPatchField<Type>& fvp = fld.boundaryField()[patchI];
274 
275  SubList<Type>(flatFld, fvp.size(), fvp.patch().start()).assign(fvp);
276  }
277 
278  // Map all faces
279  Field<Type> internalField(flatFld, mapper);
280 
281  // Trim to internal faces (note: could also have special mapper)
282  internalField.setSize
283  (
284  min
285  (
286  internalField.size(),
287  baseMesh_.nInternalFaces()
288  )
289  );
290 
291 
292  // Create the patchFields by remote mapping
293  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
294  // Note: patchFields still on mesh, not baseMesh
295 
296  PtrList<fvsPatchField<Type> > patchFields(fld.mesh().boundary().size());
297 
298  const typename GeometricField
299  <
300  Type,
303  >::GeometricBoundaryField& bfld = fld.boundaryField();
304 
305  forAll(bfld, patchI)
306  {
307  if (patchFaceMaps_.set(patchI))
308  {
309  // Clone local patch field
310  patchFields.set(patchI, bfld[patchI].clone());
311 
313  (
314  labelUList::null(),
315  patchFaceMaps_[patchI]
316  );
317 
318  // Map into local copy
319  patchFields[patchI].autoMap(mapper);
320  }
321  }
322 
323 
324  PtrList<fvsPatchField<Type> > basePatchFields
325  (
326  baseMesh_.boundary().size()
327  );
328 
329  // Clone the patchFields onto the base patches. This is just to reset
330  // the reference to the patch, size and content stay the same.
331  forAll(patchFields, patchI)
332  {
333  if (patchFields.set(patchI))
334  {
335  const fvPatch& basePatch = baseMesh_.boundary()[patchI];
336 
337  const fvsPatchField<Type>& pfld = patchFields[patchI];
338 
339  labelList dummyMap(identity(pfld.size()));
340  directFvPatchFieldMapper dummyMapper(dummyMap);
341 
342  basePatchFields.set
343  (
344  patchI,
346  (
347  pfld,
348  basePatch,
350  dummyMapper
351  )
352  );
353  }
354  }
355 
356  // Add some empty patches on remaining patches (tbd.probably processor
357  // patches)
358  forAll(basePatchFields, patchI)
359  {
360  if (patchI >= patchFields.size() || !patchFields.set(patchI))
361  {
362  basePatchFields.set
363  (
364  patchI,
366  (
368  baseMesh_.boundary()[patchI],
370  )
371  );
372  }
373  }
374 
375  // Construct a volField
376  IOobject baseIO
377  (
378  fld.name(),
379  baseMesh_.time().timeName(),
380  fld.local(),
381  baseMesh_,
382  IOobject::NO_READ,
383  IOobject::NO_WRITE
384  );
385 
387  (
389  (
390  baseIO,
391  baseMesh_,
392  fld.dimensions(),
394  basePatchFields
395  )
396  );
397 }
398 
399 
400 template<class Type>
403 (
404  const IOobject& fieldIoObject
405 ) const
406 {
407  // Read the field
409  (
410  fieldIoObject,
411  procMesh_
412  );
413 
414  return reconstructFvSurfaceField(fld);
415 }
416 
417 
418 template<class Type>
420 (
421  const IOobjectList& objects,
422  const HashSet<word>& selectedFields
423 ) const
424 {
425  const word& fieldClassName = DimensionedField<Type, volMesh>::typeName;
426 
427  IOobjectList fields = objects.lookupClass(fieldClassName);
428 
429  if (fields.size())
430  {
431  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
432 
433  forAllConstIter(IOobjectList, fields, fieldIter)
434  {
435  if
436  (
437  selectedFields.empty()
438  || selectedFields.found(fieldIter()->name())
439  )
440  {
441  Info<< " " << fieldIter()->name() << endl;
442 
444  (
445  reconstructFvVolumeInternalField<Type>(*fieldIter())
446  );
447 
448  if (isWriteProc_)
449  {
450  tfld().write();
451  }
452  }
453  }
454  Info<< endl;
455  }
456 }
457 
458 
459 template<class Type>
461 (
462  const IOobjectList& objects,
463  const HashSet<word>& selectedFields
464 ) const
465 {
466  const word& fieldClassName =
468 
469  IOobjectList fields = objects.lookupClass(fieldClassName);
470 
471  if (fields.size())
472  {
473  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
474 
475  forAllConstIter(IOobjectList, fields, fieldIter)
476  {
477  const word& name = fieldIter()->name();
478 
479  if
480  (
481  (selectedFields.empty() || selectedFields.found(name))
482  && name != "cellDist"
483  )
484  {
485  Info<< " " << name << endl;
486 
488  (
489  reconstructFvVolumeField<Type>(*fieldIter())
490  );
491  if (isWriteProc_)
492  {
493  tfld().write();
494  }
495  }
496  }
497  Info<< endl;
498  }
499 }
500 
501 
502 template<class Type>
504 (
505  const IOobjectList& objects,
506  const HashSet<word>& selectedFields
507 ) const
508 {
509  const word& fieldClassName =
511 
512  IOobjectList fields = objects.lookupClass(fieldClassName);
513 
514  if (fields.size())
515  {
516  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
517 
518  forAllConstIter(IOobjectList, fields, fieldIter)
519  {
520  if
521  (
522  selectedFields.empty()
523  || selectedFields.found(fieldIter()->name())
524  )
525  {
526  Info<< " " << fieldIter()->name() << endl;
527 
529  (
530  reconstructFvSurfaceField<Type>(*fieldIter())
531  );
532  if (isWriteProc_)
533  {
534  tfld().write();
535  }
536  }
537  }
538  Info<< endl;
539  }
540 }
541 
542 
543 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::distributedUnallocatedDirectFieldMapper
FieldMapper with direct mapping from remote quantities.
Definition: distributedUnallocatedDirectFieldMapper.H:43
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
processorFvPatch.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalField
tmp< DimensionedField< Type, volMesh > > reconstructFvVolumeInternalField(const DimensionedField< Type, volMesh > &) const
Reconstruct volume internal field.
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:47
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::surfaceMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: surfaceMesh.H:47
Foam::emptyFvsPatchField
Foam::emptyFvsPatchField.
Definition: emptyFvsPatchField.H:50
Foam::fvsPatchField< Type >
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
emptyFvsPatchField.H
Foam::HashSet
A HashTable with keys but without contents.
Definition: HashSet.H:59
parFvFieldReconstructor.H
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
mapDistributePolyMesh.H
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::Field< Type >
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalFields
void reconstructFvVolumeInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected volume internal fields.
Definition: parFvFieldReconstructorReconstructFields.C:420
Foam::HashTable< nil, word, string::hash >::empty
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
Foam::Info
messageStream Info
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
directFvPatchFieldMapper.H
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::emptyFvPatchField
This boundary condition provides an 'empty' condition for reduced dimensions cases,...
Definition: emptyFvPatchField.H:66
Foam::parFvFieldReconstructor::reconstructFvSurfaceField
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructFvSurfaceField(const GeometricField< Type, fvsPatchField, surfaceMesh > &) const
Reconstruct surface field.
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::HashTable< nil, word, string::hash >::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
emptyFvPatch.H
Foam::parFvFieldReconstructor::reconstructFvSurfaceFields
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected surface fields.
Definition: parFvFieldReconstructorReconstructFields.C:504
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
distributedUnallocatedDirectFieldMapper.H
internalField
conserve internalField()+
emptyFvPatchField.H
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
fvPatchFields.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::parFvFieldReconstructor::reconstructFvVolumeField
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructFvVolumeField(const GeometricField< Type, fvPatchField, volMesh > &fld) const
Reconstruct volume field.
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
distributedUnallocatedDirectFvPatchFieldMapper.H
Foam::directFvPatchFieldMapper
direct fvPatchFieldMapper
Definition: directFvPatchFieldMapper.H:45
PtrList.H
Foam::parFvFieldReconstructor::reconstructFvVolumeFields
void reconstructFvVolumeFields(const IOobjectList &objects, const HashSet< word > &selectedFields) const
Read, reconstruct and write all/selected volume fields.
Definition: parFvFieldReconstructorReconstructFields.C:461
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::distributedUnallocatedDirectFvPatchFieldMapper
FieldMapper with direct mapping from remote quantities.
Definition: distributedUnallocatedDirectFvPatchFieldMapper.H:45
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51