pointFieldReconstructorTemplates.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"
28 #include "reverseFieldMapper.H"
29 #include "setSizeFieldMapper.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class FieldType>
34 bool Foam::pointFieldReconstructor::reconstructs
35 (
36  const IOobjectList& objects,
37  const HashSet<word>& selectedFields
38 )
39 {
40  IOobjectList fields = objects.lookupClass(FieldType::typeName);
41 
42  if (fields.size() && selectedFields.empty())
43  {
44  return true;
45  }
46 
47  forAllConstIter(IOobjectList, fields, fieldIter)
48  {
49  if (selectedFields.found(fieldIter()->name()))
50  {
51  return true;
52  }
53  }
54 
55  return false;
56 }
57 
58 
59 template<class Type>
61 Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
62 {
63  // Read the field for all the processors
64  PtrList<PointField<Type>> procFields(procMeshes_.size());
65 
66  forAll(procMeshes_, proci)
67  {
68  procFields.set
69  (
70  proci,
71  new PointField<Type>
72  (
73  IOobject
74  (
75  fieldIoObject.name(),
76  procMeshes_[proci].time().name(),
77  procMeshes_[proci],
80  ),
81  pointMesh::New(procMeshes_[proci])
82  )
83  );
84  }
85 
86  // Create the internalField
87  Field<Type> internalField(completeMesh_.size());
88 
89  // Create the patch fields
90  PtrList<pointPatchField<Type>> patchFields(completeMesh_.boundary().size());
91 
92  forAll(procMeshes_, proci)
93  {
94  const PointField<Type>&
95  procField = procFields[proci];
96 
97  // Get processor-to-global addressing for use in rmap
98  const labelList& procToGlobalAddr = pointProcAddressing_[proci];
99 
100  // Set the cell values in the reconstructed field
101  internalField.rmap
102  (
103  procField.primitiveField(),
104  procToGlobalAddr
105  );
106 
107  // Set the boundary patch values in the reconstructed field
108  forAll(procField.boundaryField(), patchi)
109  {
110  // Get patch index of the original patch
111  const label curBPatch =
112  patchi < completeMesh_.boundary().size() ? patchi : -1;
113 
114  // check if the boundary patch is not a processor patch
115  if (curBPatch != -1)
116  {
117  if (!patchFields(curBPatch))
118  {
119  patchFields.set
120  (
121  curBPatch,
123  (
124  procField.boundaryField()[patchi],
125  completeMesh_.boundary()[curBPatch],
127  setSizeFieldMapper
128  (
129  completeMesh_.boundary()[curBPatch].size()
130  )
131  )
132  );
133  }
134 
135  patchFields[curBPatch].map
136  (
137  procField.boundaryField()[patchi],
138  reverseFieldMapper
139  (
140  patchPointAddressing_[proci][patchi]
141  )
142  );
143  }
144  }
145  }
146 
147  // Construct and return the field
148  return tmp<PointField<Type>>
149  (
150  new PointField<Type>
151  (
152  IOobject
153  (
154  fieldIoObject.name(),
155  completeMesh_().time().name(),
156  completeMesh_(),
159  ),
160  completeMesh_,
161  procFields[0].dimensions(),
162  internalField,
163  patchFields
164  )
165  );
166 }
167 
168 
169 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170 
171 template<class Type>
173 (
174  const IOobjectList& objects,
175  const HashSet<word>& selectedFields
176 )
177 {
178  word fieldClassName
179  (
181  );
182 
183  IOobjectList fields = objects.lookupClass(fieldClassName);
184 
185  if (fields.size())
186  {
187  Info<< nl << " Reconstructing " << fieldClassName << "s"
188  << nl << endl;
189 
190  forAllConstIter(IOobjectList, fields, fieldIter)
191  {
192  if
193  (
194  !selectedFields.size()
195  || selectedFields.found(fieldIter()->name())
196  )
197  {
198  Info<< " " << fieldIter()->name() << endl;
199 
200  reconstructField<Type>(*fieldIter())().write();
201  }
202  }
203  }
204 }
205 
206 
207 // ************************************************************************* //
#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 pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
static const char *const typeName
Definition: Field.H:106
void reconstructFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Reconstruct and write all fields.
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:228
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
static const char nl
Definition: Ostream.H:266
objects