1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
26 This function calculates the list of patch edges, defined on the list of
27 points supporting the patch. The edges are ordered:
28 - 0..nInternalEdges-1 : upper triangular order
29 - nInternalEdges.. : boundary edges (no particular order)
31 Other patch addressing information is also calculated:
32 - faceFaces with neighbour faces in ascending order
33 - edgeFaces with ascending face order
34 - faceEdges sorted according to edges of a face
36 \*---------------------------------------------------------------------------*/
38 #include "PrimitivePatch.H"
39 #include "DynamicList.H"
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 template<class> class FaceList,
52 Foam::PrimitivePatch<Face, FaceList, PointField, PointType>::
53 calcAddressing() const
57 Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
58 << "calcAddressing() : calculating patch addressing"
62 if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
64 // it is considered an error to attempt to recalculate
65 // if already allocated
68 "PrimitivePatch<Face, FaceList, PointField, PointType>::"
70 ) << "addressing already calculated"
74 // get reference to localFaces
75 const List<Face>& locFcs = localFaces();
77 // get reference to pointFaces
78 const labelListList& pf = pointFaces();
80 // Guess the max number of edges and neighbours for a face
82 forAll (locFcs, faceI)
84 maxEdges += locFcs[faceI].size();
87 // create the lists for the various results. (resized on completion)
88 edgesPtr_ = new edgeList(maxEdges);
89 edgeList& edges = *edgesPtr_;
91 edgeFacesPtr_ = new labelListList(maxEdges);
92 labelListList& edgeFaces = *edgeFacesPtr_;
94 // faceFaces created using a dynamic list. Cannot guess size because
95 // of multiple connections
96 List<DynamicList<label> > ff(locFcs.size());
98 faceEdgesPtr_ = new labelListList(locFcs.size());
99 labelListList& faceEdges = *faceEdgesPtr_;
101 // count the number of face neighbours
102 labelList noFaceFaces(locFcs.size());
104 // initialise the lists of subshapes for each face to avoid duplication
105 edgeListList faceIntoEdges(locFcs.size());
107 forAll (locFcs, faceI)
109 faceIntoEdges[faceI] = locFcs[faceI].edges();
111 labelList& curFaceEdges = faceEdges[faceI];
112 curFaceEdges.setSize(faceIntoEdges[faceI].size());
114 forAll (curFaceEdges, faceEdgeI)
116 curFaceEdges[faceEdgeI] = -1;
120 // This algorithm will produce a separated list of edges, internal edges
121 // starting from 0 and boundary edges starting from the top and
128 // Note that faceIntoEdges is sorted acc. to local vertex numbering
129 // in face (i.e. curEdges[0] is edge between f[0] and f[1])
131 // For all local faces ...
132 forAll (locFcs, faceI)
134 // Get reference to vertices of current face and corresponding edges.
135 const Face& curF = locFcs[faceI];
136 const edgeList& curEdges = faceIntoEdges[faceI];
138 // Record the neighbour face. Multiple connectivity allowed
139 List<DynamicList<label> > neiFaces(curF.size());
140 List<DynamicList<label> > edgeOfNeiFace(curF.size());
142 label nNeighbours = 0;
145 forAll (curEdges, edgeI)
147 // If the edge is already detected, skip
148 if (faceEdges[faceI][edgeI] >= 0) continue;
152 // Set reference to the current edge
153 const edge& e = curEdges[edgeI];
155 // Collect neighbours for the current face vertex.
157 const labelList& nbrFaces = pf[e.start()];
159 forAll (nbrFaces, nbrFaceI)
161 // set reference to the current neighbour
162 label curNei = nbrFaces[nbrFaceI];
164 // Reject neighbours with the lower label
167 // get the reference to subshapes of the neighbour
168 const edgeList& searchEdges = faceIntoEdges[curNei];
170 forAll (searchEdges, neiEdgeI)
172 if (searchEdges[neiEdgeI] == e)
177 neiFaces[edgeI].append(curNei);
178 edgeOfNeiFace[edgeI].append(neiEdgeI);
180 // Record faceFaces both ways
181 ff[faceI].append(curNei);
182 ff[curNei].append(faceI);
184 // Keep searching due to multiple connectivity
188 } // End of neighbouring faces
192 // Register another detected internal edge
195 } // End of current edges
197 // Add the edges in increasing number of neighbours.
198 // Note: for multiply connected surfaces, the lower index neighbour for
199 // an edge will come first.
201 // Add the faces in the increasing order of neighbours
202 for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
204 // Find the lowest neighbour which is still valid
206 label minNei = locFcs.size();
208 forAll (neiFaces, nfI)
210 if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
213 minNei = neiFaces[nfI][0];
219 // Add the face to the list of faces
220 edges[nEdges] = curEdges[nextNei];
222 // Set face-edge and face-neighbour-edge to current face label
223 faceEdges[faceI][nextNei] = nEdges;
225 DynamicList<label>& cnf = neiFaces[nextNei];
226 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
228 // Set edge-face addressing
229 labelList& curEf = edgeFaces[nEdges];
230 curEf.setSize(cnf.size() + 1);
235 faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
237 curEf[cnfI + 1] = cnf[cnfI];
240 // Stop the neighbour from being used again
244 // Increment number of faces counter
251 "PrimitivePatch<Face, FaceList, PointField, PointType>::"
253 ) << "Error in internal edge insertion"
254 << abort(FatalError);
259 nInternalEdges_ = nEdges;
263 forAll (faceEdges, faceI)
265 labelList& curEdges = faceEdges[faceI];
267 forAll (curEdges, edgeI)
269 if (curEdges[edgeI] < 0)
271 // Grab edge and faceEdge
272 edges[nEdges] = faceIntoEdges[faceI][edgeI];
273 curEdges[edgeI] = nEdges;
276 labelList& curEf = edgeFaces[nEdges];
286 edges.setSize(nEdges);
289 edgeFaces.setSize(nEdges);
292 faceFacesPtr_ = new labelListList(locFcs.size());
293 labelListList& faceFaces = *faceFacesPtr_;
295 forAll (faceFaces, faceI)
297 faceFaces[faceI].transfer(ff[faceI]);
303 Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
304 << "calcAddressing() : finished calculating patch addressing"
310 // ************************************************************************* //