fvFieldReconstructorReconstructFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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 
370  {
371  // Assume all scalar surfaceFields are oriented flux fields
372  const labelList& faceMap = faceProcAddressing_[proci];
373 
374  // Correctly oriented copy of internal field
375  Field<Type> procInternalField(procField.primitiveField());
376 
377  // Addressing into original field
378  // It is necessary to create a copy of the addressing array to
379  // take care of the face direction offset trick.
380  labelList curAddr(procInternalField.size());
381 
382  forAll(procInternalField, i)
383  {
384  curAddr[i] = mag(faceMap[i]) - 1;
385  if (faceMap[i] < 0)
386  {
387  procInternalField[i] = -procInternalField[i];
388  }
389  }
390 
391  // Map
392  internalField.rmap(procInternalField, curAddr);
393  }
394  else
395  {
396  // Map
397  internalField.rmap
398  (
399  procField.primitiveField(),
400  mag(labelField(faceProcAddressing_[proci])) - 1
401  );
402  }
403 
404  // Set the boundary patch values in the reconstructed field
405  forAll(boundaryProcAddressing_[proci], patchi)
406  {
407  // Get patch index of the original patch
408  const label curBPatch = boundaryProcAddressing_[proci][patchi];
409 
410  // Get addressing slice for this patch
411  const labelList::subList cp =
412  procMeshes_[proci].boundary()[patchi].patchSlice
413  (
414  faceProcAddressing_[proci]
415  );
416 
417  // check if the boundary patch is not a processor patch
418  if (curBPatch >= 0)
419  {
420  // Regular patch. Fast looping
421 
422  if (!patchFields(curBPatch))
423  {
424  patchFields.set
425  (
426  curBPatch,
428  (
429  procField.boundaryField()[patchi],
430  mesh_.boundary()[curBPatch],
433  (
434  mesh_.boundary()[curBPatch].size()
435  )
436  )
437  );
438  }
439 
440  const label curPatchStart =
441  mesh_.boundaryMesh()[curBPatch].start();
442 
443  labelList reverseAddressing(cp.size());
444 
445  forAll(cp, facei)
446  {
447  // Subtract one to take into account offsets for
448  // face direction.
449  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
450  }
451 
452  patchFields[curBPatch].rmap
453  (
454  procField.boundaryField()[patchi],
455  reverseAddressing
456  );
457  }
458  else
459  {
460  const Field<Type>& curProcPatch =
461  procField.boundaryField()[patchi];
462 
463  // In processor patches, there's a mix of internal faces (some
464  // of them turned) and possible cyclics. Slow loop
465  forAll(cp, facei)
466  {
467  label curF = cp[facei] - 1;
468 
469  // Is the face turned the right side round
470  if (curF >= 0)
471  {
472  // Is the face on the boundary?
473  if (curF >= mesh_.nInternalFaces())
474  {
475  label curBPatch =
476  mesh_.boundaryMesh().whichPatch(curF);
477 
478  if (!patchFields(curBPatch))
479  {
480  patchFields.set
481  (
482  curBPatch,
484  (
485  mesh_.boundary()[curBPatch].type(),
486  mesh_.boundary()[curBPatch],
488  ::null()
489  )
490  );
491  }
492 
493  // add the face
494  label curPatchFace =
495  mesh_.boundaryMesh()
496  [curBPatch].whichFace(curF);
497 
498  patchFields[curBPatch][curPatchFace] =
499  curProcPatch[facei];
500  }
501  else
502  {
503  // Internal face
504  internalField[curF] = curProcPatch[facei];
505  }
506  }
507  }
508  }
509  }
510  }
511 
512  forAll(mesh_.boundary(), patchi)
513  {
514  // add empty patches
515  if
516  (
517  isType<emptyFvPatch>(mesh_.boundary()[patchi])
518  && !patchFields(patchi)
519  )
520  {
521  patchFields.set
522  (
523  patchi,
525  (
527  mesh_.boundary()[patchi],
529  )
530  );
531  }
532  }
533 
534 
535  // Now construct and write the field
536  // setting the internalField and patchFields
538  (
540  (
541  fieldIoObject,
542  mesh_,
543  procFields[0].dimensions(),
544  internalField,
545  patchFields
546  )
547  );
548 }
549 
550 
551 template<class Type>
554 (
555  const IOobject& fieldIoObject
556 ) const
557 {
558  // Read the field for all the processors
560  (
561  procMeshes_.size()
562  );
563 
564  forAll(procMeshes_, proci)
565  {
566  procFields.set
567  (
568  proci,
570  (
571  IOobject
572  (
573  fieldIoObject.name(),
574  procMeshes_[proci].time().timeName(),
575  procMeshes_[proci],
576  IOobject::MUST_READ,
577  IOobject::NO_WRITE
578  ),
579  procMeshes_[proci]
580  )
581  );
582  }
583 
584  return reconstructFvSurfaceField
585  (
586  IOobject
587  (
588  fieldIoObject.name(),
589  mesh_.time().timeName(),
590  mesh_,
591  IOobject::NO_READ,
592  IOobject::NO_WRITE
593  ),
594  procFields
595  );
596 }
597 
598 
599 template<class Type>
601 (
602  const IOobjectList& objects,
603  const HashSet<word>& selectedFields
604 )
605 {
606  const word& fieldClassName = DimensionedField<Type, volMesh>::typeName;
607 
608  IOobjectList fields = objects.lookupClass(fieldClassName);
609 
610  if (fields.size())
611  {
612  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
613 
614  forAllConstIter(IOobjectList, fields, fieldIter)
615  {
616  if
617  (
618  selectedFields.empty()
619  || selectedFields.found(fieldIter()->name())
620  )
621  {
622  Info<< " " << fieldIter()->name() << endl;
623 
624  reconstructFvVolumeInternalField<Type>(*fieldIter())().write();
625 
626  nReconstructed_++;
627  }
628  }
629  Info<< endl;
630  }
631 }
632 
633 
634 template<class Type>
636 (
637  const IOobjectList& objects,
638  const HashSet<word>& selectedFields
639 )
640 {
641  const word& fieldClassName =
643 
644  IOobjectList fields = objects.lookupClass(fieldClassName);
645 
646  if (fields.size())
647  {
648  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
649 
650  forAllConstIter(IOobjectList, fields, fieldIter)
651  {
652  if
653  (
654  selectedFields.empty()
655  || selectedFields.found(fieldIter()->name())
656  )
657  {
658  Info<< " " << fieldIter()->name() << endl;
659 
660  reconstructFvVolumeField<Type>(*fieldIter())().write();
661 
662  nReconstructed_++;
663  }
664  }
665  Info<< endl;
666  }
667 }
668 
669 
670 template<class Type>
672 (
673  const IOobjectList& objects,
674  const HashSet<word>& selectedFields
675 )
676 {
677  const word& fieldClassName =
679 
680  IOobjectList fields = objects.lookupClass(fieldClassName);
681 
682  if (fields.size())
683  {
684  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
685 
686  forAllConstIter(IOobjectList, fields, fieldIter)
687  {
688  if
689  (
690  selectedFields.empty()
691  || selectedFields.found(fieldIter()->name())
692  )
693  {
694  Info<< " " << fieldIter()->name() << endl;
695 
696  reconstructFvSurfaceField<Type>(*fieldIter())().write();
697 
698  nReconstructed_++;
699  }
700  }
701  Info<< endl;
702  }
703 }
704 
705 
706 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
Definition: IOobject.H:295
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructFvVolumeField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
Reconstruct volume field.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void reconstructFvSurfaceFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected surface fields.
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Traits class for primitives.
Definition: pTraits.H:50
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
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
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
Pre-declare SubField and related Field type.
Definition: Field.H:56
A class for handling words, derived from string.
Definition: word.H:59
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const volScalarField & cp
const Mesh & mesh() const
Return mesh.
Mapper for sizing only - does not do any actual mapping.
const Field< Type > & field() const
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:192
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:381
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:70
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.
tmp< DimensionedField< Type, volMesh > > reconstructFvVolumeInternalField(const IOobject &fieldIoObject, const PtrList< DimensionedField< Type, volMesh >> &procFields) const
Reconstruct volume internal field.
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:53
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructFvSurfaceField(const IOobject &fieldIoObject, const PtrList< GeometricField< Type, fvsPatchField, surfaceMesh >> &) const
Reconstruct surface field.