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 "emptyFvPatch.H"
33 #include "fvPatchFields.H"
34 #include "fvsPatchFields.H"
35 
36 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37 
38 template<class Type>
41 (
42  const IOobject& fieldObject,
43  const PtrList<DimensionedField<Type, volMesh>>& procFields
44 ) const
45 {
46  // Create the internalField
47  Field<Type> internalField(mesh_.nCells());
48 
49  forAll(procMeshes_, proci)
50  {
51  const DimensionedField<Type, volMesh>& procField = procFields[proci];
52 
53  // Set the cell values in the reconstructed field
54  internalField.rmap
55  (
56  procField.field(),
57  cellProcAddressing_[proci]
58  );
59  }
60 
62  (
63  fieldObject,
64  mesh_,
65  procFields[0].dimensions(),
66  internalField
67  );
68 
69  tfield.ref().oriented() = procFields[0].oriented();
70 
71  return tfield;
72 }
73 
74 
75 template<class Type>
78 (
79  const IOobject& fieldObject,
81 ) const
82 {
83  // Create the internalField
84  Field<Type> internalField(mesh_.nCells());
85 
86  // Create the patch fields
87  PtrList<fvPatchField<Type>> patchFields(mesh_.boundary().size());
88 
89  forAll(procFields, proci)
90  {
91  const auto& procField = procFields[proci];
92 
93  // Set the cell values in the reconstructed field
94  internalField.rmap
95  (
96  procField.primitiveField(),
97  cellProcAddressing_[proci]
98  );
99 
100  // Set the boundary patch values in the reconstructed field
101  forAll(boundaryProcAddressing_[proci], patchi)
102  {
103  // Get patch index of the original patch
104  const label curBPatch = boundaryProcAddressing_[proci][patchi];
105 
106  // Get addressing slice for this patch
107  const labelList::subList cp =
108  procField.mesh().boundary()[patchi].patchSlice
109  (
110  faceProcAddressing_[proci]
111  );
112 
113  // check if the boundary patch is not a processor patch
114  if (curBPatch >= 0)
115  {
116  // Regular patch. Fast looping
117 
118  if (!patchFields.set(curBPatch))
119  {
120  patchFields.set
121  (
122  curBPatch,
124  (
125  procField.boundaryField()[patchi],
126  mesh_.boundary()[curBPatch],
128  fvPatchFieldReconstructor
129  (
130  mesh_.boundary()[curBPatch].size()
131  )
132  )
133  );
134  }
135 
136  const label curPatchStart =
137  mesh_.boundaryMesh()[curBPatch].start();
138 
139  labelList reverseAddressing(cp.size());
140 
141  forAll(cp, facei)
142  {
143  // Check
144  if (cp[facei] <= 0)
145  {
147  << "Processor " << proci
148  << " patch "
149  << procField.mesh().boundary()[patchi].name()
150  << " face " << facei
151  << " originates from reversed face since "
152  << cp[facei]
153  << exit(FatalError);
154  }
155 
156  // Subtract one to take into account offsets for
157  // face direction.
158  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
159  }
160 
161 
162  patchFields[curBPatch].rmap
163  (
164  procField.boundaryField()[patchi],
165  reverseAddressing
166  );
167  }
168  else
169  {
170  const Field<Type>& curProcPatch =
171  procField.boundaryField()[patchi];
172 
173  // In processor patches, there's a mix of internal faces (some
174  // of them turned) and possible cyclics. Slow loop
175  forAll(cp, facei)
176  {
177  // Subtract one to take into account offsets for
178  // face direction.
179  label curF = cp[facei] - 1;
180 
181  // Is the face on the boundary?
182  if (curF >= mesh_.nInternalFaces())
183  {
184  label curBPatch = mesh_.boundaryMesh().whichPatch(curF);
185 
186  if (!patchFields.set(curBPatch))
187  {
188  patchFields.set
189  (
190  curBPatch,
192  (
193  mesh_.boundary()[curBPatch].type(),
194  mesh_.boundary()[curBPatch],
196  )
197  );
198  }
199 
200  // add the face
201  label curPatchFace =
202  mesh_.boundaryMesh()
203  [curBPatch].whichFace(curF);
204 
205  patchFields[curBPatch][curPatchFace] =
206  curProcPatch[facei];
207  }
208  }
209  }
210  }
211  }
212 
213  forAll(mesh_.boundary(), patchi)
214  {
215  // add empty patches
216  if
217  (
218  isType<emptyFvPatch>(mesh_.boundary()[patchi])
219  && !patchFields.set(patchi)
220  )
221  {
222  patchFields.set
223  (
224  patchi,
226  (
228  mesh_.boundary()[patchi],
230  )
231  );
232  }
233  }
234 
235 
236  // Now construct and write the field
237  // setting the internalField and patchFields
238  auto tfield = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
239  (
240  fieldObject,
241  mesh_,
242  procFields[0].dimensions(),
243  internalField,
244  patchFields
245  );
246 
247  tfield.ref().oriented() = procFields[0].oriented();
248 
249  return tfield;
250 }
251 
252 
253 template<class Type>
256 (
257  const IOobject& fieldObject,
259 ) const
260 {
261  // Create the internalField
262  Field<Type> internalField(mesh_.nInternalFaces());
263 
264  // Create the patch fields
265  PtrList<fvsPatchField<Type>> patchFields(mesh_.boundary().size());
266 
267 
268  forAll(procMeshes_, proci)
269  {
270  const auto& procField = procFields[proci];
271 
272  // Set the face values in the reconstructed field
273 
274  // It is necessary to create a copy of the addressing array to
275  // take care of the face direction offset trick.
276  //
277  {
278  const labelList& faceMap = faceProcAddressing_[proci];
279 
280  // Correctly oriented copy of internal field
281  Field<Type> procInternalField(procField.primitiveField());
282  // Addressing into original field
283  labelList curAddr(procInternalField.size());
284 
285  forAll(procInternalField, addrI)
286  {
287  curAddr[addrI] = mag(faceMap[addrI])-1;
288  if (faceMap[addrI] < 0)
289  {
290  procInternalField[addrI] = -procInternalField[addrI];
291  }
292  }
293 
294  // Map
295  internalField.rmap(procInternalField, curAddr);
296  }
297 
298  // Set the boundary patch values in the reconstructed field
299  forAll(boundaryProcAddressing_[proci], patchi)
300  {
301  // Get patch index of the original patch
302  const label curBPatch = boundaryProcAddressing_[proci][patchi];
303 
304  // Get addressing slice for this patch
305  const labelList::subList cp =
306  procMeshes_[proci].boundary()[patchi].patchSlice
307  (
308  faceProcAddressing_[proci]
309  );
310 
311  // check if the boundary patch is not a processor patch
312  if (curBPatch >= 0)
313  {
314  // Regular patch. Fast looping
315 
316  if (!patchFields.set(curBPatch))
317  {
318  patchFields.set
319  (
320  curBPatch,
322  (
323  procField.boundaryField()[patchi],
324  mesh_.boundary()[curBPatch],
326  fvPatchFieldReconstructor
327  (
328  mesh_.boundary()[curBPatch].size()
329  )
330  )
331  );
332  }
333 
334  const label curPatchStart =
335  mesh_.boundaryMesh()[curBPatch].start();
336 
337  labelList reverseAddressing(cp.size());
338 
339  forAll(cp, facei)
340  {
341  // Subtract one to take into account offsets for
342  // face direction.
343  reverseAddressing[facei] = cp[facei] - 1 - curPatchStart;
344  }
345 
346  patchFields[curBPatch].rmap
347  (
348  procField.boundaryField()[patchi],
349  reverseAddressing
350  );
351  }
352  else
353  {
354  const Field<Type>& curProcPatch =
355  procField.boundaryField()[patchi];
356 
357  // In processor patches, there's a mix of internal faces (some
358  // of them turned) and possible cyclics. Slow loop
359  forAll(cp, facei)
360  {
361  label curF = cp[facei] - 1;
362 
363  // Is the face turned the right side round
364  if (curF >= 0)
365  {
366  // Is the face on the boundary?
367  if (curF >= mesh_.nInternalFaces())
368  {
369  label curBPatch =
370  mesh_.boundaryMesh().whichPatch(curF);
371 
372  if (!patchFields.set(curBPatch))
373  {
374  patchFields.set
375  (
376  curBPatch,
378  (
379  mesh_.boundary()[curBPatch].type(),
380  mesh_.boundary()[curBPatch],
382  )
383  );
384  }
385 
386  // add the face
387  label curPatchFace =
388  mesh_.boundaryMesh()
389  [curBPatch].whichFace(curF);
390 
391  patchFields[curBPatch][curPatchFace] =
392  curProcPatch[facei];
393  }
394  else
395  {
396  // Internal face
397  internalField[curF] = curProcPatch[facei];
398  }
399  }
400  }
401  }
402  }
403  }
404 
405  forAll(mesh_.boundary(), patchi)
406  {
407  // add empty patches
408  if
409  (
410  isType<emptyFvPatch>(mesh_.boundary()[patchi])
411  && !patchFields.set(patchi)
412  )
413  {
414  patchFields.set
415  (
416  patchi,
418  (
420  mesh_.boundary()[patchi],
422  )
423  );
424  }
425  }
426 
427 
428  // Now construct and write the field
429  // setting the internalField and patchFields
430  auto tfield = tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
431  (
432  fieldObject,
433  mesh_,
434  procFields[0].dimensions(),
435  internalField,
436  patchFields
437  );
438 
439  tfield.ref().oriented() = procFields[0].oriented();
440 
441  return tfield;
442 }
443 
444 
445 template<class Type>
448 (
449  const IOobject& fieldObject
450 ) const
451 {
452  // Read the field for all the processors
454  (
455  procMeshes_.size()
456  );
457 
458  forAll(procMeshes_, proci)
459  {
460  procFields.set
461  (
462  proci,
464  (
465  IOobject
466  (
467  fieldObject.name(),
468  procMeshes_[proci].thisDb().time().timeName(),
469  procMeshes_[proci].thisDb(),
472  ),
473  procMeshes_[proci]
474  )
475  );
476  }
477 
478  return reconstructField
479  (
480  IOobject
481  (
482  fieldObject.name(),
483  mesh_.thisDb().time().timeName(),
484  mesh_.thisDb(),
487  ),
488  procFields
489  );
490 }
491 
492 
493 template<class Type>
496 (
497  const IOobject& fieldObject
498 ) const
499 {
500  // Read the field for all the processors
502  (
503  procMeshes_.size()
504  );
505 
506  forAll(procMeshes_, proci)
507  {
508  procFields.set
509  (
510  proci,
512  (
513  IOobject
514  (
515  fieldObject.name(),
516  procMeshes_[proci].thisDb().time().timeName(),
517  procMeshes_[proci].thisDb(),
520  ),
521  procMeshes_[proci]
522  )
523  );
524  }
525 
526  return reconstructField
527  (
528  IOobject
529  (
530  fieldObject.name(),
531  mesh_.thisDb().time().timeName(),
532  mesh_.thisDb(),
535  ),
536  procFields
537  );
538 }
539 
540 
541 template<class Type>
544 (
545  const IOobject& fieldObject
546 ) const
547 {
548  // Read the field for all the processors
550  (
551  procMeshes_.size()
552  );
553 
554  forAll(procMeshes_, proci)
555  {
556  procFields.set
557  (
558  proci,
560  (
561  IOobject
562  (
563  fieldObject.name(),
564  procMeshes_[proci].thisDb().time().timeName(),
565  procMeshes_[proci].thisDb(),
568  ),
569  procMeshes_[proci]
570  )
571  );
572  }
573 
574  return reconstructField
575  (
576  IOobject
577  (
578  fieldObject.name(),
579  mesh_.thisDb().time().timeName(),
580  mesh_.thisDb(),
583  ),
584  procFields
585  );
586 }
587 
588 
589 template<class Type>
591 (
592  const UPtrList<const IOobject>& fieldObjects
593 )
594 {
595  typedef DimensionedField<Type, volMesh> fieldType;
596 
597  label nFields = 0;
598 
599  for (const IOobject& io : fieldObjects)
600  {
601  if (io.isHeaderClass<fieldType>())
602  {
603  if (verbose_)
604  {
605  if (!nFields)
606  {
607  Info<< " Reconstructing "
608  << fieldType::typeName << "s\n" << nl;
609  }
610  Info<< " " << io.name() << endl;
611  }
612  ++nFields;
613 
614  reconstructInternalField<Type>(io)().write();
615  ++nReconstructed_;
616  }
617  }
618 
619  if (verbose_ && nFields) Info<< endl;
620  return nFields;
621 }
622 
623 
624 template<class Type>
626 (
627  const UPtrList<const IOobject>& fieldObjects
628 )
629 {
631 
632  label nFields = 0;
633 
634  for (const IOobject& io : fieldObjects)
635  {
636  if (io.isHeaderClass<fieldType>())
637  {
638  if (verbose_)
639  {
640  if (!nFields)
641  {
642  Info<< " Reconstructing "
643  << fieldType::typeName << "s\n" << nl;
644  }
645  Info<< " " << io.name() << endl;
646  }
647  ++nFields;
648 
649  reconstructVolumeField<Type>(io)().write();
650  ++nReconstructed_;
651  }
652  }
653 
654  if (verbose_ && nFields) Info<< endl;
655  return nFields;
656 }
657 
658 
659 template<class Type>
661 (
662  const UPtrList<const IOobject>& fieldObjects
663 )
664 {
666 
667  label nFields = 0;
668 
669  for (const IOobject& io : fieldObjects)
670  {
671  if (io.isHeaderClass<fieldType>())
672  {
673  if (verbose_)
674  {
675  if (!nFields)
676  {
677  Info<< " Reconstructing "
678  << fieldType::typeName << "s\n" << nl;
679  }
680  Info<< " " << io.name() << endl;
681  }
682  ++nFields;
683 
684 
685  reconstructSurfaceField<Type>(io)().write();
686  ++nReconstructed_;
687  }
688  }
689 
690  if (verbose_ && nFields) Info<< endl;
691  return nFields;
692 }
693 
694 
695 template<class Type>
697 (
698  const IOobjectList& objects,
699  const wordRes& selectedFields
700 )
701 {
702  typedef DimensionedField<Type, volMesh> fieldType;
703 
704  return reconstructInternalFields<Type>
705  (
706  (
707  selectedFields.empty()
708  ? objects.csorted<fieldType>()
709  : objects.csorted<fieldType>(selectedFields)
710  )
711  );
712 }
713 
714 
715 template<class Type>
717 (
718  const IOobjectList& objects,
719  const wordRes& selectedFields
720 )
721 {
723 
724  return reconstructVolumeFields<Type>
725  (
726  (
727  selectedFields.empty()
728  ? objects.csorted<fieldType>()
729  : objects.csorted<fieldType>(selectedFields)
730  )
731  );
732 }
733 
734 
735 template<class Type>
737 (
738  const IOobjectList& objects,
739  const wordRes& selectedFields
740 )
741 {
743 
744  return reconstructSurfaceFields<Type>
745  (
746  (
747  selectedFields.empty()
748  ? objects.csorted<fieldType>()
749  : objects.csorted<fieldType>(selectedFields)
750  )
751  );
752 }
753 
754 
755 // ************************************************************************* //
static const word & emptyType() noexcept
The type name for empty patch fields.
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...
static const word & emptyType() noexcept
The type name for empty patch fields.
Definition: fvPatchField.H:196
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
static const DimensionedField< Type, GeoMesh > & null() noexcept
Return a null DimensionedField (reference to a nullObject).
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:675
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.
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
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.
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...
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:180
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