PrimitivePatchAddressing.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-2021 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
27Description
28 This function calculates the list of patch edges, defined on the list of
29 points supporting the patch. The edges are ordered:
30 - 0..nInternalEdges-1 : upper triangular order
31 - nInternalEdges.. : boundary edges (no particular order)
32
33 Other patch addressing information is also calculated:
34 - faceFaces with neighbour faces in ascending order
35 - edgeFaces with ascending face order
36 - faceEdges sorted according to edges of a face
37
38\*---------------------------------------------------------------------------*/
39
40#include "PrimitivePatch.H"
41#include "DynamicList.H"
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45template<class FaceList, class PointField>
46void
48{
49 DebugInFunction << "Calculating patch addressing" << nl;
50
51 if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
52 {
53 // An error to recalculate if already allocated
55 << "addressing already calculated"
56 << abort(FatalError);
57 }
58
59 // get reference to localFaces
60 const List<face_type>& locFcs = localFaces();
61
62 // get reference to pointFaces
63 const labelListList& pf = pointFaces();
64
65 // Guess the max number of edges/neighbours for a face
66 label maxEdges = 0;
67 for (const auto& f : locFcs)
68 {
69 maxEdges += f.nEdges();
70 }
71
72 // create the lists for the various results. (resized on completion)
73 edgesPtr_.reset(new edgeList(maxEdges));
74 auto& edges = *edgesPtr_;
75
76 edgeFacesPtr_.reset(new labelListList(maxEdges));
77 auto& edgeFaces = *edgeFacesPtr_;
78
79 faceEdgesPtr_.reset(new labelListList(locFcs.size()));
80 auto& faceEdges = *faceEdgesPtr_;
81
82 // initialise the lists of subshapes for each face to avoid duplication
83 edgeListList faceIntoEdges(locFcs.size());
84
85 forAll(locFcs, facei)
86 {
87 faceIntoEdges[facei] = locFcs[facei].edges();
88 faceEdges[facei].resize(faceIntoEdges[facei].size(), -1);
89 }
90
91 // faceFaces created using a dynamic list. Cannot guess size because
92 // of multiple connections
93 List<DynamicList<label>> ff(locFcs.size());
94
95
96 // This algorithm will produce a separated list of edges, internal edges
97 // starting from 0 and boundary edges starting from the top and
98 // growing down.
99
100 label nEdges = 0;
101
102 // Note that faceIntoEdges is sorted acc. to local vertex numbering
103 // in face (i.e. curEdges[0] is edge between f[0] and f[1])
104
105 // For all local faces ...
106 forAll(locFcs, facei)
107 {
108 // Get reference to vertices of current face and corresponding edges.
109 const face_type& curF = locFcs[facei];
110 const edgeList& curEdges = faceIntoEdges[facei];
111
112 // Record the neighbour face. Multiple connectivity allowed
113 List<DynamicList<label>> neiFaces(curF.size());
114 List<DynamicList<label>> edgeOfNeiFace(curF.size());
115
116 label nNeighbours = 0;
117
118 // For all edges ...
119 forAll(curEdges, edgeI)
120 {
121 // If the edge is already detected, skip
122 if (faceEdges[facei][edgeI] >= 0) continue;
123
124 // Set reference to the current edge
125 const edge& e = curEdges[edgeI];
126
127 // Collect neighbours for the current face vertex.
128
129 const labelList& nbrFaces = pf[e.start()];
130
131 bool found = false;
132
133 forAll(nbrFaces, nbrFacei)
134 {
135 // set reference to the current neighbour
136 label curNei = nbrFaces[nbrFacei];
137
138 // Reject neighbours with the lower label
139 if (curNei > facei)
140 {
141 // get the reference to subshapes of the neighbour
142 const edgeList& searchEdges = faceIntoEdges[curNei];
143
144 forAll(searchEdges, neiEdgeI)
145 {
146 if (searchEdges[neiEdgeI] == e)
147 {
148 // Match
149 found = true;
150
151 neiFaces[edgeI].append(curNei);
152 edgeOfNeiFace[edgeI].append(neiEdgeI);
153
154 // Record faceFaces both ways
155 ff[facei].append(curNei);
156 ff[curNei].append(facei);
157
158 // Keep searching due to multiple connectivity
159 }
160 }
161 }
162 } // End of neighbouring faces
163
164 if (found)
165 {
166 // Register another detected internal edge
167 nNeighbours++;
168 }
169 } // End of current edges
170
171 // Add the edges in increasing number of neighbours.
172 // Note: for multiply connected surfaces, the lower index neighbour for
173 // an edge will come first.
174
175 // Add the faces in the increasing order of neighbours
176 for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
177 {
178 // Find the lowest neighbour which is still valid
179 label nextNei = -1;
180 label minNei = locFcs.size();
181
182 forAll(neiFaces, nfI)
183 {
184 if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
185 {
186 nextNei = nfI;
187 minNei = neiFaces[nfI][0];
188 }
189 }
190
191 if (nextNei > -1)
192 {
193 // Add the face to the list of faces
194 edges[nEdges] = curEdges[nextNei];
195
196 // Set face-edge and face-neighbour-edge to current face label
197 faceEdges[facei][nextNei] = nEdges;
198
199 DynamicList<label>& cnf = neiFaces[nextNei];
200 DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
201
202 // Set edge-face addressing
203 labelList& curEf = edgeFaces[nEdges];
204 curEf.resize(cnf.size() + 1);
205 curEf[0] = facei;
206
207 forAll(cnf, cnfI)
208 {
209 faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
210
211 curEf[cnfI + 1] = cnf[cnfI];
212 }
213
214 // Stop the neighbour from being used again
215 cnf.clear();
216 eonf.clear();
217
218 // Increment number of faces counter
219 nEdges++;
220 }
221 else
222 {
224 << "Error in internal edge insertion"
225 << abort(FatalError);
226 }
227 }
228 }
229
230 nInternalEdges_ = nEdges;
231
232 // Do boundary faces
233
234 forAll(faceEdges, facei)
235 {
236 labelList& curEdges = faceEdges[facei];
237
238 forAll(curEdges, edgeI)
239 {
240 if (curEdges[edgeI] < 0)
241 {
242 // Grab edge and faceEdge
243 edges[nEdges] = faceIntoEdges[facei][edgeI];
244 curEdges[edgeI] = nEdges;
245
246 // Add edgeFace
247 labelList& curEf = edgeFaces[nEdges];
248 curEf.resize(1);
249 curEf[0] = facei;
250
251 nEdges++;
252 }
253 }
254 }
255
256 // edges
257 edges.resize(nEdges);
258
259 // edgeFaces list
260 edgeFaces.resize(nEdges);
261
262 // faceFaces list
263 faceFacesPtr_.reset(new labelListList(locFcs.size()));
264 auto& faceFaces = *faceFacesPtr_;
265
266 forAll(faceFaces, facei)
267 {
268 faceFaces[facei].transfer(ff[facei]);
269 }
270
271
272 DebugInFunction << "Calculated patch addressing" << nl;
273}
274
275
276// ************************************************************************* //
bool found
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A list of faces which address into the list of points.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
#define DebugInFunction
Report an information message using Foam::Info.
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
List< label > labelList
A List of labels.
Definition: List.H:66
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
List< edgeList > edgeListList
A List of edgeList.
Definition: edgeList.H:65
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333