PatchToolsSortEdges.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2013 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PatchTools.H"
29 #include "SortableList.H"
30 #include "transform.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template
35 <
36  class Face,
37  template<class> class FaceList,
38  class PointField,
39  class PointType
40 >
41 
44 (
46 )
47 {
48  const edgeList& edges = p.edges();
49  const labelListList& edgeFaces = p.edgeFaces();
50  const List<Face>& localFaces = p.localFaces();
51  const Field<PointType>& localPoints = p.localPoints();
52 
53  // create the lists for the various results. (resized on completion)
54  labelListList sortedEdgeFaces(edgeFaces.size());
55 
56  forAll(edgeFaces, edgeI)
57  {
58  const labelList& faceNbs = edgeFaces[edgeI];
59 
60  if (faceNbs.size() > 2)
61  {
62  // Get point on edge and normalized direction of edge (= e2 base
63  // of our coordinate system)
64  const edge& e = edges[edgeI];
65 
66  const point& edgePt = localPoints[e.start()];
67 
68  const vector e2 = e.unitVec(localPoints);
69 
70  // Get the vertex on 0th face that forms a vector with the first
71  // edge point that has the largest angle with the edge
72  const Face& f0 = localFaces[faceNbs[0]];
73 
74  scalar maxAngle = GREAT;
75  vector maxAngleEdgeDir(vector::max);
76 
77  forAll(f0, fpI)
78  {
79  if (f0[fpI] != e.start())
80  {
81  const vector faceEdgeDir =
83  (
84  localPoints[f0[fpI]] - edgePt
85  );
86 
87  const scalar angle = e2 & faceEdgeDir;
88 
89  if (mag(angle) < maxAngle)
90  {
91  maxAngle = angle;
92  maxAngleEdgeDir = faceEdgeDir;
93  }
94  }
95  }
96 
97  // Get vector normal both to e2 and to edge from opposite vertex
98  // to edge (will be x-axis of our coordinate system)
99  const vector e0 = normalised(e2 ^ maxAngleEdgeDir);
100 
101  // Get y-axis of coordinate system
102  const vector e1 = e2 ^ e0;
103 
104  SortableList<scalar> faceAngles(faceNbs.size());
105 
106  // e0 is reference so angle is 0
107  faceAngles[0] = 0;
108 
109  for (label nbI = 1; nbI < faceNbs.size(); nbI++)
110  {
111  // Get the vertex on face that forms a vector with the first
112  // edge point that has the largest angle with the edge
113  const Face& f = localFaces[faceNbs[nbI]];
114 
115  maxAngle = GREAT;
116  maxAngleEdgeDir = vector::max;
117 
118  forAll(f, fpI)
119  {
120  if (f[fpI] != e.start())
121  {
122  const vector faceEdgeDir =
123  normalised
124  (
125  localPoints[f[fpI]] - edgePt
126  );
127 
128  const scalar angle = e2 & faceEdgeDir;
129 
130  if (mag(angle) < maxAngle)
131  {
132  maxAngle = angle;
133  maxAngleEdgeDir = faceEdgeDir;
134  }
135  }
136  }
137 
138  const vector vec = normalised(e2 ^ maxAngleEdgeDir);
139 
140  faceAngles[nbI] = pseudoAngle
141  (
142  e0,
143  e1,
144  vec
145  );
146  }
147 
148  faceAngles.sort();
149 
150  sortedEdgeFaces[edgeI] = labelUIndList
151  (
152  faceNbs,
153  faceAngles.indices()
154  );
155  }
156  else
157  {
158  // No need to sort. Just copy.
159  sortedEdgeFaces[edgeI] = faceNbs;
160  }
161  }
162 
163  return sortedEdgeFaces;
164 }
165 
166 
167 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::SortableList::sort
void sort()
Forward (stable) sort the list (if changed after construction).
Definition: SortableList.C:138
PatchTools.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:62
Foam::Field< PointType >
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:66
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:424
Foam::PatchTools::sortedEdgeFaces
static labelListList sortedEdgeFaces(const PrimitivePatch< Face, FaceList, PointField, PointType > &)
Return edge-face addressing sorted by angle around the edge.
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:496
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:102
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:113
transform.H
3D tensor transformation operations.
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:59