pointFieldReconstructorTemplates.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 
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
36 (
37  const IOobject& fieldObject,
39 ) const
40 {
42 
43  // Create the internalField
44  Field<Type> internalField(mesh_.size());
45 
46  // Create the patch fields
47  PtrList<pointPatchField<Type>> patchFields(mesh_.boundary().size());
48 
49 
50  forAll(procMeshes_, proci)
51  {
52  const auto& procField = procFields[proci];
53 
54  // Get processor-to-global addressing for use in rmap
55  const labelList& procToGlobalAddr = pointProcAddressing_[proci];
56 
57  // Set the cell values in the reconstructed field
58  internalField.rmap
59  (
60  procField.primitiveField(),
61  procToGlobalAddr
62  );
63 
64  // Set the boundary patch values in the reconstructed field
65  forAll(boundaryProcAddressing_[proci], patchi)
66  {
67  // Get patch index of the original patch
68  const label curBPatch = boundaryProcAddressing_[proci][patchi];
69 
70  // check if the boundary patch is not a processor patch
71  if (curBPatch >= 0)
72  {
73  if (!patchFields.set(curBPatch))
74  {
75  patchFields.set
76  (
77  curBPatch,
79  (
80  procField.boundaryField()[patchi],
81  mesh_.boundary()[curBPatch],
84  (
85  mesh_.boundary()[curBPatch].size()
86  )
87  )
88  );
89  }
90 
91  patchFields[curBPatch].rmap
92  (
93  procField.boundaryField()[patchi],
94  patchPointAddressing_[proci][patchi]
95  );
96  }
97  }
98  }
99 
100  // Construct and write the field
101  // setting the internalField and patchFields
102  return tmp<fieldType>::New
103  (
104  IOobject
105  (
106  fieldObject.name(),
107  mesh_.thisDb().time().timeName(),
108  mesh_.thisDb(),
111  ),
112  mesh_,
113  procFields[0].dimensions(),
114  internalField,
115  patchFields
116  );
117 }
118 
119 
120 template<class Type>
123 (
124  const IOobject& fieldObject
125 )
126 {
128 
129  // Read the field for all the processors
130  PtrList<fieldType> procFields(procMeshes_.size());
131 
132  forAll(procMeshes_, proci)
133  {
134  procFields.set
135  (
136  proci,
137  new fieldType
138  (
139  IOobject
140  (
141  fieldObject.name(),
142  procMeshes_[proci].thisDb().time().timeName(),
143  procMeshes_[proci].thisDb(),
146  ),
147  procMeshes_[proci]
148  )
149  );
150  }
151 
152  return reconstructField
153  (
154  IOobject
155  (
156  fieldObject.name(),
157  mesh_.thisDb().time().timeName(),
158  mesh_.thisDb(),
161  ),
162  procFields
163  );
164 }
165 
166 
167 template<class Type>
169 (
170  const UPtrList<const IOobject>& fieldObjects
171 )
172 {
174 
175  label nFields = 0;
176 
177  for (const IOobject& io : fieldObjects)
178  {
179  if (io.isHeaderClass<fieldType>())
180  {
181  if (verbose_)
182  {
183  if (!nFields)
184  {
185  Info<< " Reconstructing "
186  << fieldType::typeName << "s\n" << nl;
187  }
188  Info<< " " << io.name() << endl;
189  }
190  ++nFields;
191 
192  reconstructPointField<Type>(io)().write();
193  ++nReconstructed_;
194  }
195  }
196 
197  if (verbose_ && nFields) Info<< endl;
198  return nFields;
199 }
200 
201 
202 template<class Type>
204 (
205  const IOobjectList& objects,
206  const wordRes& selectedFields
207 )
208 {
210 
211  return reconstructPointFields<Type>
212  (
213  (
214  selectedFields.empty()
215  ? objects.csorted<fieldType>()
216  : objects.csorted<fieldType>(selectedFields)
217  )
218  );
219 }
220 
221 
222 // ************************************************************************* //
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
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:666
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static label size(const Mesh &mesh)
Return size. Number of points.
Definition: pointMesh.H:114
Ignore writing from objectRegistry::writeObject()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< GeometricField< Type, pointPatchField, pointMesh > > reconstructField(const IOobject &fieldObject, const PtrList< GeometricField< Type, pointPatchField, pointMesh >> &) const
Reconstruct field.
Abstract base class for point-mesh patch fields.
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
Generic templated field type.
Definition: Field.H:62
const Time & time() const noexcept
Return time registry.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:154
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
tmp< GeometricField< Type, pointPatchField, pointMesh > > reconstructPointField(const IOobject &fieldObject)
Read and reconstruct point field.
label reconstructPointFields(const UPtrList< const IOobject > &fieldObjects)
Reconstruct and write specified point fields.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
const pointBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: pointMesh.H:130
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
UPtrList< const IOobject > csorted() const
The sorted list of IOobjects with headerClassName == Type::typeName.
bool isHeaderClass() const
Check if headerClassName() equals Type::typeName.
Definition: IOobjectI.H:258