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-2023 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 #include "fvMesh.H"
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<class Type>
34 Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
35 {
36  // Read the field for all the processors
37  PtrList<PointField<Type>> procFields
38  (
39  procMeshes_.size()
40  );
41 
42  forAll(procMeshes_, proci)
43  {
44  procFields.set
45  (
46  proci,
47  new PointField<Type>
48  (
49  IOobject
50  (
51  fieldIoObject.name(),
52  procMeshes_[proci].time().name(),
53  procMeshes_[proci],
56  ),
57  pointMesh::New(procMeshes_[proci])
58  )
59  );
60  }
61 
62 
63  // Create the internalField
64  Field<Type> internalField(completeMesh_.size());
65 
66  // Create the patch fields
67  PtrList<pointPatchField<Type>> patchFields(completeMesh_.boundary().size());
68 
69 
70  forAll(procMeshes_, proci)
71  {
72  const PointField<Type>&
73  procField = procFields[proci];
74 
75  // Get processor-to-global addressing for use in rmap
76  const labelList& procToGlobalAddr = pointProcAddressing_[proci];
77 
78  // Set the cell values in the reconstructed field
79  internalField.rmap
80  (
81  procField.primitiveField(),
82  procToGlobalAddr
83  );
84 
85  // Set the boundary patch values in the reconstructed field
86  forAll(procField.boundaryField(), patchi)
87  {
88  // Get patch index of the original patch
89  const label curBPatch =
90  patchi < completeMesh_.boundary().size() ? patchi : -1;
91 
92  // check if the boundary patch is not a processor patch
93  if (curBPatch != -1)
94  {
95  if (!patchFields(curBPatch))
96  {
97  patchFields.set
98  (
99  curBPatch,
101  (
102  procField.boundaryField()[patchi],
103  completeMesh_.boundary()[curBPatch],
105  setSizePointPatchFieldMapper
106  (
107  completeMesh_.boundary()[curBPatch].size()
108  )
109  )
110  );
111  }
112 
113  patchFields[curBPatch].map
114  (
115  procField.boundaryField()[patchi],
116  reversePointPatchFieldMapper
117  (
118  patchPointAddressing_[proci][patchi]
119  )
120  );
121  }
122  }
123  }
124 
125  // Construct and write the field
126  // setting the internalField and patchFields
127  return tmp<PointField<Type>>
128  (
129  new PointField<Type>
130  (
131  IOobject
132  (
133  fieldIoObject.name(),
134  completeMesh_().time().name(),
135  completeMesh_(),
138  ),
139  completeMesh_,
140  procFields[0].dimensions(),
141  internalField,
142  patchFields
143  )
144  );
145 }
146 
147 
148 // Reconstruct and write all point fields
149 template<class Type>
151 (
152  const IOobjectList& objects,
153  const HashSet<word>& selectedFields
154 )
155 {
156  word fieldClassName
157  (
159  );
160 
161  IOobjectList fields = objects.lookupClass(fieldClassName);
162 
163  if (fields.size())
164  {
165  Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
166 
167  forAllConstIter(IOobjectList, fields, fieldIter)
168  {
169  if
170  (
171  !selectedFields.size()
172  || selectedFields.found(fieldIter()->name())
173  )
174  {
175  Info<< " " << fieldIter()->name() << endl;
176 
177  reconstructField<Type>(*fieldIter())().write();
178 
179  nReconstructed_++;
180  }
181  }
182 
183  Info<< endl;
184  }
185 }
186 
187 
188 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
Definition: Field.H:105
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
tmp< PointField< Type > > reconstructField(const IOobject &fieldIoObject)
Reconstruct field.
void reconstructFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Reconstruct and write all fields.
label size() const
Return number of points.
Definition: pointMesh.H:92
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:104
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: tmp.H:55
label patchi
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
objects