1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "PatchTools.H"
27 #include "SortableList.H"
28 #include "transform.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 template<class> class FaceList,
41 Foam::PatchTools::sortedEdgeFaces
43 const PrimitivePatch<Face, FaceList, PointField, PointType>& p
46 const edgeList& edges = p.edges();
47 const labelListList& edgeFaces = p.edgeFaces();
48 const List<Face>& localFaces = p.localFaces();
49 const Field<PointType>& localPoints = p.localPoints();
51 // create the lists for the various results. (resized on completion)
52 labelListList sortedEdgeFaces(edgeFaces.size());
54 forAll (edgeFaces, edgeI)
56 const labelList& faceNbs = edgeFaces[edgeI];
58 if (faceNbs.size() > 2)
60 // Get point on edge and normalized direction of edge (= e2 base
61 // of our coordinate system)
62 const edge& e = edges[edgeI];
64 const point& edgePt = localPoints[e.start()];
66 vector e2 = e.vec(localPoints);
67 e2 /= mag(e2) + VSMALL;
69 // Get opposite vertex for 0th face
70 const Face& f = localFaces[faceNbs[0]];
72 label fp0 = findIndex(f, e[0]);
73 label fp1 = f.fcIndex(fp0);
74 label vertI = (f[fp1] != e[1] ? f[fp1] : f.fcIndex(fp1));
76 // Get vector normal both to e2 and to edge from opposite vertex
77 // to edge (will be x-axis of our coordinate system)
78 vector e0 = e2 ^ (localPoints[vertI] - edgePt);
79 e0 /= mag(e0) + VSMALL;
81 // Get y-axis of coordinate system
84 SortableList<scalar> faceAngles(faceNbs.size());
86 // e0 is reference so angle is 0
89 for (label nbI = 1; nbI < faceNbs.size(); nbI++)
91 // Get opposite vertex
92 const Face& f = localFaces[faceNbs[nbI]];
93 label fp0 = findIndex(f, e[0]);
94 label fp1 = f.fcIndex(fp0);
95 label vertI = (f[fp1] != e[1] ? f[fp1] : f.fcIndex(fp1));
97 vector vec = e2 ^ (localPoints[vertI] - edgePt);
98 vec /= mag(vec) + VSMALL;
100 faceAngles[nbI] = pseudoAngle
110 sortedEdgeFaces[edgeI] = UIndirectList<label>
118 // No need to sort. Just copy.
119 sortedEdgeFaces[edgeI] = faceNbs;
123 return sortedEdgeFaces;
127 // ************************************************************************* //