pointFieldReconstructor.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-2015 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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#include "pointFields.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
35
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
40(
41 const pointMesh& mesh,
42 const PtrList<pointMesh>& procMeshes,
43 const PtrList<labelIOList>& pointProcAddressing,
44 const PtrList<labelIOList>& boundaryProcAddressing
45)
46:
47 mesh_(mesh),
48 procMeshes_(procMeshes),
49 pointProcAddressing_(pointProcAddressing),
50 boundaryProcAddressing_(boundaryProcAddressing),
51 patchPointAddressing_(procMeshes.size()),
52 nReconstructed_(0)
53{
54 // Inverse-addressing of the patch point labels.
55 labelList pointMap(mesh_.size(), -1);
56
57 // Create the pointPatch addressing
58 forAll(procMeshes_, proci)
59 {
60 const pointMesh& procMesh = procMeshes_[proci];
61
62 patchPointAddressing_[proci].setSize(procMesh.boundary().size());
63
64 forAll(procMesh.boundary(), patchi)
65 {
66 if (boundaryProcAddressing_[proci][patchi] >= 0)
67 {
68 labelList& procPatchAddr = patchPointAddressing_[proci][patchi];
69 procPatchAddr.setSize(procMesh.boundary()[patchi].size(), -1);
70
71 const labelList& patchPointLabels =
72 mesh_.boundary()[boundaryProcAddressing_[proci][patchi]]
73 .meshPoints();
74
75 // Create the inverse-addressing of the patch point labels.
76 forAll(patchPointLabels, pointi)
77 {
78 pointMap[patchPointLabels[pointi]] = pointi;
79 }
80
81 const labelList& procPatchPoints =
82 procMesh.boundary()[patchi].meshPoints();
83
84 forAll(procPatchPoints, pointi)
85 {
86 procPatchAddr[pointi] =
87 pointMap
88 [
89 pointProcAddressing_[proci][procPatchPoints[pointi]]
90 ];
91 }
92
93 if (procPatchAddr.size() && min(procPatchAddr) < 0)
94 {
96 << "Incomplete patch point addressing"
97 << abort(FatalError);
98 }
99 }
100 }
101 }
102}
103
104
105// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106
108(
109 const IOobjectList& objects,
110 const wordRes& selected
111)
112{
113 label nTotal = 0;
114
115 nTotal += reconstructPointFields<scalar>(objects, selected);
116 nTotal += reconstructPointFields<vector>(objects, selected);
117 nTotal += reconstructPointFields<symmTensor>(objects, selected);
118 nTotal += reconstructPointFields<sphericalTensor>(objects, selected);
119 nTotal += reconstructPointFields<tensor>(objects, selected);
120
121 return nTotal;
122}
123
124
125// ************************************************************************* //
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Point field reconstructor.
static int verbose_
Output verbosity when writing.
label reconstructAllFields(const IOobjectList &objects, const wordRes &selectedFields=wordRes())
Reconstruct and write all known field types.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
static label size(const Mesh &mesh)
Return size. Number of points.
Definition: pointMesh.H:102
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:114
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333