fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / OpenFOAM / meshes / primitiveMesh / PrimitivePatch / PrimitivePatchAddressing.C
blob72a1606979c7f78b948b066ed5534e6479ff6805
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
44 template
46     class Face,
47     template<class> class FaceList,
48     class PointField,
49     class PointType
51 void
52 Foam::PrimitivePatch<Face, FaceList, PointField, PointType>::
53 calcAddressing() const
55     if (debug)
56     {
57         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
58             << "calcAddressing() : calculating patch addressing"
59             << endl;
60     }
62     if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
63     {
64         // it is considered an error to attempt to recalculate
65         // if already allocated
66         FatalErrorIn
67         (
68             "PrimitivePatch<Face, FaceList, PointField, PointType>::"
69             "calcAddressing()"
70         )   << "addressing already calculated"
71             << abort(FatalError);
72     }
74     if (this->size() == 0)
75     {
76         // No faces in patch.  All lists are empty
78         edgesPtr_ = new edgeList(0);
79         faceFacesPtr_ = new labelListList(0);
80         edgeFacesPtr_ = new labelListList(0);
81         faceEdgesPtr_ = new labelListList(0);
82         nInternalEdges_ = 0;
84         return;
85     }
87     // Get reference to localFaces
88     const List<Face>& locFcs = localFaces();
90     // Get reference to pointFaces
91     const labelListList& pf = pointFaces();
93     // Guess the max number of edges and neighbours for a face
94     label maxEdges = 0;
95     forAll (locFcs, faceI)
96     {
97         maxEdges += locFcs[faceI].size();
98     }
100     // Create the lists for the various results. (resized on completion)
101     edgesPtr_ = new edgeList(maxEdges);
102     edgeList& edges = *edgesPtr_;
104     edgeFacesPtr_ = new labelListList(maxEdges);
105     labelListList& edgeFaces = *edgeFacesPtr_;
107     // faceFaces created using a dynamic list.  Cannot guess size because
108     // of multiple connections
109     List<DynamicList<label> > ff(locFcs.size());
111     faceEdgesPtr_ = new labelListList(locFcs.size());
112     labelListList& faceEdges = *faceEdgesPtr_;
114     // Count the number of face neighbours
115     labelList noFaceFaces(locFcs.size());
117     // initialise the lists of subshapes for each face to avoid duplication
118     edgeListList faceIntoEdges(locFcs.size());
120     forAll (locFcs, faceI)
121     {
122         faceIntoEdges[faceI] = locFcs[faceI].edges();
124         labelList& curFaceEdges = faceEdges[faceI];
125         curFaceEdges.setSize(faceIntoEdges[faceI].size());
127         forAll (curFaceEdges, faceEdgeI)
128         {
129             curFaceEdges[faceEdgeI] = -1;
130         }
131     }
133     // This algorithm will produce a separated list of edges, internal edges
134     // starting from 0 and boundary edges starting from the top and
135     // growing down.
137     label nEdges = 0;
139     bool found = false;
141     // Note that faceIntoEdges is sorted acc. to local vertex numbering
142     // in face (i.e. curEdges[0] is edge between f[0] and f[1])
144     // For all local faces ...
145     forAll (locFcs, faceI)
146     {
147         // Get reference to vertices of current face and corresponding edges.
148         const Face& curF = locFcs[faceI];
149         const edgeList& curEdges = faceIntoEdges[faceI];
151         // Record the neighbour face.  Multiple connectivity allowed
152         List<DynamicList<label> > neiFaces(curF.size());
153         List<DynamicList<label> > edgeOfNeiFace(curF.size());
155         label nNeighbours = 0;
157         // For all edges ...
158         forAll (curEdges, edgeI)
159         {
160             // If the edge is already detected, skip
161             if (faceEdges[faceI][edgeI] >= 0) continue;
163             found = false;
165             // Set reference to the current edge
166             const edge& e = curEdges[edgeI];
168             // Collect neighbours for the current face vertex.
170             const labelList& nbrFaces = pf[e.start()];
172             forAll (nbrFaces, nbrFaceI)
173             {
174                 // set reference to the current neighbour
175                 label curNei = nbrFaces[nbrFaceI];
177                 // Reject neighbours with the lower label
178                 if (curNei > faceI)
179                 {
180                     // get the reference to subshapes of the neighbour
181                     const edgeList& searchEdges = faceIntoEdges[curNei];
183                     forAll (searchEdges, neiEdgeI)
184                     {
185                         if (searchEdges[neiEdgeI] == e)
186                         {
187                             // Match
188                             found = true;
190                             neiFaces[edgeI].append(curNei);
191                             edgeOfNeiFace[edgeI].append(neiEdgeI);
193                             // Record faceFaces both ways
194                             ff[faceI].append(curNei);
195                             ff[curNei].append(faceI);
197                             // Keep searching due to multiple connectivity
198                         }
199                     }
200                 }
201             } // End of neighbouring faces
203             if (found)
204             {
205                 // Register another detected internal edge
206                 nNeighbours++;
207             }
208         } // End of current edges
210         // Add the edges in increasing number of neighbours.
211         // Note: for multiply connected surfaces, the lower index neighbour for
212         // an edge will come first.
214         // Add the faces in the increasing order of neighbours
215         for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
216         {
217             // Find the lowest neighbour which is still valid
218             label nextNei = -1;
219             label minNei = locFcs.size();
221             forAll (neiFaces, nfI)
222             {
223                 if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
224                 {
225                     nextNei = nfI;
226                     minNei = neiFaces[nfI][0];
227                 }
228             }
230             if (nextNei > -1)
231             {
232                 // Add the face to the list of faces
233                 edges[nEdges] = curEdges[nextNei];
235                 // Set face-edge and face-neighbour-edge to current face label
236                 faceEdges[faceI][nextNei] = nEdges;
238                 DynamicList<label>& cnf = neiFaces[nextNei];
239                 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
241                 // Set edge-face addressing
242                 labelList& curEf = edgeFaces[nEdges];
243                 curEf.setSize(cnf.size() + 1);
244                 curEf[0] = faceI;
246                 forAll (cnf, cnfI)
247                 {
248                     faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
250                     curEf[cnfI + 1] = cnf[cnfI];
251                 }
253                 // Stop the neighbour from being used again
254                 cnf.clear();
255                 eonf.clear();
257                 // Increment number of faces counter
258                 nEdges++;
259             }
260             else
261             {
262                 FatalErrorIn
263                 (
264                     "PrimitivePatch<Face, FaceList, PointField, PointType>::"
265                     "calcAddressing()"
266                 )   << "Error in internal edge insertion"
267                     << abort(FatalError);
268             }
269         }
270     }
272     nInternalEdges_ = nEdges;
274     // Do boundary faces
276     forAll (faceEdges, faceI)
277     {
278         labelList& curEdges = faceEdges[faceI];
280         forAll (curEdges, edgeI)
281         {
282             if (curEdges[edgeI] < 0)
283             {
284                 // Grab edge and faceEdge
285                 edges[nEdges] = faceIntoEdges[faceI][edgeI];
286                 curEdges[edgeI] = nEdges;
288                 // Add edgeFace
289                 labelList& curEf = edgeFaces[nEdges];
290                 curEf.setSize(1);
291                 curEf[0] = faceI;
293                 nEdges++;
294             }
295         }
296     }
298     // edges
299     edges.setSize(nEdges);
301     // edgeFaces list
302     edgeFaces.setSize(nEdges);
304     // faceFaces list
305     faceFacesPtr_ = new labelListList(locFcs.size());
306     labelListList& faceFaces = *faceFacesPtr_;
308     forAll (faceFaces, faceI)
309     {
310         faceFaces[faceI].transfer(ff[faceI]);
311     }
314     if (debug)
315     {
316         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
317             << "calcAddressing() : finished calculating patch addressing"
318             << endl;
319     }
323 // ************************************************************************* //