pointFieldReconstructorReconstructFields.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-2022 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 
28 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29 
30 template<class Type>
32 Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
33 {
34  // Read the field for all the processors
35  PtrList<GeometricField<Type, pointPatchField, pointMesh>> procFields
36  (
37  procMeshes_.size()
38  );
39 
40  forAll(procMeshes_, proci)
41  {
42  procFields.set
43  (
44  proci,
45  new GeometricField<Type, pointPatchField, pointMesh>
46  (
47  IOobject
48  (
49  fieldIoObject.name(),
50  procMeshes_[proci]().time().timeName(),
51  procMeshes_[proci](),
54  ),
55  procMeshes_[proci]
56  )
57  );
58  }
59 
60 
61  // Create the internalField
62  Field<Type> internalField(completeMesh_.size());
63 
64  // Create the patch fields
65  PtrList<pointPatchField<Type>> patchFields(completeMesh_.boundary().size());
66 
67 
68  forAll(procMeshes_, proci)
69  {
70  const GeometricField<Type, pointPatchField, pointMesh>&
71  procField = procFields[proci];
72 
73  // Get processor-to-global addressing for use in rmap
74  const labelList& procToGlobalAddr = pointProcAddressing_[proci];
75 
76  // Set the cell values in the reconstructed field
77  internalField.rmap
78  (
79  procField.primitiveField(),
80  procToGlobalAddr
81  );
82 
83  // Set the boundary patch values in the reconstructed field
84  forAll(procField.boundaryField(), patchi)
85  {
86  // Get patch index of the original patch
87  const label curBPatch =
88  patchi < completeMesh_.boundary().size() ? patchi : -1;
89 
90  // check if the boundary patch is not a processor patch
91  if (curBPatch != -1)
92  {
93  if (!patchFields(curBPatch))
94  {
95  patchFields.set(
96  curBPatch,
98  (
99  procField.boundaryField()[patchi],
100  completeMesh_.boundary()[curBPatch],
102  pointPatchFieldReconstructor
103  (
104  completeMesh_.boundary()[curBPatch].size()
105  )
106  )
107  );
108  }
109 
110  patchFields[curBPatch].rmap
111  (
112  procField.boundaryField()[patchi],
113  patchPointAddressing_[proci][patchi]
114  );
115  }
116  }
117  }
118 
119  // Construct and write the field
120  // setting the internalField and patchFields
121  return tmp<GeometricField<Type, pointPatchField, pointMesh>>
122  (
123  new GeometricField<Type, pointPatchField, pointMesh>
124  (
125  IOobject
126  (
127  fieldIoObject.name(),
128  completeMesh_().time().timeName(),
129  completeMesh_(),
132  ),
133  completeMesh_,
134  procFields[0].dimensions(),
135  internalField,
136  patchFields
137  )
138  );
139 }
140 
141 
142 // Reconstruct and write all point fields
143 template<class Type>
145 (
146  const IOobjectList& objects,
147  const HashSet<word>& selectedFields
148 )
149 {
150  word fieldClassName
151  (
153  );
154 
155  IOobjectList fields = objects.lookupClass(fieldClassName);
156 
157  if (fields.size())
158  {
159  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
160 
161  forAllConstIter(IOobjectList, fields, fieldIter)
162  {
163  if
164  (
165  !selectedFields.size()
166  || selectedFields.found(fieldIter()->name())
167  )
168  {
169  Info<< " " << fieldIter()->name() << endl;
170 
171  reconstructField<Type>(*fieldIter())().write();
172  }
173  }
174 
175  Info<< endl;
176  }
177 }
178 
179 
180 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static const char *const typeName
Definition: Field.H:105
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
void reconstructFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Reconstruct and write all fields.
List< label > labelList
A List of labels.
Definition: labelList.H:56
label size() const
Return number of points.
Definition: pointMesh.H:85
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 UPtrList.
Definition: UPtrListI.H:29
label patchi
const Time & time() const
Return time.
Definition: IOobject.C:318
messageStream Info
static autoPtr< pointPatchField< Type > > New(const word &, const pointPatch &, const DimensionedField< Type, pointMesh > &)
Return a pointer to a new patchField created on freestore given.
A class for managing temporary objects.
Definition: PtrList.H:53
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:97
tmp< GeometricField< Type, pointPatchField, pointMesh > > reconstructField(const IOobject &fieldIoObject)
Reconstruct field.