PrimitivePatchAddressing.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  This function calculates the list of patch edges, defined on the list of
26  points supporting the patch. The edges are ordered:
27  - 0..nInternalEdges-1 : upper triangular order
28  - nInternalEdges.. : boundary edges (no particular order)
29 
30  Other patch addressing information is also calculated:
31  - faceFaces with neighbour faces in ascending order
32  - edgeFaces with ascending face order
33  - faceEdges sorted according to edges of a face
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "PrimitivePatch.H"
38 #include "DynamicList.H"
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 template<class FaceList, class PointField>
45 {
46  if (debug)
47  {
48  InfoInFunction << "calculating patch addressing" << endl;
49  }
50 
51  if (edgesPtr_ || faceFacesPtr_ || edgeFacesPtr_ || faceEdgesPtr_)
52  {
53  // it is considered an error to attempt to recalculate
54  // if already allocated
56  << "addressing already calculated"
57  << abort(FatalError);
58  }
59 
60  // get reference to localFaces
61  const List<FaceType>& locFcs = localFaces();
62 
63  // get reference to pointFaces
64  const labelListList& pf = pointFaces();
65 
66  // Guess the max number of edges and neighbours for a face
67  label maxEdges = 0;
68  forAll(locFcs, facei)
69  {
70  maxEdges += locFcs[facei].size();
71  }
72 
73  // create the lists for the various results. (resized on completion)
74  edgesPtr_ = new edgeList(maxEdges);
75  edgeList& edges = *edgesPtr_;
76 
77  edgeFacesPtr_ = new labelListList(maxEdges);
78  labelListList& edgeFaces = *edgeFacesPtr_;
79 
80  // faceFaces created using a dynamic list. Cannot guess size because
81  // of multiple connections
82  List<DynamicList<label>> ff(locFcs.size());
83 
84  faceEdgesPtr_ = new labelListList(locFcs.size());
85  labelListList& faceEdges = *faceEdgesPtr_;
86 
87  // count the number of face neighbours
88  labelList noFaceFaces(locFcs.size());
89 
90  // initialise the lists of subshapes for each face to avoid duplication
91  edgeListList faceIntoEdges(locFcs.size());
92 
93  forAll(locFcs, facei)
94  {
95  faceIntoEdges[facei] = locFcs[facei].edges();
96 
97  labelList& curFaceEdges = faceEdges[facei];
98  curFaceEdges.setSize(faceIntoEdges[facei].size());
99 
100  forAll(curFaceEdges, faceEdgeI)
101  {
102  curFaceEdges[faceEdgeI] = -1;
103  }
104  }
105 
106  // This algorithm will produce a separated list of edges, internal edges
107  // starting from 0 and boundary edges starting from the top and
108  // growing down.
109 
110  label nEdges = 0;
111 
112  bool found = false;
113 
114  // Note that faceIntoEdges is sorted acc. to local vertex numbering
115  // in face (i.e. curEdges[0] is edge between f[0] and f[1])
116 
117  // For all local faces ...
118  forAll(locFcs, facei)
119  {
120  // Get reference to vertices of current face and corresponding edges.
121  const FaceType& curF = locFcs[facei];
122  const edgeList& curEdges = faceIntoEdges[facei];
123 
124  // Record the neighbour face. Multiple connectivity allowed
125  List<DynamicList<label>> neiFaces(curF.size());
126  List<DynamicList<label>> edgeOfNeiFace(curF.size());
127 
128  label nNeighbours = 0;
129 
130  // For all edges ...
131  forAll(curEdges, edgeI)
132  {
133  // If the edge is already detected, skip
134  if (faceEdges[facei][edgeI] >= 0) continue;
135 
136  found = false;
137 
138  // Set reference to the current edge
139  const edge& e = curEdges[edgeI];
140 
141  // Collect neighbours for the current face vertex.
142 
143  const labelList& nbrFaces = pf[e.start()];
144 
145  forAll(nbrFaces, nbrFacei)
146  {
147  // set reference to the current neighbour
148  label curNei = nbrFaces[nbrFacei];
149 
150  // Reject neighbours with the lower label
151  if (curNei > facei)
152  {
153  // get the reference to subshapes of the neighbour
154  const edgeList& searchEdges = faceIntoEdges[curNei];
155 
156  forAll(searchEdges, neiEdgeI)
157  {
158  if (searchEdges[neiEdgeI] == e)
159  {
160  // Match
161  found = true;
162 
163  neiFaces[edgeI].append(curNei);
164  edgeOfNeiFace[edgeI].append(neiEdgeI);
165 
166  // Record faceFaces both ways
167  ff[facei].append(curNei);
168  ff[curNei].append(facei);
169 
170  // Keep searching due to multiple connectivity
171  }
172  }
173  }
174  } // End of neighbouring faces
175 
176  if (found)
177  {
178  // Register another detected internal edge
179  nNeighbours++;
180  }
181  } // End of current edges
182 
183  // Add the edges in increasing number of neighbours.
184  // Note: for multiply connected surfaces, the lower index neighbour for
185  // an edge will come first.
186 
187  // Add the faces in the increasing order of neighbours
188  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
189  {
190  // Find the lowest neighbour which is still valid
191  label nextNei = -1;
192  label minNei = locFcs.size();
193 
194  forAll(neiFaces, nfI)
195  {
196  if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
197  {
198  nextNei = nfI;
199  minNei = neiFaces[nfI][0];
200  }
201  }
202 
203  if (nextNei > -1)
204  {
205  // Add the face to the list of faces
206  edges[nEdges] = curEdges[nextNei];
207 
208  // Set face-edge and face-neighbour-edge to current face label
209  faceEdges[facei][nextNei] = nEdges;
210 
211  DynamicList<label>& cnf = neiFaces[nextNei];
212  DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
213 
214  // Set edge-face addressing
215  labelList& curEf = edgeFaces[nEdges];
216  curEf.setSize(cnf.size() + 1);
217  curEf[0] = facei;
218 
219  forAll(cnf, cnfI)
220  {
221  faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
222 
223  curEf[cnfI + 1] = cnf[cnfI];
224  }
225 
226  // Stop the neighbour from being used again
227  cnf.clear();
228  eonf.clear();
229 
230  // Increment number of faces counter
231  nEdges++;
232  }
233  else
234  {
236  << "Error in internal edge insertion"
237  << abort(FatalError);
238  }
239  }
240  }
241 
242  nInternalEdges_ = nEdges;
243 
244  // Do boundary faces
245 
246  forAll(faceEdges, facei)
247  {
248  labelList& curEdges = faceEdges[facei];
249 
250  forAll(curEdges, edgeI)
251  {
252  if (curEdges[edgeI] < 0)
253  {
254  // Grab edge and faceEdge
255  edges[nEdges] = faceIntoEdges[facei][edgeI];
256  curEdges[edgeI] = nEdges;
257 
258  // Add edgeFace
259  labelList& curEf = edgeFaces[nEdges];
260  curEf.setSize(1);
261  curEf[0] = facei;
262 
263  nEdges++;
264  }
265  }
266  }
267 
268  // edges
269  edges.setSize(nEdges);
270 
271  // edgeFaces list
272  edgeFaces.setSize(nEdges);
273 
274  // faceFaces list
275  faceFacesPtr_ = new labelListList(locFcs.size());
276  labelListList& faceFaces = *faceFacesPtr_;
277 
278  forAll(faceFaces, facei)
279  {
280  faceFaces[facei].transfer(ff[facei]);
281  }
282 
283 
284  if (debug)
285  {
286  InfoInFunction << "Finished calculating patch addressing" << endl;
287  }
288 }
289 
290 
291 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< edgeList > edgeListList
Definition: edgeList.H:39
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A list of faces which address into the list of points.
List< edge > edgeList
Definition: edgeList.H:38
List< label > labelList
A List of labels.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of List.
Definition: List.C:281
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
#define InfoInFunction
Report an information message using Foam::Info.