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