FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

PrimitivePatchAddressing.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Description
00025     This function calculates the list of patch edges, defined on the list of
00026     points supporting the patch. The edges are ordered:
00027     - 0..nInternalEdges-1 : upper triangular order
00028     - nInternalEdges..    : boundary edges (no particular order)
00029 
00030     Other patch addressing information is also calculated:
00031     - faceFaces with neighbour faces in ascending order
00032     - edgeFaces with ascending face order
00033     - faceEdges sorted according to edges of a face
00034 
00035 \*---------------------------------------------------------------------------*/
00036 
00037 #include "PrimitivePatch_.H"
00038 #include <OpenFOAM/DynamicList.H>
00039 
00040 
00041 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00042 
00043 template
00044 <
00045     class Face,
00046     template<class> class FaceList,
00047     class PointField,
00048     class PointType
00049 >
00050 void
00051 Foam::PrimitivePatch<Face, FaceList, PointField, PointType>::
00052 calcAddressing() const
00053 {
00054     if (debug)
00055     {
00056         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
00057             << "calcAddressing() : calculating patch addressing"
00058             << endl;
00059     }
00060 
00061     if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
00062     {
00063         // it is considered an error to attempt to recalculate
00064         // if already allocated
00065         FatalErrorIn
00066         (
00067             "PrimitivePatch<Face, FaceList, PointField, PointType>::"
00068             "calcAddressing()"
00069         )   << "addressing already calculated"
00070             << abort(FatalError);
00071     }
00072 
00073     // get reference to localFaces
00074     const List<Face>& locFcs = localFaces();
00075 
00076     // get reference to pointFaces
00077     const labelListList& pf = pointFaces();
00078 
00079     // Guess the max number of edges and neighbours for a face
00080     label maxEdges = 0;
00081     forAll (locFcs, faceI)
00082     {
00083         maxEdges += locFcs[faceI].size();
00084     }
00085 
00086     // create the lists for the various results. (resized on completion)
00087     edgesPtr_ = new edgeList(maxEdges);
00088     edgeList& edges = *edgesPtr_;
00089 
00090     edgeFacesPtr_ = new labelListList(maxEdges);
00091     labelListList& edgeFaces = *edgeFacesPtr_;
00092 
00093     // faceFaces created using a dynamic list.  Cannot guess size because
00094     // of multiple connections
00095     List<DynamicList<label> > ff(locFcs.size());
00096 
00097     faceEdgesPtr_ = new labelListList(locFcs.size());
00098     labelListList& faceEdges = *faceEdgesPtr_;
00099 
00100     // count the number of face neighbours
00101     labelList noFaceFaces(locFcs.size());
00102 
00103     // initialise the lists of subshapes for each face to avoid duplication
00104     edgeListList faceIntoEdges(locFcs.size());
00105 
00106     forAll (locFcs, faceI)
00107     {
00108         faceIntoEdges[faceI] = locFcs[faceI].edges();
00109 
00110         labelList& curFaceEdges = faceEdges[faceI];
00111         curFaceEdges.setSize(faceIntoEdges[faceI].size());
00112 
00113         forAll (curFaceEdges, faceEdgeI)
00114         {
00115             curFaceEdges[faceEdgeI] = -1;
00116         }
00117     }
00118 
00119     // This algorithm will produce a separated list of edges, internal edges
00120     // starting from 0 and boundary edges starting from the top and
00121     // growing down.
00122 
00123     label nEdges = 0;
00124 
00125     bool found = false;
00126 
00127     // Note that faceIntoEdges is sorted acc. to local vertex numbering
00128     // in face (i.e. curEdges[0] is edge between f[0] and f[1])
00129 
00130     // For all local faces ...
00131     forAll (locFcs, faceI)
00132     {
00133         // Get reference to vertices of current face and corresponding edges.
00134         const Face& curF = locFcs[faceI];
00135         const edgeList& curEdges = faceIntoEdges[faceI];
00136 
00137         // Record the neighbour face.  Multiple connectivity allowed
00138         List<DynamicList<label> > neiFaces(curF.size());
00139         List<DynamicList<label> > edgeOfNeiFace(curF.size());
00140 
00141         label nNeighbours = 0;
00142 
00143         // For all edges ...
00144         forAll (curEdges, edgeI)
00145         {
00146             // If the edge is already detected, skip
00147             if (faceEdges[faceI][edgeI] >= 0) continue;
00148 
00149             found = false;
00150 
00151             // Set reference to the current edge
00152             const edge& e = curEdges[edgeI];
00153 
00154             // Collect neighbours for the current face vertex.
00155 
00156             const labelList& nbrFaces = pf[e.start()];
00157 
00158             forAll (nbrFaces, nbrFaceI)
00159             {
00160                 // set reference to the current neighbour
00161                 label curNei = nbrFaces[nbrFaceI];
00162 
00163                 // Reject neighbours with the lower label
00164                 if (curNei > faceI)
00165                 {
00166                     // get the reference to subshapes of the neighbour
00167                     const edgeList& searchEdges = faceIntoEdges[curNei];
00168 
00169                     forAll (searchEdges, neiEdgeI)
00170                     {
00171                         if (searchEdges[neiEdgeI] == e)
00172                         {
00173                             // Match
00174                             found = true;
00175 
00176                             neiFaces[edgeI].append(curNei);
00177                             edgeOfNeiFace[edgeI].append(neiEdgeI);
00178 
00179                             // Record faceFaces both ways
00180                             ff[faceI].append(curNei);
00181                             ff[curNei].append(faceI);
00182 
00183                             // Keep searching due to multiple connectivity
00184                         }
00185                     }
00186                 }
00187             } // End of neighbouring faces
00188 
00189             if (found)
00190             {
00191                 // Register another detected internal edge
00192                 nNeighbours++;
00193             }
00194         } // End of current edges
00195 
00196         // Add the edges in increasing number of neighbours.
00197         // Note: for multiply connected surfaces, the lower index neighbour for
00198         // an edge will come first.
00199 
00200         // Add the faces in the increasing order of neighbours
00201         for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
00202         {
00203             // Find the lowest neighbour which is still valid
00204             label nextNei = -1;
00205             label minNei = locFcs.size();
00206 
00207             forAll (neiFaces, nfI)
00208             {
00209                 if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
00210                 {
00211                     nextNei = nfI;
00212                     minNei = neiFaces[nfI][0];
00213                 }
00214             }
00215 
00216             if (nextNei > -1)
00217             {
00218                 // Add the face to the list of faces
00219                 edges[nEdges] = curEdges[nextNei];
00220 
00221                 // Set face-edge and face-neighbour-edge to current face label
00222                 faceEdges[faceI][nextNei] = nEdges;
00223 
00224                 DynamicList<label>& cnf = neiFaces[nextNei];
00225                 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
00226 
00227                 // Set edge-face addressing
00228                 labelList& curEf = edgeFaces[nEdges];
00229                 curEf.setSize(cnf.size() + 1);
00230                 curEf[0] = faceI;
00231 
00232                 forAll (cnf, cnfI)
00233                 {
00234                     faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
00235 
00236                     curEf[cnfI + 1] = cnf[cnfI];
00237                 }
00238 
00239                 // Stop the neighbour from being used again
00240                 cnf.clear();
00241                 eonf.clear();
00242 
00243                 // Increment number of faces counter
00244                 nEdges++;
00245             }
00246             else
00247             {
00248                 FatalErrorIn
00249                 (
00250                     "PrimitivePatch<Face, FaceList, PointField, PointType>::"
00251                     "calcAddressing()"
00252                 )   << "Error in internal edge insertion"
00253                     << abort(FatalError);
00254             }
00255         }
00256     }
00257 
00258     nInternalEdges_ = nEdges;
00259 
00260     // Do boundary faces
00261 
00262     forAll (faceEdges, faceI)
00263     {
00264         labelList& curEdges = faceEdges[faceI];
00265 
00266         forAll (curEdges, edgeI)
00267         {
00268             if (curEdges[edgeI] < 0)
00269             {
00270                 // Grab edge and faceEdge
00271                 edges[nEdges] = faceIntoEdges[faceI][edgeI];
00272                 curEdges[edgeI] = nEdges;
00273 
00274                 // Add edgeFace
00275                 labelList& curEf = edgeFaces[nEdges];
00276                 curEf.setSize(1);
00277                 curEf[0] = faceI;
00278 
00279                 nEdges++;
00280             }
00281         }
00282     }
00283 
00284     // edges
00285     edges.setSize(nEdges);
00286 
00287     // edgeFaces list
00288     edgeFaces.setSize(nEdges);
00289 
00290     // faceFaces list
00291     faceFacesPtr_ = new labelListList(locFcs.size());
00292     labelListList& faceFaces = *faceFacesPtr_;
00293 
00294     forAll (faceFaces, faceI)
00295     {
00296         faceFaces[faceI].transfer(ff[faceI]);
00297     }
00298 
00299 
00300     if (debug)
00301     {
00302         Info<< "PrimitivePatch<Face, FaceList, PointField, PointType>::"
00303             << "calcAddressing() : finished calculating patch addressing"
00304             << endl;
00305     }
00306 }
00307 
00308 
00309 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines