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-2025 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
61 template<class Type>
64 {
65  // Read the field for all the processors
66  PtrList<PointField<Type>> procFields(procMeshes_.size());
67 
68  forAll(procMeshes_, proci)
69  {
70  procFields.set
71  (
72  proci,
74  (
75  IOobject
76  (
77  fieldIoObject.name(),
78  procMeshes_[proci].time().name(),
79  procMeshes_[proci],
82  ),
83  pointMesh::New(procMeshes_[proci])
84  )
85  );
86  }
87 
88  // Create the internalField
89  Field<Type> internalField(completeMesh_.size());
90 
91  // Create the patch fields
92  PtrList<pointPatchField<Type>> patchFields(completeMesh_.boundary().size());
93 
94  forAll(procMeshes_, proci)
95  {
96  const PointField<Type>&
97  procField = procFields[proci];
98 
99  // Get processor-to-global addressing for use in rmap
100  const labelList& procToGlobalAddr = pointProcAddressing_[proci];
101 
102  // Set the cell values in the reconstructed field
103  internalField.rmap
104  (
105  procField.primitiveField(),
106  procToGlobalAddr
107  );
108 
109  // Set the boundary patch values in the reconstructed field
110  forAll(procField.boundaryField(), patchi)
111  {
112  // Get patch index of the original patch
113  const label curBPatch =
114  patchi < completeMesh_.boundary().size() ? patchi : -1;
115 
116  // check if the boundary patch is not a processor patch
117  if (curBPatch != -1)
118  {
119  if (!patchFields(curBPatch))
120  {
121  patchFields.set
122  (
123  curBPatch,
125  (
126  procField.boundaryField()[patchi],
127  completeMesh_.boundary()[curBPatch],
130  (
131  completeMesh_.boundary()[curBPatch].size()
132  )
133  )
134  );
135  }
136 
137  patchFields[curBPatch].map
138  (
139  procField.boundaryField()[patchi],
141  (
142  patchPointAddressing_[proci][patchi]
143  )
144  );
145  }
146  }
147  }
148 
149  // Construct and return the field
150  return tmp<PointField<Type>>
151  (
152  new PointField<Type>
153  (
154  IOobject
155  (
156  fieldIoObject.name(),
157  completeMesh_().time().name(),
158  completeMesh_(),
161  ),
162  completeMesh_,
163  procFields[0].dimensions(),
164  internalField,
165  patchFields
166  )
167  );
168 }
169 
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:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Pre-declare SubField and related Field type.
Definition: Field.H:83
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:386
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
A HashTable with keys but without contents.
Definition: HashSet.H:62
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:307
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
tmp< PointField< Type > > reconstructField(const IOobject &fieldIoObject)
Read and reconstruct a field.
void reconstructFields(const IOobjectList &objects, const HashSet< word > &selectedFields)
Read, reconstruct and write all/selected point fields.
Abstract base class for point-mesh patch fields.
Reverse field mapper.
Mapper which sets the field size. It does not actually map values.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
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:234
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
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:258
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
static const char nl
Definition: Ostream.H:267
objects