patchToPoly2DMesh.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 | Copyright (C) 2011-2016 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 \*---------------------------------------------------------------------------*/
25 
26 #include "patchToPoly2DMesh.H"
27 #include "PatchTools.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::patchToPoly2DMesh::flipFaceOrder()
32 {
33  const edgeList& edges = patch_.edges();
34  const faceList& localFaces = patch_.localFaces();
35  const labelList& meshPoints = patch_.meshPoints();
36 
37  Info<< "Flipping face order if necessary." << endl;
38  forAll(edges, edgeI)
39  {
40  const edge& e = edges[edgeI];
41 
42  faces_[edgeI].setSize(2);
43 
44  label edgeOwner = owner_[edgeI];
45 
46  const face& f = localFaces[edgeOwner];
47 
48  label fp = findIndex(f, e[0]);
49 
50  if (f.nextLabel(fp) != e[1])
51  {
52  Info<< "Flipping face " << faces_[edgeI] << endl;
53  faces_[edgeI][0] = meshPoints[e[1]];
54  faces_[edgeI][1] = meshPoints[e[0]];
55  }
56  else
57  {
58  faces_[edgeI][0] = meshPoints[e[0]];
59  faces_[edgeI][1] = meshPoints[e[1]];
60  }
61  }
62 }
63 
64 
65 void Foam::patchToPoly2DMesh::createNeighbours()
66 {
67  const edgeList& edges = patch_.edges();
68  const labelListList& edgeFaces = patch_.edgeFaces();
69 
70  Info<< "Calculating neighbours." << endl;
71  forAll(edges, edgeI)
72  {
73  const labelList& eFaces = edgeFaces[edgeI];
74  if (eFaces.size() == 2)
75  {
76  if (owner_[edgeI] == eFaces[0])
77  {
78  neighbour_[edgeI] = eFaces[1];
79  }
80  else
81  {
82  neighbour_[edgeI] = eFaces[0];
83  }
84  }
85  else if (eFaces.size() == 1)
86  {
87  continue;
88  }
89  else
90  {
92  << abort(FatalError);
93  }
94  }
95 }
96 
97 
98 Foam::labelList Foam::patchToPoly2DMesh::internalFaceOrder()
99 {
100  const labelListList& faceEdges = patch_.faceEdges();
101 
102  labelList oldToNew(owner_.size(), -1);
103 
104  label newFacei = 0;
105 
106  forAll(faceEdges, facei)
107  {
108  const labelList& fEdges = faceEdges[facei];
109  // Neighbouring faces
110  SortableList<label> nbr(fEdges.size(), -1);
111 
112  forAll(fEdges, feI)
113  {
114  if (fEdges[feI] < neighbour_.size())
115  {
116  // Internal edge. Get the face on other side.
117 
118  label nbrFacei = neighbour_[fEdges[feI]];
119 
120  if (nbrFacei == facei)
121  {
122  nbrFacei = owner_[fEdges[feI]];
123  }
124 
125  if (facei < nbrFacei)
126  {
127  // facei is master
128  nbr[feI] = nbrFacei;
129  }
130  }
131  }
132 
133  nbr.sort();
134 
135  forAll(nbr, i)
136  {
137  if (nbr[i] != -1)
138  {
139  oldToNew[fEdges[nbr.indices()[i]]] = newFacei++;
140  }
141  }
142  }
143 
144  return oldToNew;
145 }
146 
147 
148 void Foam::patchToPoly2DMesh::addPatchFacesToFaces()
149 {
150  const labelList& meshPoints = patch_.meshPoints();
151 
152  label offset = patch_.nInternalEdges();
153  face f(2);
154 
155  forAll(patchNames_, patchi)
156  {
157  forAllConstIter(EdgeMap<label>, mapEdgesRegion_, eIter)
158  {
159  if (eIter() == patchi)
160  {
161  f[0] = meshPoints[eIter.key().start()];
162  f[1] = meshPoints[eIter.key().end()];
163  faces_[offset++] = f;
164  }
165  }
166  }
167 
168  f.clear();
169 }
170 
171 
172 void Foam::patchToPoly2DMesh::addPatchFacesToOwner()
173 {
174  const label nInternalEdges = patch_.nInternalEdges();
175  const faceList& faces = patch_.faces();
176  const label nExternalEdges = patch_.edges().size() - nInternalEdges;
177  const labelList& meshPoints = patch_.meshPoints();
178 
179  // Reorder patch faces on owner list.
180  labelList newOwner = owner_;
181 
182  label nMatched = 0;
183 
184  for
185  (
186  label bFacei = nInternalEdges;
187  bFacei < faces_.size();
188  ++bFacei
189  )
190  {
191  const face& e = faces_[bFacei];
192 
193  bool matched = false;
194 
195  for
196  (
197  label bEdgeI = nInternalEdges;
198  bEdgeI < faces_.size();
199  ++bEdgeI
200  )
201  {
202  if
203  (
204  e[0] == meshPoints[patch_.edges()[bEdgeI][0]]
205  && e[1] == meshPoints[patch_.edges()[bEdgeI][1]]
206  )
207  {
208  const face& f = faces[owner_[bEdgeI]];
209 
210  label fp = findIndex(f, e[0]);
211 
212  newOwner[bFacei] = owner_[bEdgeI];
213 
214  if (f.nextLabel(fp) != e[1])
215  {
216  Info<< "Flipping" << endl;
217 
218  faces_[bFacei][0] = e[1];
219  faces_[bFacei][1] = e[0];
220  }
221 
222  nMatched++;
223 
224  matched = true;
225  }
226  else if
227  (
228  e[0] == meshPoints[patch_.edges()[bEdgeI][1]]
229  && e[1] == meshPoints[patch_.edges()[bEdgeI][0]]
230  )
231  {
232  Info<< "Warning: Wrong orientation." << endl;
233  nMatched++;
234  matched = true;
235  }
236  }
237  if (!matched)
238  {
239  Info<< "No match for edge." << endl;
240  }
241  }
242 
243  if (nMatched != nExternalEdges)
244  {
245  Info<< "Number of matched edges, " << nMatched
246  << ", does not match number of external edges, "
247  << nExternalEdges << endl;
248  }
249 
250  owner_ = newOwner.xfer();
251 }
252 
253 
254 void Foam::patchToPoly2DMesh::createPolyMeshComponents()
255 {
256  flipFaceOrder();
257 
258  createNeighbours();
259 
260  // New function for returning a map of old faces to new faces.
261  labelList oldToNew = internalFaceOrder();
262 
263  inplaceReorder(oldToNew, faces_);
264  inplaceReorder(oldToNew, owner_);
265  inplaceReorder(oldToNew, neighbour_);
266 
267  // Add patches.
268  addPatchFacesToFaces();
269 
270  addPatchFacesToOwner();
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
275 
276 Foam::patchToPoly2DMesh::patchToPoly2DMesh
277 (
278  const MeshedSurface<face>& patch,
279  const wordList& patchNames,
280  const labelList& patchSizes,
281  const EdgeMap<label>& mapEdgesRegion
282 )
283 :
284  patch_(patch),
285  patchNames_(patchNames),
286  patchSizes_(patchSizes),
287  patchStarts_(patchNames.size(), 0),
288  mapEdgesRegion_(mapEdgesRegion),
289  points_(patch.points()),
290  faces_(patch.nEdges()),
291  owner_(PatchTools::edgeOwner(patch)),
292  neighbour_(patch.nInternalEdges())
293 {}
294 
295 
296 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
297 
299 {}
300 
301 
302 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
303 
305 {
306  for (label edgeI = 0; edgeI < patch_.nInternalEdges(); edgeI++)
307  {
308  if (patch_.edgeFaces()[edgeI].size() != 2)
309  {
311  << "internal edge:" << edgeI
312  << " patch.edgeFaces()[edgeI]:" << patch_.edgeFaces()[edgeI]
313  << abort(FatalError);
314  }
315  }
316 
317  for
318  (
319  label edgeI = patch_.nInternalEdges();
320  edgeI < patch_.nEdges();
321  edgeI++
322  )
323  {
324  if (patch_.edgeFaces()[edgeI].size() != 1)
325  {
327  << "boundary edge:" << edgeI
328  << " patch.edgeFaces()[edgeI]:" << patch_.edgeFaces()[edgeI]
329  << abort(FatalError);
330  }
331  }
332 
333  createPolyMeshComponents();
334 
335  label startFace = patch_.nInternalEdges();
336  forAll(patchNames_, patchi)
337  {
338  patchStarts_[patchi] = startFace;
339  startFace += patchSizes_[patchi];
340  }
341 }
342 
343 
344 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
const wordList & patchNames() const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
List< face > faceList
Definition: faceListFwd.H:43
const List< Face > & faces() const
Return const access to the faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const labelList & patchSizes() const
List< edge > edgeList
Definition: edgeList.H:38
const labelListList & faceEdges() const
Return face-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
~patchToPoly2DMesh()
Destructor.
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
labelList f(nPoints)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
List< word > wordList
A List of words.
Definition: fileName.H:54
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
void createMesh()
Create the mesh.
label nInternalEdges() const
Number of internal edges.
messageStream Info
const labelListList & edgeFaces() const
Return edge-face addressing.