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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 Foam::pointFieldReconstructor::pointFieldReconstructor
33 (
34  const pointMesh& mesh,
35  const PtrList<pointMesh>& procMeshes,
36  const PtrList<labelIOList>& pointProcAddressing,
37  const PtrList<labelIOList>& boundaryProcAddressing
38 )
39 :
40  mesh_(mesh),
41  procMeshes_(procMeshes),
42  pointProcAddressing_(pointProcAddressing),
43  boundaryProcAddressing_(boundaryProcAddressing),
44  patchPointAddressing_(procMeshes.size()),
45  nReconstructed_(0)
46 {
47  // Inverse-addressing of the patch point labels.
48  labelList pointMap(mesh_.size(), -1);
49 
50  // Create the pointPatch addressing
51  forAll(procMeshes_, proci)
52  {
53  const pointMesh& procMesh = procMeshes_[proci];
54 
55  patchPointAddressing_[proci].setSize(procMesh.boundary().size());
56 
57  forAll(procMesh.boundary(), patchi)
58  {
59  if (boundaryProcAddressing_[proci][patchi] >= 0)
60  {
61  labelList& procPatchAddr = patchPointAddressing_[proci][patchi];
62  procPatchAddr.setSize(procMesh.boundary()[patchi].size(), -1);
63 
64  const labelList& patchPointLabels =
65  mesh_.boundary()[boundaryProcAddressing_[proci][patchi]]
66  .meshPoints();
67 
68  // Create the inverse-addressing of the patch point labels.
69  forAll(patchPointLabels, pointi)
70  {
71  pointMap[patchPointLabels[pointi]] = pointi;
72  }
73 
74  const labelList& procPatchPoints =
75  procMesh.boundary()[patchi].meshPoints();
76 
77  forAll(procPatchPoints, pointi)
78  {
79  procPatchAddr[pointi] =
80  pointMap
81  [
82  pointProcAddressing_[proci][procPatchPoints[pointi]]
83  ];
84  }
85 
86  if (procPatchAddr.size() && min(procPatchAddr) < 0)
87  {
89  << "Incomplete patch point addressing"
90  << abort(FatalError);
91  }
92  }
93  }
94  }
95 }
96 
97 
98 // ************************************************************************* //
Foam::pointMesh::boundary
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:109
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
pointFieldReconstructor.H
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::List< label >