Loading...
Searching...
No Matches
PatchToolsSortPoints.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-2016 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
29#include "PatchTools.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class FaceList, class PointField>
36(
38)
39{
40 // Now order the edges of each point according to whether they share a
41 // face
42 const labelListList& pointEdges = p.pointEdges();
43 const edgeList& edges = p.edges();
44 const labelListList& edgeFaces = p.edgeFaces();
45 const labelListList& faceEdges = p.faceEdges();
46
47 // create the lists for the various results. (resized on completion)
49
50 DynamicList<label> newEdgeList;
51
52 forAll(pointEdges, pointi)
53 {
54 const labelList& pEdges = pointEdges[pointi];
55
56 label nPointEdges = pEdges.size();
57
58 label edgeI = pEdges[0];
59
60 label prevFacei = edgeFaces[edgeI][0];
61
62 newEdgeList.clear();
63 newEdgeList.setCapacity(nPointEdges);
64
65 label nVisitedEdges = 0;
66
67 do
68 {
69 newEdgeList.append(edgeI);
70
71 // Cross edge to next face
72 const labelList& eFaces = edgeFaces[edgeI];
73
74 if (eFaces.size() != 2)
75 {
76 break;
77 }
78
79 label facei = eFaces[0];
80 if (facei == prevFacei)
81 {
82 facei = eFaces[1];
83 }
84
85 // Cross face to next edge
86 const labelList& fEdges = faceEdges[facei];
87
88 forAll(fEdges, feI)
89 {
90 const label nextEdgeI = fEdges[feI];
91 const edge& nextEdge = edges[nextEdgeI];
92
93 if
94 (
95 nextEdgeI != edgeI
96 && (nextEdge.start() == pointi || nextEdge.end() == pointi)
97 )
98 {
99 edgeI = nextEdgeI;
100 break;
101 }
102 }
103
104 prevFacei = facei;
105
106 nVisitedEdges++;
107 if (nVisitedEdges > nPointEdges)
108 {
110 << "Unable to order pointEdges as the face connections "
111 << "are not circular" << nl
112 << " Original pointEdges = " << pEdges << nl
113 << " New pointEdges = " << newEdgeList
114 << endl;
115
116 newEdgeList = pEdges;
117
118 break;
119 }
120
121 } while (edgeI != pEdges[0]);
122
123 if (newEdgeList.size() == nPointEdges)
124 {
125 forAll(pEdges, eI)
126 {
127 if (!newEdgeList.found(pEdges[eI]))
128 {
130 << "Cannot find all original edges in the new list"
131 << nl
132 << " Original pointEdges = " << pEdges << nl
133 << " New pointEdges = " << newEdgeList
134 << endl;
135
136 newEdgeList = pEdges;
137
138 break;
139 }
140 }
141
142 sortedPointEdges[pointi] = newEdgeList;
143 }
144 }
145
146 return sortedPointEdges;
147}
148
149
150// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
void setCapacity(const label len)
Alter the size of the underlying storage.
static labelListList sortedPointEdges(const PrimitivePatch< FaceList, PointField > &)
Return point-edge addressing sorted by order around the point.
A list of faces which address into the list of points.
bool found(const T &val, label pos=0) const
Same as contains().
Definition UList.H:983
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
label & end() noexcept
The end (second/last) vertex label.
Definition edge.H:165
label start() const noexcept
The start (first) vertex label.
Definition edge.H:155
volScalarField & p
#define WarningInFunction
Report a warning using Foam::Warning.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299