PatchToolsSortEdges.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
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 -------------------------------------------------------------------------------
8 License
9  This file is part of foam-extend.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PatchTools.H"
27 #include "SortableList.H"
28 #include "transform.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template
33 <
34  class Face,
35  template<class> class FaceList,
36  class PointField,
37  class PointType
38 >
39 
42 (
44 )
45 {
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();
50 
51  // create the lists for the various results. (resized on completion)
52  labelListList sortedEdgeFaces(edgeFaces.size());
53 
54  forAll (edgeFaces, edgeI)
55  {
56  const labelList& faceNbs = edgeFaces[edgeI];
57 
58  if (faceNbs.size() > 2)
59  {
60  // Get point on edge and normalized direction of edge (= e2 base
61  // of our coordinate system)
62  const edge& e = edges[edgeI];
63 
64  const point& edgePt = localPoints[e.start()];
65 
66  vector e2 = e.vec(localPoints);
67  e2 /= mag(e2) + VSMALL;
68 
69  // Get opposite vertex for 0th face
70  const Face& f = localFaces[faceNbs[0]];
71 
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));
75 
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;
80 
81  // Get y-axis of coordinate system
82  vector e1 = e2 ^ e0;
83 
84  SortableList<scalar> faceAngles(faceNbs.size());
85 
86  // e0 is reference so angle is 0
87  faceAngles[0] = 0;
88 
89  for (label nbI = 1; nbI < faceNbs.size(); nbI++)
90  {
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));
96 
97  vector vec = e2 ^ (localPoints[vertI] - edgePt);
98  vec /= mag(vec) + VSMALL;
99 
100  faceAngles[nbI] = pseudoAngle
101  (
102  e0,
103  e1,
104  vec
105  );
106  }
107 
108  faceAngles.sort();
109 
110  sortedEdgeFaces[edgeI] = UIndirectList<label>
111  (
112  faceNbs,
113  faceAngles.indices()
114  );
115  }
116  else
117  {
118  // No need to sort. Just copy.
119  sortedEdgeFaces[edgeI] = faceNbs;
120  }
121  }
122 
123  return sortedEdgeFaces;
124 }
125 
126 
127 // ************************************************************************* //
p
p
Definition: pEqn.H:62
Foam::SortableList::sort
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:95
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
SortableList.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
PatchTools.H
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:65
Foam::pseudoAngle
scalar pseudoAngle(const vector &e0, const vector &e1, const vector &vec)
Estimate angle of vec in coordinate system (e0, e1, e0^e1).
Definition: transform.H:223
Foam::PatchTools::sortedEdgeFaces
static labelListList sortedEdgeFaces(const PrimitivePatch< Face, FaceList, PointField, PointType > &)
Return edge-face addressing sorted by angle around the edge.
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:90
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
transform.H
3D tensor transformation operations.
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88