pointFieldDecomposer.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 | Copyright (C) 2011-2015 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 
26 #include "pointFieldDecomposer.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const pointPatch& completeMeshPatch,
33  const pointPatch& procMeshPatch,
34  const labelList& directAddr
35 )
36 :
37  pointPatchFieldMapperPatchRef
38  (
39  completeMeshPatch,
40  procMeshPatch
41  ),
42  directAddressing_(procMeshPatch.size(), -1),
43  hasUnmapped_(false)
44 {
45  // Create the inverse-addressing of the patch point labels.
46  labelList pointMap(completeMeshPatch.boundaryMesh().mesh().size(), -1);
47 
48  const labelList& completeMeshPatchPoints = completeMeshPatch.meshPoints();
49 
50  forAll(completeMeshPatchPoints, pointi)
51  {
52  pointMap[completeMeshPatchPoints[pointi]] = pointi;
53  }
54 
55  // Use the inverse point addressing to create the addressing table for this
56  // patch
57  const labelList& procMeshPatchPoints = procMeshPatch.meshPoints();
58 
59  forAll(procMeshPatchPoints, pointi)
60  {
61  directAddressing_[pointi] =
62  pointMap[directAddr[procMeshPatchPoints[pointi]]];
63  }
64 
65  // Check that all the patch point addresses are set
66  if (directAddressing_.size() && min(directAddressing_) < 0)
67  {
68  hasUnmapped_ = true;
69 
71  << "Incomplete patch point addressing"
72  << abort(FatalError);
73  }
74 }
75 
76 
77 Foam::pointFieldDecomposer::pointFieldDecomposer
78 (
79  const pointMesh& completeMesh,
80  const pointMesh& procMesh,
81  const labelList& pointAddressing,
82  const labelList& boundaryAddressing
83 )
84 :
85  completeMesh_(completeMesh),
86  procMesh_(procMesh),
87  pointAddressing_(pointAddressing),
88  boundaryAddressing_(boundaryAddressing),
89  patchFieldDecomposerPtrs_
90  (
91  procMesh_.boundary().size(),
92  static_cast<patchFieldDecomposer*>(NULL)
93  )
94 {
95  forAll(boundaryAddressing_, patchi)
96  {
97  if (boundaryAddressing_[patchi] >= 0)
98  {
99  patchFieldDecomposerPtrs_[patchi] = new patchFieldDecomposer
100  (
101  completeMesh_.boundary()[boundaryAddressing_[patchi]],
102  procMesh_.boundary()[patchi],
103  pointAddressing_
104  );
105  }
106  }
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
113 {
114  forAll(patchFieldDecomposerPtrs_, patchi)
115  {
116  if (patchFieldDecomposerPtrs_[patchi])
117  {
118  delete patchFieldDecomposerPtrs_[patchi];
119  }
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< label > labelList
A List of labels.
Definition: labelList.H:56
faceListList boundary(nPatches)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label patchi
patchFieldDecomposer(const pointPatch &completeMeshPatch, const pointPatch &procMeshPatch, const labelList &directAddr)
Construct given addressing.
~pointFieldDecomposer()
Destructor.