fvFieldReconstructorTemplates.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018-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 "fvFieldReconstructor.H"
30 #include "Time.H"
31 #include "PtrList.H"
32 #include "fvPatchFields.H"
33 #include "emptyFvPatch.H"
34 #include "emptyFvPatchField.H"
36 
37 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
38 
39 template<class Type>
42 (
43  const IOobject& fieldObject,
44  const PtrList<DimensionedField<Type, volMesh>>& procFields
45 ) const
46 {
47  // Create the internalField
48  Field<Type> internalField(mesh_.nCells());
49 
50  forAll(procMeshes_, proci)
51  {
52  const DimensionedField<Type, volMesh>& procField = procFields[proci];
53 
54  // Set the cell values in the reconstructed field
55  internalField.rmap
56  (
57  procField.field(),
58  cellProcAddressing_[proci]
59  );
60  }
61 
63  (
64  fieldObject,
65  mesh_,
66  procFields[0].dimensions(),
67  internalField
68  );
69 
70  tfield.ref().oriented() = procFields[0].oriented();
71 
72  return tfield;
73 }
74 
75 
76 template<class Type>
79 (
80  const IOobject& fieldObject,
82 ) const
83 {
84  // Create the internalField
85  Field<Type> internalField(mesh_.nCells());
86 
87  // Create the patch fields
88  PtrList<fvPatchField<Type>> patchFields(mesh_.boundary().size());
89 
90  forAll(procFields, proci)
91  {
92  const auto& procField = procFields[proci];
93 
94  // Set the cell values in the reconstructed field
95  internalField.rmap
96  (
97  procField.primitiveField(),
98  cellProcAddressing_[proci]
99  );
100 
101  // Set the boundary patch values in the reconstructed field
102  forAll(boundaryProcAddressing_[proci], patchi)
103  {
104  // Get patch index of the original patch
105  const label curBPatch = boundaryProcAddressing_[proci][patchi];
106 
107  // Get addressing slice for this patch
108  const labelList::subList cp =
109  procField.mesh().boundary()[patchi].patchSlice
110  (
111  faceProcAddressing_[proci]
112  );
113 
114  // check if the boundary patch is not a processor patch
115  if (curBPatch >= 0)
116  {
117  // Regular patch. Fast looping
118 
119  if (!patchFields.set(curBPatch))
120  {
121  patchFields.set
122  (
123  curBPatch,
125  (
126  procField.boundaryField()[patchi],
127  mesh_.boundary()[curBPatch],
129  fvPatchFieldReconstructor
130  (
131  mesh_.boundary()[curBPatch].size()
132  )
133  )
134  );
135  }
136 
137  const label curPatchStart =
138  mesh_.boundaryMesh()[curBPatch].start();
139 
140  labelList reverseAddressing(cp.size());
141 
142  forAll(cp, facei)
143  {
144  // Check
145  if (cp[facei] <= 0)
146  {
148  << "Processor " << proci
149  << " patch "
150  << procField.mesh().boundary()[patchi].name()
151  << " face " << facei
152  << " originates from reversed face since "
153  << cp[facei]
154  << exit(FatalError);
155  }
156 
157  // Subtract one to take into account offsets for
158  // face direction.
159  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
160  }
161 
162 
163  patchFields[curBPatch].rmap
164  (
165  procField.boundaryField()[patchi],
166  reverseAddressing
167  );
168  }
169  else
170  {
171  const Field<Type>& curProcPatch =
172  procField.boundaryField()[patchi];
173 
174  // In processor patches, there's a mix of internal faces (some
175  // of them turned) and possible cyclics. Slow loop
176  forAll(cp, facei)
177  {
178  // Subtract one to take into account offsets for
179  // face direction.
180  label curF = cp[facei] - 1;
181 
182  // Is the face on the boundary?
183  if (curF >= mesh_.nInternalFaces())
184  {
185  label curBPatch = mesh_.boundaryMesh().whichPatch(curF);
186 
187  if (!patchFields.set(curBPatch))
188  {
189  patchFields.set
190  (
191  curBPatch,
193  (
194  mesh_.boundary()[curBPatch].type(),
195  mesh_.boundary()[curBPatch],
197  )
198  );
199  }
200 
201  // add the face
202  label curPatchFace =
203  mesh_.boundaryMesh()
204  [curBPatch].whichFace(curF);
205 
206  patchFields[curBPatch][curPatchFace] =
207  curProcPatch[facei];
208  }
209  }
210  }
211  }
212  }
213 
214  forAll(mesh_.boundary(), patchi)
215  {
216  // add empty patches
217  if
218  (
219  isType<emptyFvPatch>(mesh_.boundary()[patchi])
220  && !patchFields.set(patchi)
221  )
222  {
223  patchFields.set
224  (
225  patchi,
227  (
229  mesh_.boundary()[patchi],
231  )
232  );
233  }
234  }
235 
236 
237  // Now construct and write the field
238  // setting the internalField and patchFields
239  auto tfield = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
240  (
241  fieldObject,
242  mesh_,
243  procFields[0].dimensions(),
244  internalField,
245  patchFields
246  );
247 
248  tfield.ref().oriented() = procFields[0].oriented();
249 
250  return tfield;
251 }
252 
253 
254 template<class Type>
257 (
258  const IOobject& fieldObject,
260 ) const
261 {
262  // Create the internalField
263  Field<Type> internalField(mesh_.nInternalFaces());
264 
265  // Create the patch fields
266  PtrList<fvsPatchField<Type>> patchFields(mesh_.boundary().size());
267 
268 
269  forAll(procMeshes_, proci)
270  {
271  const auto& procField = procFields[proci];
272 
273  // Set the face values in the reconstructed field
274 
275  // It is necessary to create a copy of the addressing array to
276  // take care of the face direction offset trick.
277  //
278  {
279  const labelList& faceMap = faceProcAddressing_[proci];
280 
281  // Correctly oriented copy of internal field
282  Field<Type> procInternalField(procField.primitiveField());
283  // Addressing into original field
284  labelList curAddr(procInternalField.size());
285 
286  forAll(procInternalField, addrI)
287  {
288  curAddr[addrI] = mag(faceMap[addrI])-1;
289  if (faceMap[addrI] < 0)
290  {
291  procInternalField[addrI] = -procInternalField[addrI];
292  }
293  }
294 
295  // Map
296  internalField.rmap(procInternalField, curAddr);
297  }
298 
299  // Set the boundary patch values in the reconstructed field
300  forAll(boundaryProcAddressing_[proci], patchi)
301  {
302  // Get patch index of the original patch
303  const label curBPatch = boundaryProcAddressing_[proci][patchi];
304 
305  // Get addressing slice for this patch
306  const labelList::subList cp =
307  procMeshes_[proci].boundary()[patchi].patchSlice
308  (
309  faceProcAddressing_[proci]
310  );
311 
312  // check if the boundary patch is not a processor patch
313  if (curBPatch >= 0)
314  {
315  // Regular patch. Fast looping
316 
317  if (!patchFields.set(curBPatch))
318  {
319  patchFields.set
320  (
321  curBPatch,
323  (
324  procField.boundaryField()[patchi],
325  mesh_.boundary()[curBPatch],
327  fvPatchFieldReconstructor
328  (
329  mesh_.boundary()[curBPatch].size()
330  )
331  )
332  );
333  }
334 
335  const label curPatchStart =
336  mesh_.boundaryMesh()[curBPatch].start();
337 
338  labelList reverseAddressing(cp.size());
339 
340  forAll(cp, facei)
341  {
342  // Subtract one to take into account offsets for
343  // face direction.
344  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
345  }
346 
347  patchFields[curBPatch].rmap
348  (
349  procField.boundaryField()[patchi],
350  reverseAddressing
351  );
352  }
353  else
354  {
355  const Field<Type>& curProcPatch =
356  procField.boundaryField()[patchi];
357 
358  // In processor patches, there's a mix of internal faces (some
359  // of them turned) and possible cyclics. Slow loop
360  forAll(cp, facei)
361  {
362  label curF = cp[facei] - 1;
363 
364  // Is the face turned the right side round
365  if (curF >= 0)
366  {
367  // Is the face on the boundary?
368  if (curF >= mesh_.nInternalFaces())
369  {
370  label curBPatch =
371  mesh_.boundaryMesh().whichPatch(curF);
372 
373  if (!patchFields.set(curBPatch))
374  {
375  patchFields.set
376  (
377  curBPatch,
379  (
380  mesh_.boundary()[curBPatch].type(),
381  mesh_.boundary()[curBPatch],
383  ::null()
384  )
385  );
386  }
387 
388  // add the face
389  label curPatchFace =
390  mesh_.boundaryMesh()
391  [curBPatch].whichFace(curF);
392 
393  patchFields[curBPatch][curPatchFace] =
394  curProcPatch[facei];
395  }
396  else
397  {
398  // Internal face
399  internalField[curF] = curProcPatch[facei];
400  }
401  }
402  }
403  }
404  }
405  }
406 
407  forAll(mesh_.boundary(), patchi)
408  {
409  // add empty patches
410  if
411  (
412  isType<emptyFvPatch>(mesh_.boundary()[patchi])
413  && !patchFields.set(patchi)
414  )
415  {
416  patchFields.set
417  (
418  patchi,
420  (
422  mesh_.boundary()[patchi],
424  )
425  );
426  }
427  }
428 
429 
430  // Now construct and write the field
431  // setting the internalField and patchFields
432  auto tfield = tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
433  (
434  fieldObject,
435  mesh_,
436  procFields[0].dimensions(),
437  internalField,
438  patchFields
439  );
440 
441  tfield.ref().oriented() = procFields[0].oriented();
442 
443  return tfield;
444 }
445 
446 
447 template<class Type>
450 (
451  const IOobject& fieldObject
452 ) const
453 {
454  // Read the field for all the processors
456  (
457  procMeshes_.size()
458  );
459 
460  forAll(procMeshes_, proci)
461  {
462  procFields.set
463  (
464  proci,
466  (
467  IOobject
468  (
469  fieldObject.name(),
470  procMeshes_[proci].thisDb().time().timeName(),
471  procMeshes_[proci].thisDb(),
474  ),
475  procMeshes_[proci]
476  )
477  );
478  }
479 
480  return reconstructField
481  (
482  IOobject
483  (
484  fieldObject.name(),
485  mesh_.thisDb().time().timeName(),
486  mesh_.thisDb(),
489  ),
490  procFields
491  );
492 }
493 
494 
495 template<class Type>
498 (
499  const IOobject& fieldObject
500 ) const
501 {
502  // Read the field for all the processors
504  (
505  procMeshes_.size()
506  );
507 
508  forAll(procMeshes_, proci)
509  {
510  procFields.set
511  (
512  proci,
514  (
515  IOobject
516  (
517  fieldObject.name(),
518  procMeshes_[proci].thisDb().time().timeName(),
519  procMeshes_[proci].thisDb(),
522  ),
523  procMeshes_[proci]
524  )
525  );
526  }
527 
528  return reconstructField
529  (
530  IOobject
531  (
532  fieldObject.name(),
533  mesh_.thisDb().time().timeName(),
534  mesh_.thisDb(),
537  ),
538  procFields
539  );
540 }
541 
542 
543 template<class Type>
546 (
547  const IOobject& fieldObject
548 ) const
549 {
550  // Read the field for all the processors
552  (
553  procMeshes_.size()
554  );
555 
556  forAll(procMeshes_, proci)
557  {
558  procFields.set
559  (
560  proci,
562  (
563  IOobject
564  (
565  fieldObject.name(),
566  procMeshes_[proci].thisDb().time().timeName(),
567  procMeshes_[proci].thisDb(),
570  ),
571  procMeshes_[proci]
572  )
573  );
574  }
575 
576  return reconstructField
577  (
578  IOobject
579  (
580  fieldObject.name(),
581  mesh_.thisDb().time().timeName(),
582  mesh_.thisDb(),
585  ),
586  procFields
587  );
588 }
589 
590 
591 template<class Type>
593 (
594  const UPtrList<const IOobject>& fieldObjects
595 )
596 {
597  typedef DimensionedField<Type, volMesh> fieldType;
598 
599  label nFields = 0;
600 
601  for (const IOobject& io : fieldObjects)
602  {
603  if (io.isHeaderClass<fieldType>())
604  {
605  if (verbose_)
606  {
607  if (!nFields)
608  {
609  Info<< " Reconstructing "
610  << fieldType::typeName << "s\n" << nl;
611  }
612  Info<< " " << io.name() << endl;
613  }
614  ++nFields;
615 
616  reconstructInternalField<Type>(io)().write();
617  ++nReconstructed_;
618  }
619  }
620 
621  if (verbose_ && nFields) Info<< endl;
622  return nFields;
623 }
624 
625 
626 template<class Type>
628 (
629  const UPtrList<const IOobject>& fieldObjects
630 )
631 {
633 
634  label nFields = 0;
635 
636  for (const IOobject& io : fieldObjects)
637  {
638  if (io.isHeaderClass<fieldType>())
639  {
640  if (verbose_)
641  {
642  if (!nFields)
643  {
644  Info<< " Reconstructing "
645  << fieldType::typeName << "s\n" << nl;
646  }
647  Info<< " " << io.name() << endl;
648  }
649  ++nFields;
650 
651  reconstructVolumeField<Type>(io)().write();
652  ++nReconstructed_;
653  }
654  }
655 
656  if (verbose_ && nFields) Info<< endl;
657  return nFields;
658 }
659 
660 
661 template<class Type>
663 (
664  const UPtrList<const IOobject>& fieldObjects
665 )
666 {
668 
669  label nFields = 0;
670 
671  for (const IOobject& io : fieldObjects)
672  {
673  if (io.isHeaderClass<fieldType>())
674  {
675  if (verbose_)
676  {
677  if (!nFields)
678  {
679  Info<< " Reconstructing "
680  << fieldType::typeName << "s\n" << nl;
681  }
682  Info<< " " << io.name() << endl;
683  }
684  ++nFields;
685 
686 
687  reconstructSurfaceField<Type>(io)().write();
688  ++nReconstructed_;
689  }
690  }
691 
692  if (verbose_ && nFields) Info<< endl;
693  return nFields;
694 }
695 
696 
697 template<class Type>
699 (
700  const IOobjectList& objects,
701  const wordRes& selectedFields
702 )
703 {
704  typedef DimensionedField<Type, volMesh> fieldType;
705 
706  return reconstructInternalFields<Type>
707  (
708  (
709  selectedFields.empty()
710  ? objects.csorted<fieldType>()
711  : objects.csorted<fieldType>(selectedFields)
712  )
713  );
714 }
715 
716 
717 template<class Type>
719 (
720  const IOobjectList& objects,
721  const wordRes& selectedFields
722 )
723 {
725 
726  return reconstructVolumeFields<Type>
727  (
728  (
729  selectedFields.empty()
730  ? objects.csorted<fieldType>()
731  : objects.csorted<fieldType>(selectedFields)
732  )
733  );
734 }
735 
736 
737 template<class Type>
739 (
740  const IOobjectList& objects,
741  const wordRes& selectedFields
742 )
743 {
745 
746  return reconstructSurfaceFields<Type>
747  (
748  (
749  selectedFields.empty()
750  ? objects.csorted<fieldType>()
751  : objects.csorted<fieldType>(selectedFields)
752  )
753  );
754 }
755 
756 
757 // ************************************************************************* //
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > reconstructVolumeField(const IOobject &fieldObject) const
Read and reconstruct volume field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > reconstructSurfaceField(const IOobject &fieldObject) const
Read and reconstruct surface field.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const char *const typeName
Typename for Field.
Definition: Field.H:86
A List obtained as a section of another List.
Definition: SubList.H:50
tmp< DimensionedField< Type, volMesh > > reconstructInternalField(const IOobject &fieldObject) const
Read and reconstruct volume internal field.
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
Generic templated field type.
Definition: Field.H:62
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
label reconstructVolumeFields(const UPtrList< const IOobject > &fieldObjects)
Read, reconstruct and write specified volume fields.
label reconstructInternalFields(const UPtrList< const IOobject > &fieldObjects)
Read, reconstruct and write specified volume internal fields.
const Field< Type > & field() const noexcept
Return const-reference to the field values.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:529
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
tmp< DimensionedField< Type, volMesh > > reconstructField(const IOobject &fieldObject, const PtrList< DimensionedField< Type, volMesh >> &procFields) const
Reconstruct volume internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label reconstructSurfaceFields(const UPtrList< const IOobject > &fieldObjects)
Read, reconstruct and write specified surface fields.
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
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.
bool isHeaderClass() const
Check if headerClassName() equals Type::typeName.
Definition: IOobjectI.H:258