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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Description
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 
45 template<class FaceList, class PointField>
46 void
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 and neighbours for a face
66  label maxEdges = 0;
67  for (const auto& f : locFcs)
68  {
69  maxEdges += f.size();
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  // faceFaces created using a dynamic list. Cannot guess size because
80  // of multiple connections
81  List<DynamicList<label>> ff(locFcs.size());
82 
83  faceEdgesPtr_.reset(new labelListList(locFcs.size()));
84  auto& faceEdges = *faceEdgesPtr_;
85 
86  // count the number of face neighbours
87  labelList noFaceFaces(locFcs.size());
88 
89  // initialise the lists of subshapes for each face to avoid duplication
90  edgeListList faceIntoEdges(locFcs.size());
91 
92  forAll(locFcs, facei)
93  {
94  faceIntoEdges[facei] = locFcs[facei].edges();
95 
96  labelList& curFaceEdges = faceEdges[facei];
97  curFaceEdges.setSize(faceIntoEdges[facei].size());
98 
99  forAll(curFaceEdges, faceEdgeI)
100  {
101  curFaceEdges[faceEdgeI] = -1;
102  }
103  }
104 
105  // This algorithm will produce a separated list of edges, internal edges
106  // starting from 0 and boundary edges starting from the top and
107  // growing down.
108 
109  label nEdges = 0;
110 
111  bool found = false;
112 
113  // Note that faceIntoEdges is sorted acc. to local vertex numbering
114  // in face (i.e. curEdges[0] is edge between f[0] and f[1])
115 
116  // For all local faces ...
117  forAll(locFcs, facei)
118  {
119  // Get reference to vertices of current face and corresponding edges.
120  const face_type& curF = locFcs[facei];
121  const edgeList& curEdges = faceIntoEdges[facei];
122 
123  // Record the neighbour face. Multiple connectivity allowed
124  List<DynamicList<label>> neiFaces(curF.size());
125  List<DynamicList<label>> edgeOfNeiFace(curF.size());
126 
127  label nNeighbours = 0;
128 
129  // For all edges ...
130  forAll(curEdges, edgeI)
131  {
132  // If the edge is already detected, skip
133  if (faceEdges[facei][edgeI] >= 0) continue;
134 
135  found = false;
136 
137  // Set reference to the current edge
138  const edge& e = curEdges[edgeI];
139 
140  // Collect neighbours for the current face vertex.
141 
142  const labelList& nbrFaces = pf[e.start()];
143 
144  forAll(nbrFaces, nbrFacei)
145  {
146  // set reference to the current neighbour
147  label curNei = nbrFaces[nbrFacei];
148 
149  // Reject neighbours with the lower label
150  if (curNei > facei)
151  {
152  // get the reference to subshapes of the neighbour
153  const edgeList& searchEdges = faceIntoEdges[curNei];
154 
155  forAll(searchEdges, neiEdgeI)
156  {
157  if (searchEdges[neiEdgeI] == e)
158  {
159  // Match
160  found = true;
161 
162  neiFaces[edgeI].append(curNei);
163  edgeOfNeiFace[edgeI].append(neiEdgeI);
164 
165  // Record faceFaces both ways
166  ff[facei].append(curNei);
167  ff[curNei].append(facei);
168 
169  // Keep searching due to multiple connectivity
170  }
171  }
172  }
173  } // End of neighbouring faces
174 
175  if (found)
176  {
177  // Register another detected internal edge
178  nNeighbours++;
179  }
180  } // End of current edges
181 
182  // Add the edges in increasing number of neighbours.
183  // Note: for multiply connected surfaces, the lower index neighbour for
184  // an edge will come first.
185 
186  // Add the faces in the increasing order of neighbours
187  for (label neiSearch = 0; neiSearch < nNeighbours; neiSearch++)
188  {
189  // Find the lowest neighbour which is still valid
190  label nextNei = -1;
191  label minNei = locFcs.size();
192 
193  forAll(neiFaces, nfI)
194  {
195  if (neiFaces[nfI].size() && neiFaces[nfI][0] < minNei)
196  {
197  nextNei = nfI;
198  minNei = neiFaces[nfI][0];
199  }
200  }
201 
202  if (nextNei > -1)
203  {
204  // Add the face to the list of faces
205  edges[nEdges] = curEdges[nextNei];
206 
207  // Set face-edge and face-neighbour-edge to current face label
208  faceEdges[facei][nextNei] = nEdges;
209 
210  DynamicList<label>& cnf = neiFaces[nextNei];
211  DynamicList<label>& eonf = edgeOfNeiFace[nextNei];
212 
213  // Set edge-face addressing
214  labelList& curEf = edgeFaces[nEdges];
215  curEf.setSize(cnf.size() + 1);
216  curEf[0] = facei;
217 
218  forAll(cnf, cnfI)
219  {
220  faceEdges[cnf[cnfI]][eonf[cnfI]] = nEdges;
221 
222  curEf[cnfI + 1] = cnf[cnfI];
223  }
224 
225  // Stop the neighbour from being used again
226  cnf.clear();
227  eonf.clear();
228 
229  // Increment number of faces counter
230  nEdges++;
231  }
232  else
233  {
235  << "Error in internal edge insertion"
236  << abort(FatalError);
237  }
238  }
239  }
240 
241  nInternalEdges_ = nEdges;
242 
243  // Do boundary faces
244 
245  forAll(faceEdges, facei)
246  {
247  labelList& curEdges = faceEdges[facei];
248 
249  forAll(curEdges, edgeI)
250  {
251  if (curEdges[edgeI] < 0)
252  {
253  // Grab edge and faceEdge
254  edges[nEdges] = faceIntoEdges[facei][edgeI];
255  curEdges[edgeI] = nEdges;
256 
257  // Add edgeFace
258  labelList& curEf = edgeFaces[nEdges];
259  curEf.setSize(1);
260  curEf[0] = facei;
261 
262  nEdges++;
263  }
264  }
265  }
266 
267  // edges
268  edges.setSize(nEdges);
269 
270  // edgeFaces list
271  edgeFaces.setSize(nEdges);
272 
273  // faceFaces list
274  faceFacesPtr_.reset(new labelListList(locFcs.size()));
275  auto& faceFaces = *faceFacesPtr_;
276 
277  forAll(faceFaces, facei)
278  {
279  faceFaces[facei].transfer(ff[facei]);
280  }
281 
282 
283  DebugInFunction << "Calculated patch addressing" << nl;
284 }
285 
286 
287 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
PrimitivePatch.H
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::fv::ff
const FieldField< fvPatchField, Type > & ff(const FieldField< fvPatchField, Type > &bf)
Definition: CrankNicolsonDdtScheme.C:275
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
f
labelList f(nPoints)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
DynamicList.H
Foam::edgeListList
List< edgeList > edgeListList
A List of edgeList.
Definition: edgeList.H:65
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:85