pointFieldReconstructor.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 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const pointMesh& completeMesh,
34  const PtrList<fvMesh>& procMeshes,
35  const labelListList& pointProcAddressing
36 )
37 :
38  completeMesh_(completeMesh),
39  procMeshes_(procMeshes),
40  pointProcAddressing_(pointProcAddressing),
41  patchPointAddressing_(procMeshes.size())
42 {
43  // Inverse-addressing of the patch point labels.
44  labelList pointMap(completeMesh_.size(), -1);
45 
46  // Create the pointPatch addressing
47  forAll(procMeshes_, proci)
48  {
49  const pointMesh& procMesh = pointMesh::New(procMeshes_[proci]);
50 
51  patchPointAddressing_[proci].setSize(procMesh.boundary().size());
52 
53  forAll(procMesh.boundary(), patchi)
54  {
55  if (patchi < completeMesh_.boundary().size())
56  {
57  labelList& procPatchAddr = patchPointAddressing_[proci][patchi];
58  procPatchAddr.setSize(procMesh.boundary()[patchi].size(), -1);
59 
60  const labelList& patchPointLabels =
61  completeMesh_.boundary()[patchi].meshPoints();
62 
63  // Create the inverse-addressing of the patch point labels.
64  forAll(patchPointLabels, pointi)
65  {
66  pointMap[patchPointLabels[pointi]] = pointi;
67  }
68 
69  const labelList& procPatchPoints =
70  procMesh.boundary()[patchi].meshPoints();
71 
72  forAll(procPatchPoints, pointi)
73  {
74  procPatchAddr[pointi] =
75  pointMap
76  [
77  pointProcAddressing_[proci][procPatchPoints[pointi]]
78  ];
79  }
80 
81  if (procPatchAddr.size() && min(procPatchAddr) < 0)
82  {
84  << "Incomplete patch point addressing"
85  << abort(FatalError);
86  }
87  }
88  }
89  }
90 }
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 bool Foam::pointFieldReconstructor::reconstructs
96 (
97  const IOobjectList& objects,
98  const HashSet<word>& selectedFields
99 )
100 {
101  bool result = false;
102 
103  #define DO_POINT_FIELDS_TYPE(Type, nullArg) \
104  result = result \
105  || reconstructs<PointField<Type>>(objects, selectedFields);
106  FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
107  #undef DO_POINT_FIELDS_TYPE
108 
109  return result;
110 }
111 
112 
113 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void setSize(const label)
Reset size of List.
Definition: List.C:281
pointFieldReconstructor(const pointMesh &mesh, const PtrList< fvMesh > &procMeshes, const labelListList &pointProcAddressing)
Construct from components.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
objects