1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
25 This function calculates the list of patch edges, defined on the list of
26 points supporting the patch. The edges are ordered:
27 - 0..nInternalEdges-1 : upper triangular order
28 - nInternalEdges.. : boundary edges (no particular order)
30 Other patch addressing information is also calculated:
31 - faceFaces with neighbour faces in ascending order
32 - edgeFaces with ascending face order
33 - faceEdges sorted according to edges of a face
35 \*---------------------------------------------------------------------------*/
37 #include "PrimitivePatch.H"
38 #include "DynamicList.H"
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 template<class> class FaceList,
51 Foam::PrimitivePatch<Face, FaceList, PointField, PointType>::
52 calcAddressing() const
56 Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
57 << "calcAddressing() : calculating patch addressing"
61 if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
63 // it is considered an error to attempt to recalculate
64 // if already allocated
67 "PrimitivePatch<Face, FaceList, PointField, PointType>::"
69 ) << "addressing already calculated"
73 // get reference to localFaces
74 const List<Face>& locFcs = localFaces();
76 // get reference to pointFaces
77 const labelListList& pf = pointFaces();
79 // Guess the max number of edges and neighbours for a face
83 maxEdges += locFcs[faceI].size();
86 // create the lists for the various results. (resized on completion)
87 edgesPtr_ = new edgeList(maxEdges);
88 edgeList& edges = *edgesPtr_;
90 edgeFacesPtr_ = new labelListList(maxEdges);
91 labelListList& edgeFaces = *edgeFacesPtr_;
93 // faceFaces created using a dynamic list. Cannot guess size because
94 // of multiple connections
95 List<DynamicList<label> > ff(locFcs.size());
97 faceEdgesPtr_ = new labelListList(locFcs.size());
98 labelListList& faceEdges = *faceEdgesPtr_;
100 // count the number of face neighbours
101 labelList noFaceFaces(locFcs.size());
103 // initialise the lists of subshapes for each face to avoid duplication
104 edgeListList faceIntoEdges(locFcs.size());
106 forAll(locFcs, faceI)
108 faceIntoEdges[faceI] = locFcs[faceI].edges();
110 labelList& curFaceEdges = faceEdges[faceI];
111 curFaceEdges.setSize(faceIntoEdges[faceI].size());
113 forAll(curFaceEdges, faceEdgeI)
115 curFaceEdges[faceEdgeI] = -1;
119 // This algorithm will produce a separated list of edges, internal edges
120 // starting from 0 and boundary edges starting from the top and
127 // Note that faceIntoEdges is sorted acc. to local vertex numbering
128 // in face (i.e. curEdges[0] is edge between f[0] and f[1])
130 // For all local faces ...
131 forAll(locFcs, faceI)
133 // Get reference to vertices of current face and corresponding edges.
134 const Face& curF = locFcs[faceI];
135 const edgeList& curEdges = faceIntoEdges[faceI];
137 // Record the neighbour face. Multiple connectivity allowed
138 List<DynamicList<label> > neiFaces(curF.size());
139 List<DynamicList<label> > edgeOfNeiFace(curF.size());
141 label nNeighbours = 0;
144 forAll(curEdges, edgeI)
146 // If the edge is already detected, skip
147 if (faceEdges[faceI][edgeI] >= 0) continue;
151 // Set reference to the current edge
152 const edge& e = curEdges[edgeI];
154 // Collect neighbours for the current face vertex.
156 const labelList& nbrFaces = pf[e.start()];
158 forAll(nbrFaces, nbrFaceI)
160 // set reference to the current neighbour
161 label curNei = nbrFaces[nbrFaceI];
163 // Reject neighbours with the lower label
166 // get the reference to subshapes of the neighbour
167 const edgeList& searchEdges = faceIntoEdges[curNei];
169 forAll(searchEdges, neiEdgeI)
171 if (searchEdges[neiEdgeI] == e)
176 neiFaces[edgeI].append(curNei);
177 edgeOfNeiFace[edgeI].append(neiEdgeI);
179 // Record faceFaces both ways
180 ff[faceI].append(curNei);
181 ff[curNei].append(faceI);
183 // Keep searching due to multiple connectivity
187 } // End of neighbouring faces
191 // Register another detected internal edge
194 } // End of current edges
196 // Add the edges in increasing number of neighbours.
197 // Note: for multiply connected surfaces, the lower index neighbour for
198 // an edge will come first.
200 // Add the faces in the increasing order of neighbours
201 for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
203 // Find the lowest neighbour which is still valid
205 label minNei = locFcs.size();
207 forAll(neiFaces, nfI)
209 if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
212 minNei = neiFaces[nfI][0];
218 // Add the face to the list of faces
219 edges[nEdges] = curEdges[nextNei];
221 // Set face-edge and face-neighbour-edge to current face label
222 faceEdges[faceI][nextNei] = nEdges;
224 DynamicList<label>& cnf = neiFaces[nextNei];
225 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
227 // Set edge-face addressing
228 labelList& curEf = edgeFaces[nEdges];
229 curEf.setSize(cnf.size() + 1);
234 faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
236 curEf[cnfI + 1] = cnf[cnfI];
239 // Stop the neighbour from being used again
243 // Increment number of faces counter
250 "PrimitivePatch<Face, FaceList, PointField, PointType>::"
252 ) << "Error in internal edge insertion"
253 << abort(FatalError);
258 nInternalEdges_ = nEdges;
262 forAll(faceEdges, faceI)
264 labelList& curEdges = faceEdges[faceI];
266 forAll(curEdges, edgeI)
268 if (curEdges[edgeI] < 0)
270 // Grab edge and faceEdge
271 edges[nEdges] = faceIntoEdges[faceI][edgeI];
272 curEdges[edgeI] = nEdges;
275 labelList& curEf = edgeFaces[nEdges];
285 edges.setSize(nEdges);
288 edgeFaces.setSize(nEdges);
291 faceFacesPtr_ = new labelListList(locFcs.size());
292 labelListList& faceFaces = *faceFacesPtr_;
294 forAll(faceFaces, faceI)
296 faceFaces[faceI].transfer(ff[faceI]);
302 Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
303 << "calcAddressing() : finished calculating patch addressing"
309 // ************************************************************************* //