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