initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / OpenFOAM / meshes / primitiveMesh / PrimitivePatch / PrimitivePatchAddressing.C
blob825af4cb706757fe69614ceeb264cd00208e164d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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
25 Description
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
49 template
51     class Face,
52     template<class> class FaceList,
53     class PointField,
54     class PointType
56 void PrimitivePatch<Face, FaceList, PointField, PointType>::calcAddressing()
57  const
59     if (debug)
60     {
61         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
62             << "calcAddressing() : calculating patch addressing"
63             << endl;
64     }
66     if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
67     {
68         // it is considered an error to attempt to recalculate
69         // if already allocated
70         FatalErrorIn
71         (
72             "PrimitivePatch<Face, FaceList, PointField, PointType>::"
73             "calcAddressing()"
74         )   << "addressing already calculated"
75             << abort(FatalError);
76     }
78     // get reference to localFaces
79     const List<Face>& locFcs = localFaces();
81     // get reference to pointFaces
82     const labelListList& pf = pointFaces();
84     // Guess the max number of edges and neighbours for a face
85     label maxEdges = 0;
86     forAll (locFcs, faceI)
87     {
88         maxEdges += locFcs[faceI].size();
89     }
91     // create the lists for the various results. (resized on completion)
92     edgesPtr_ = new edgeList(maxEdges);
93     edgeList& edges = *edgesPtr_;
95     edgeFacesPtr_ = new labelListList(maxEdges);
96     labelListList& edgeFaces = *edgeFacesPtr_;
98     // faceFaces created using a dynamic list.  Cannot guess size because
99     // of multiple connections
100     List<DynamicList<label> > ff(locFcs.size());
102     faceEdgesPtr_ = new labelListList(locFcs.size());
103     labelListList& faceEdges = *faceEdgesPtr_;
105     // count the number of face neighbours
106     labelList noFaceFaces(locFcs.size());
108     // initialise the lists of subshapes for each face to avoid duplication
109     edgeListList faceIntoEdges(locFcs.size());
111     forAll (locFcs, faceI)
112     {
113         faceIntoEdges[faceI] = locFcs[faceI].edges();
115         labelList& curFaceEdges = faceEdges[faceI];
116         curFaceEdges.setSize(faceIntoEdges[faceI].size());
118         forAll (curFaceEdges, faceEdgeI)
119         {
120             curFaceEdges[faceEdgeI] = -1;
121         }
122     }
124     // This algorithm will produce a separated list of edges, internal edges
125     // starting from 0 and boundary edges starting from the top and
126     // growing down.
128     label nEdges = 0;
130     bool found = false;
132     // Note that faceIntoEdges is sorted acc. to local vertex numbering
133     // in face (i.e. curEdges[0] is edge between f[0] and f[1])
135     // For all local faces ...
136     forAll (locFcs, faceI)
137     {
138         // Get reference to vertices of current face and corresponding edges.
139         const Face& curF = locFcs[faceI];
140         const edgeList& curEdges = faceIntoEdges[faceI];
142         // Record the neighbour face.  Multiple connectivity allowed
143         List<DynamicList<label> > neiFaces(curF.size());
144         List<DynamicList<label> > edgeOfNeiFace(curF.size());
146         label nNeighbours = 0;
148         // For all edges ...
149         forAll (curEdges, edgeI)
150         {
151             // If the edge is already detected, skip
152             if (faceEdges[faceI][edgeI] >= 0) continue;
154             found = false;
156             // Set reference to the current edge
157             const edge& e = curEdges[edgeI];
159             // Collect neighbours for the current face vertex.
161             const labelList& nbrFaces = pf[e.start()];
163             forAll (nbrFaces, nbrFaceI)
164             {
165                 // set reference to the current neighbour
166                 label curNei = nbrFaces[nbrFaceI];
168                 // Reject neighbours with the lower label
169                 if (curNei > faceI)
170                 {
171                     // get the reference to subshapes of the neighbour
172                     const edgeList& searchEdges = faceIntoEdges[curNei];
174                     forAll (searchEdges, neiEdgeI)
175                     {
176                         if (searchEdges[neiEdgeI] == e)
177                         {
178                             // Match
179                             found = true;
181                             neiFaces[edgeI].append(curNei);
182                             edgeOfNeiFace[edgeI].append(neiEdgeI);
184                             // Record faceFaces both ways
185                             ff[faceI].append(curNei);
186                             ff[curNei].append(faceI);
188                             // Keep searching due to multiple connectivity
189                         }
190                     }
191                 }
192             } // End of neighbouring faces
194             if (found)
195             {
196                 // Register another detected internal edge
197                 nNeighbours++;
198             }
199         } // End of current edges
201         // Add the edges in increasing number of neighbours.
202         // Note: for multiply connected surfaces, the lower index neighbour for
203         // an edge will come first.
205         // Add the faces in the increasing order of neighbours
206         for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
207         {
208             // Find the lowest neighbour which is still valid
209             label nextNei = -1;
210             label minNei = locFcs.size();
212             forAll (neiFaces, nfI)
213             {
214                 if (neiFaces[nfI].size() > 0 && neiFaces[nfI][0] < minNei)
215                 {
216                     nextNei = nfI;
217                     minNei = neiFaces[nfI][0];
218                 }
219             }
221             if (nextNei > -1)
222             {
223                 // Add the face to the list of faces
224                 edges[nEdges] = curEdges[nextNei];
226                 // Set face-edge and face-neighbour-edge to current face label
227                 faceEdges[faceI][nextNei] = nEdges;
229                 DynamicList<label>& cnf = neiFaces[nextNei];
230                 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
232                 // Set edge-face addressing
233                 labelList& curEf = edgeFaces[nEdges];
234                 curEf.setSize(cnf.size() + 1);
235                 curEf[0] = faceI;
237                 forAll (cnf, cnfI)
238                 {
239                     faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
241                     curEf[cnfI + 1] = cnf[cnfI];
242                 }
244                 // Stop the neighbour from being used again
245                 cnf.clear();
246                 eonf.clear();
248                 // Increment number of faces counter
249                 nEdges++;
250             }
251             else
252             {
253                 FatalErrorIn
254                 (
255                     "PrimitivePatch<Face, FaceList, PointField, PointType>::"
256                     "calcAddressing()"
257                 )   << "Error in internal edge insertion"
258                     << abort(FatalError);
259             }
260         }
261     }
263     nInternalEdges_ = nEdges;
265     // Do boundary faces
267     forAll (faceEdges, faceI)
268     {
269         labelList& curEdges = faceEdges[faceI];
271         forAll (curEdges, edgeI)
272         {
273             if (curEdges[edgeI] < 0)
274             {
275                 // Grab edge and faceEdge
276                 edges[nEdges] = faceIntoEdges[faceI][edgeI];
277                 curEdges[edgeI] = nEdges;
279                 // Add edgeFace
280                 labelList& curEf = edgeFaces[nEdges];
281                 curEf.setSize(1);
282                 curEf[0] = faceI;
284                 nEdges++;
285             }
286         }
287     }
289     // edges
290     edges.setSize(nEdges);
292     // edgeFaces list
293     edgeFaces.setSize(nEdges);
295     // faceFaces list
296     faceFacesPtr_ = new labelListList(locFcs.size());
297     labelListList& faceFaces = *faceFacesPtr_;
299     forAll (faceFaces, faceI)
300     {
301         faceFaces[faceI].transfer(ff[faceI].shrink());
302     }
305     if (debug)
306     {
307         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
308             << "calcAddressing() : finished calculating patch addressing"
309             << endl;
310     }
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 } // End namespace Foam
318 // ************************************************************************* //