fvFieldReconstructorReconstructFields.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) 2011-2016 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 
26 #include "fvFieldReconstructor.H"
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 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
39 (
40  const IOobject& fieldIoObject,
41  const PtrList<DimensionedField<Type, volMesh>>& procFields
42 ) const
43 {
44  // Create the internalField
45  Field<Type> internalField(mesh_.nCells());
46 
47  forAll(procMeshes_, proci)
48  {
49  const DimensionedField<Type, volMesh>& procField = procFields[proci];
50 
51  // Set the cell values in the reconstructed field
52  internalField.rmap
53  (
54  procField.field(),
55  cellProcAddressing_[proci]
56  );
57  }
58 
60  (
62  (
63  fieldIoObject,
64  mesh_,
65  procFields[0].dimensions(),
66  internalField
67  )
68  );
69 }
70 
71 
72 template<class Type>
75 (
76  const IOobject& fieldIoObject
77 ) const
78 {
79  // Read the field for all the processors
81  (
82  procMeshes_.size()
83  );
84 
85  forAll(procMeshes_, proci)
86  {
87  procFields.set
88  (
89  proci,
91  (
92  IOobject
93  (
94  fieldIoObject.name(),
95  procMeshes_[proci].time().timeName(),
96  procMeshes_[proci],
97  IOobject::MUST_READ,
98  IOobject::NO_WRITE
99  ),
100  procMeshes_[proci]
101  )
102  );
103  }
104 
105 
106  return reconstructFvVolumeInternalField
107  (
108  IOobject
109  (
110  fieldIoObject.name(),
111  mesh_.time().timeName(),
112  mesh_,
113  IOobject::NO_READ,
114  IOobject::NO_WRITE
115  ),
116  procFields
117  );
118 }
119 
120 
121 template<class Type>
124 (
125  const IOobject& fieldIoObject,
127 ) const
128 {
129  // Create the internalField
130  Field<Type> internalField(mesh_.nCells());
131 
132  // Create the patch fields
133  PtrList<fvPatchField<Type>> patchFields(mesh_.boundary().size());
134 
135  forAll(procFields, proci)
136  {
138  procFields[proci];
139 
140  // Set the cell values in the reconstructed field
141  internalField.rmap
142  (
143  procField.primitiveField(),
144  cellProcAddressing_[proci]
145  );
146 
147  // Set the boundary patch values in the reconstructed field
148  forAll(boundaryProcAddressing_[proci], patchi)
149  {
150  // Get patch index of the original patch
151  const label curBPatch = boundaryProcAddressing_[proci][patchi];
152 
153  // Get addressing slice for this patch
154  const labelList::subList cp =
155  procField.mesh().boundary()[patchi].patchSlice
156  (
157  faceProcAddressing_[proci]
158  );
159 
160  // check if the boundary patch is not a processor patch
161  if (curBPatch >= 0)
162  {
163  // Regular patch. Fast looping
164 
165  if (!patchFields(curBPatch))
166  {
167  patchFields.set
168  (
169  curBPatch,
171  (
172  procField.boundaryField()[patchi],
173  mesh_.boundary()[curBPatch],
176  (
177  mesh_.boundary()[curBPatch].size()
178  )
179  )
180  );
181  }
182 
183  const label curPatchStart =
184  mesh_.boundaryMesh()[curBPatch].start();
185 
186  labelList reverseAddressing(cp.size());
187 
188  forAll(cp, facei)
189  {
190  // Check
191  if (cp[facei] <= 0)
192  {
194  << "Processor " << proci
195  << " patch "
196  << procField.mesh().boundary()[patchi].name()
197  << " face " << facei
198  << " originates from reversed face since "
199  << cp[facei]
200  << exit(FatalError);
201  }
202 
203  // Subtract one to take into account offsets for
204  // face direction.
205  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
206  }
207 
208 
209  patchFields[curBPatch].rmap
210  (
211  procField.boundaryField()[patchi],
212  reverseAddressing
213  );
214  }
215  else
216  {
217  const Field<Type>& curProcPatch =
218  procField.boundaryField()[patchi];
219 
220  // In processor patches, there's a mix of internal faces (some
221  // of them turned) and possible cyclics. Slow loop
222  forAll(cp, facei)
223  {
224  // Subtract one to take into account offsets for
225  // face direction.
226  label curF = cp[facei] - 1;
227 
228  // Is the face on the boundary?
229  if (curF >= mesh_.nInternalFaces())
230  {
231  label curBPatch = mesh_.boundaryMesh().whichPatch(curF);
232 
233  if (!patchFields(curBPatch))
234  {
235  patchFields.set
236  (
237  curBPatch,
239  (
240  mesh_.boundary()[curBPatch].type(),
241  mesh_.boundary()[curBPatch],
243  )
244  );
245  }
246 
247  // add the face
248  label curPatchFace =
249  mesh_.boundaryMesh()
250  [curBPatch].whichFace(curF);
251 
252  patchFields[curBPatch][curPatchFace] =
253  curProcPatch[facei];
254  }
255  }
256  }
257  }
258  }
259 
260  forAll(mesh_.boundary(), patchi)
261  {
262  // add empty patches
263  if
264  (
265  isType<emptyFvPatch>(mesh_.boundary()[patchi])
266  && !patchFields(patchi)
267  )
268  {
269  patchFields.set
270  (
271  patchi,
273  (
275  mesh_.boundary()[patchi],
277  )
278  );
279  }
280  }
281 
282 
283  // Now construct and write the field
284  // setting the internalField and patchFields
286  (
288  (
289  fieldIoObject,
290  mesh_,
291  procFields[0].dimensions(),
292  internalField,
293  patchFields
294  )
295  );
296 }
297 
298 
299 template<class Type>
302 (
303  const IOobject& fieldIoObject
304 ) const
305 {
306  // Read the field for all the processors
308  (
309  procMeshes_.size()
310  );
311 
312  forAll(procMeshes_, proci)
313  {
314  procFields.set
315  (
316  proci,
318  (
319  IOobject
320  (
321  fieldIoObject.name(),
322  procMeshes_[proci].time().timeName(),
323  procMeshes_[proci],
324  IOobject::MUST_READ,
325  IOobject::NO_WRITE
326  ),
327  procMeshes_[proci]
328  )
329  );
330  }
331 
332  return reconstructFvVolumeField
333  (
334  IOobject
335  (
336  fieldIoObject.name(),
337  mesh_.time().timeName(),
338  mesh_,
339  IOobject::NO_READ,
340  IOobject::NO_WRITE
341  ),
342  procFields
343  );
344 }
345 
346 
347 template<class Type>
350 (
351  const IOobject& fieldIoObject,
353 ) const
354 {
355  // Create the internalField
356  Field<Type> internalField(mesh_.nInternalFaces());
357 
358  // Create the patch fields
359  PtrList<fvsPatchField<Type>> patchFields(mesh_.boundary().size());
360 
361 
362  forAll(procMeshes_, proci)
363  {
365  procFields[proci];
366 
367  // Set the face values in the reconstructed field
368 
369  // It is necessary to create a copy of the addressing array to
370  // take care of the face direction offset trick.
371  //
372  {
373  const labelList& faceMap = faceProcAddressing_[proci];
374 
375  // Correctly oriented copy of internal field
376  Field<Type> procInternalField(procField.primitiveField());
377  // Addressing into original field
378  labelList curAddr(procInternalField.size());
379 
380  forAll(procInternalField, addrI)
381  {
382  curAddr[addrI] = mag(faceMap[addrI])-1;
383  if (faceMap[addrI] < 0)
384  {
385  procInternalField[addrI] = -procInternalField[addrI];
386  }
387  }
388 
389  // Map
390  internalField.rmap(procInternalField, curAddr);
391  }
392 
393  // Set the boundary patch values in the reconstructed field
394  forAll(boundaryProcAddressing_[proci], patchi)
395  {
396  // Get patch index of the original patch
397  const label curBPatch = boundaryProcAddressing_[proci][patchi];
398 
399  // Get addressing slice for this patch
400  const labelList::subList cp =
401  procMeshes_[proci].boundary()[patchi].patchSlice
402  (
403  faceProcAddressing_[proci]
404  );
405 
406  // check if the boundary patch is not a processor patch
407  if (curBPatch >= 0)
408  {
409  // Regular patch. Fast looping
410 
411  if (!patchFields(curBPatch))
412  {
413  patchFields.set
414  (
415  curBPatch,
417  (
418  procField.boundaryField()[patchi],
419  mesh_.boundary()[curBPatch],
422  (
423  mesh_.boundary()[curBPatch].size()
424  )
425  )
426  );
427  }
428 
429  const label curPatchStart =
430  mesh_.boundaryMesh()[curBPatch].start();
431 
432  labelList reverseAddressing(cp.size());
433 
434  forAll(cp, facei)
435  {
436  // Subtract one to take into account offsets for
437  // face direction.
438  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
439  }
440 
441  patchFields[curBPatch].rmap
442  (
443  procField.boundaryField()[patchi],
444  reverseAddressing
445  );
446  }
447  else
448  {
449  const Field<Type>& curProcPatch =
450  procField.boundaryField()[patchi];
451 
452  // In processor patches, there's a mix of internal faces (some
453  // of them turned) and possible cyclics. Slow loop
454  forAll(cp, facei)
455  {
456  label curF = cp[facei] - 1;
457 
458  // Is the face turned the right side round
459  if (curF >= 0)
460  {
461  // Is the face on the boundary?
462  if (curF >= mesh_.nInternalFaces())
463  {
464  label curBPatch =
465  mesh_.boundaryMesh().whichPatch(curF);
466 
467  if (!patchFields(curBPatch))
468  {
469  patchFields.set
470  (
471  curBPatch,
473  (
474  mesh_.boundary()[curBPatch].type(),
475  mesh_.boundary()[curBPatch],
477  ::null()
478  )
479  );
480  }
481 
482  // add the face
483  label curPatchFace =
484  mesh_.boundaryMesh()
485  [curBPatch].whichFace(curF);
486 
487  patchFields[curBPatch][curPatchFace] =
488  curProcPatch[facei];
489  }
490  else
491  {
492  // Internal face
493  internalField[curF] = curProcPatch[facei];
494  }
495  }
496  }
497  }
498  }
499  }
500 
501  forAll(mesh_.boundary(), patchi)
502  {
503  // add empty patches
504  if
505  (
506  isType<emptyFvPatch>(mesh_.boundary()[patchi])
507  && !patchFields(patchi)
508  )
509  {
510  patchFields.set
511  (
512  patchi,
514  (
516  mesh_.boundary()[patchi],
518  )
519  );
520  }
521  }
522 
523 
524  // Now construct and write the field
525  // setting the internalField and patchFields
527  (
529  (
530  fieldIoObject,
531  mesh_,
532  procFields[0].dimensions(),
533  internalField,
534  patchFields
535  )
536  );
537 }
538 
539 
540 template<class Type>
543 (
544  const IOobject& fieldIoObject
545 ) const
546 {
547  // Read the field for all the processors
549  (
550  procMeshes_.size()
551  );
552 
553  forAll(procMeshes_, proci)
554  {
555  procFields.set
556  (
557  proci,
559  (
560  IOobject
561  (
562  fieldIoObject.name(),
563  procMeshes_[proci].time().timeName(),
564  procMeshes_[proci],
565  IOobject::MUST_READ,
566  IOobject::NO_WRITE
567  ),
568  procMeshes_[proci]
569  )
570  );
571  }
572 
573  return reconstructFvSurfaceField
574  (
575  IOobject
576  (
577  fieldIoObject.name(),
578  mesh_.time().timeName(),
579  mesh_,
580  IOobject::NO_READ,
581  IOobject::NO_WRITE
582  ),
583  procFields
584  );
585 }
586 
587 
588 template<class Type>
590 (
591  const IOobjectList& objects,
592  const HashSet<word>& selectedFields
593 )
594 {
595  const word& fieldClassName = DimensionedField<Type, volMesh>::typeName;
596 
597  IOobjectList fields = objects.lookupClass(fieldClassName);
598 
599  if (fields.size())
600  {
601  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
602 
603  forAllConstIter(IOobjectList, fields, fieldIter)
604  {
605  if
606  (
607  selectedFields.empty()
608  || selectedFields.found(fieldIter()->name())
609  )
610  {
611  Info<< " " << fieldIter()->name() << endl;
612 
613  reconstructFvVolumeInternalField<Type>(*fieldIter())().write();
614 
615  nReconstructed_++;
616  }
617  }
618  Info<< endl;
619  }
620 }
621 
622 
623 template<class Type>
625 (
626  const IOobjectList& objects,
627  const HashSet<word>& selectedFields
628 )
629 {
630  const word& fieldClassName =
632 
633  IOobjectList fields = objects.lookupClass(fieldClassName);
634 
635  if (fields.size())
636  {
637  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
638 
639  forAllConstIter(IOobjectList, fields, fieldIter)
640  {
641  if
642  (
643  selectedFields.empty()
644  || selectedFields.found(fieldIter()->name())
645  )
646  {
647  Info<< " " << fieldIter()->name() << endl;
648 
649  reconstructFvVolumeField<Type>(*fieldIter())().write();
650 
651  nReconstructed_++;
652  }
653  }
654  Info<< endl;
655  }
656 }
657 
658 
659 template<class Type>
661 (
662  const IOobjectList& objects,
663  const HashSet<word>& selectedFields
664 )
665 {
666  const word& fieldClassName =
668 
669  IOobjectList fields = objects.lookupClass(fieldClassName);
670 
671  if (fields.size())
672  {
673  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
674 
675  forAllConstIter(IOobjectList, fields, fieldIter)
676  {
677  if
678  (
679  selectedFields.empty()
680  || selectedFields.found(fieldIter()->name())
681  )
682  {
683  Info<< " " << fieldIter()->name() << endl;
684 
685  reconstructFvSurfaceField<Type>(*fieldIter())().write();
686 
687  nReconstructed_++;
688  }
689  }
690  Info<< endl;
691  }
692 }
693 
694 
695 // ************************************************************************* //
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructFvSurfaceField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvsPatchField, surfaceMesh >> &) const
Reconstruct surface field.
A HashTable with keys but without contents.
Definition: HashSet.H:59
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:197
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
tmp< DimensionedField< Type, volMesh > > reconstructFvVolumeInternalField(const IOobject &fieldIoObject, const PtrList< DimensionedField< Type, volMesh >> &procFields) const
Reconstruct volume internal field.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
A List obtained as a section of another List.
Definition: SubList.H:53
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructFvVolumeField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
Reconstruct volume field.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
bool empty() const
Return true if the hash table is empty.
bool found(const Key &) const
Return true if hashedEntry is found in table.
const Field< Type > & field() const
const volScalarField & cp
Mapper for sizing only - does not do any actual mapping.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:577
label patchi
Foam::emptyFvsPatchField.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const Mesh & mesh() const
Return mesh.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual Ostream & write(const token &)=0
Write next token to stream.
This boundary condition provides an &#39;empty&#39; condition for reduced dimensions cases, i.e. 1- and 2-D geometries. Apply this condition to patches whose normal is aligned to geometric directions that do not constitue solution directions.
void reconstructFvVolumeFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume fields.
A class for managing temporary objects.
Definition: PtrList.H:54
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
void reconstructFvVolumeInternalFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected volume internal fields.
const word & name() const
Return name.
Definition: IOobject.H:260