addCellLayer.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) 2021 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
30#include "polyMesh.H"
31#include "primitiveMesh.H"
32#include "polyTopoChange.H"
33#include "polyTopoChanger.H"
34#include "polyAddPoint.H"
35#include "polyAddCell.H"
36#include "polyAddFace.H"
37#include "polyModifyFace.H"
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41Foam::tmp<Foam::vectorField> Foam::layerAdditionRemoval::extrusionDir() const
42{
43 const polyMesh& mesh = topoChanger().mesh();
44 const primitiveFacePatch& masterFaceLayer =
45 mesh.faceZones()[faceZoneID_.index()]();
46
47 const pointField& points = mesh.points();
48 const labelList& mp = masterFaceLayer.meshPoints();
49
50 tmp<vectorField> textrusionDir(new vectorField(mp.size()));
51 vectorField& extrusionDir = textrusionDir.ref();
52
53 if (setLayerPairing())
54 {
55 if (debug)
56 {
57 Pout<< "void layerAdditionRemoval::extrusionDir() const "
58 << " for object " << name() << " : "
59 << "Using edges for point insertion" << endl;
60 }
61
62 // Detected a valid layer. Grab the point and face collapse mapping
63 const labelList& ptc = pointsPairing();
64
65 forAll(extrusionDir, mpI)
66 {
67 extrusionDir[mpI] = points[ptc[mpI]] - points[mp[mpI]];
68 }
69 }
70 else
71 {
72 if (debug)
73 {
74 Pout<< "void layerAdditionRemoval::extrusionDir() const "
75 << " for object " << name() << " : "
76 << "A valid layer could not be found in front of "
77 << "the addition face layer. Using face-based "
78 << "point normals for point addition"
79 << endl;
80 }
81
82 extrusionDir = minLayerThickness_*masterFaceLayer.pointNormals();
83 }
84
85 return textrusionDir;
86}
87
88
89void Foam::layerAdditionRemoval::addCellLayer
90(
91 polyTopoChange& ref
92) const
93{
94 // Insert the layer addition instructions into the topological change
95
96 // Algorithm:
97 // 1. For every point in the master zone add a new point
98 // 2. For every face in the master zone add a cell
99 // 3. Go through all the edges of the master zone. For all internal edges,
100 // add a face with the correct orientation and make it internal.
101 // For all boundary edges, find the patch it belongs to and add the face
102 // to this patch
103 // 4. Create a set of new faces on the patch image and assign them to be
104 // between the old master cells and the newly created cells
105 // 5. Modify all the faces in the patch such that they are located
106 // between the old slave cells and newly created cells
107
108 if (debug)
109 {
110 Pout<< "void layerAdditionRemoval::addCellLayer("
111 << "polyTopoChange& ref) const for object " << name() << " : "
112 << "Adding cell layer" << endl;
113 }
114
115 // Create the points
116
117 const polyMesh& mesh = topoChanger().mesh();
118 const primitiveFacePatch& masterFaceLayer =
119 mesh.faceZones()[faceZoneID_.index()]();
120
121 const pointField& points = mesh.points();
122 const labelList& mp = masterFaceLayer.meshPoints();
123
124 // Get the extrusion direction for the added points
125
126 const vectorField pointOffsets(extrusionDir());
127
128 // Add the new points
129 labelList addedPoints(mp.size());
130
131 forAll(mp, pointi)
132 {
133 // Add the point nominal thickness away from the master zone point
134 // and grab the label
135 addedPoints[pointi] =
136 ref.setAction
137 (
138 polyAddPoint
139 (
140 points[mp[pointi]] // point
141 + addDelta_*pointOffsets[pointi],
142 mp[pointi], // master point
143 -1, // zone for point
144 true // supports a cell
145 )
146 );
147 }
148
149 if (debug > 1)
150 {
151 Pout<< "mp: " << mp << " addedPoints: " << addedPoints << endl;
152 }
153
154 // Create the cells
155
156 const labelList& mc =
157 mesh.faceZones()[faceZoneID_.index()].masterCells();
158 const labelList& sc =
159 mesh.faceZones()[faceZoneID_.index()].slaveCells();
160
161 const labelList& mf = mesh.faceZones()[faceZoneID_.index()];
162 const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
163
164 const labelList& own = mesh.faceOwner();
165 const labelList& nei = mesh.faceNeighbour();
166
167 labelList addedCells(mf.size());
168
169 forAll(mf, facei)
170 {
171 label celli = mc[facei];
172 label zoneI = mesh.cellZones().whichZone(celli);
173
174 addedCells[facei] =
175 ref.setAction
176 (
177 polyAddCell
178 (
179 -1, // master point
180 -1, // master edge
181 mf[facei], // master face
182 -1, // master cell
183 zoneI // zone for cell
184 )
185 );
186 }
187
188 // Create the new faces
189
190 // Grab the local faces from the master zone
191 const faceList& zoneFaces = masterFaceLayer.localFaces();
192
193 // Create the new faces by copying the master zone. All the added
194 // faces need to point into the newly created cells, which means
195 // all the faces from the master layer are flipped. The flux flip
196 // is determined from the original master layer face and the face
197 // owner: if the master cell is equal to the face owner the flux
198 // remains the same; otherwise it is flipped
199
200 forAll(zoneFaces, facei)
201 {
202 const face oldFace = zoneFaces[facei].reverseFace();
203
204 face newFace(oldFace.size());
205
206 forAll(oldFace, pointi)
207 {
208 newFace[pointi] = addedPoints[oldFace[pointi]];
209 }
210
211 bool flipFaceFlux = false;
212
213 // Flip the face as necessary
214 if
215 (
216 !mesh.isInternalFace(mf[facei])
217 || mc[facei] == nei[mf[facei]]
218 )
219 {
220 flipFaceFlux = true;
221 }
222
223 // Add the face
224 ref.setAction
225 (
226 polyAddFace
227 (
228 newFace, // face
229 mc[facei], // owner
230 addedCells[facei], // neighbour
231 -1, // master point
232 -1, // master edge
233 mf[facei], // master face for addition
234 flipFaceFlux, // flux flip
235 -1, // patch for face
236 -1, // zone for face
237 false // face zone flip
238 )
239 );
240
241 if (debug > 1)
242 {
243 Pout<< "adding face: " << newFace
244 << " own: " << mc[facei]
245 << " nei: " << addedCells[facei]
246 << endl;
247 }
248 }
249
250 // Modify the faces from the master zone for the new neighbour
251
252 const faceList& faces = mesh.faces();
253
254 // Pout<< "mfFlip: " << mfFlip << endl;
255
256 forAll(mf, facei)
257 {
258 const label curfaceID = mf[facei];
259
260 // If the face is internal, modify its owner to be the newly
261 // created cell. No flip is necessary
262 if (!mesh.isInternalFace(curfaceID))
263 {
264 ref.setAction
265 (
266 polyModifyFace
267 (
268 faces[curfaceID], // modified face
269 curfaceID, // label of face being modified
270 addedCells[facei], // owner
271 -1, // neighbour
272 false, // face flip
273 mesh.boundaryMesh().whichPatch(curfaceID),// patch for face
274 false, // remove from zone
275 faceZoneID_.index(), // zone for face
276 mfFlip[facei] // face flip in zone
277 )
278 );
279
280 if (debug > 1)
281 {
282 Pout<< "Modifying a boundary face. Face: " << curfaceID
283 << " flip: " << mfFlip[facei]
284 << endl;
285 }
286 }
287
288 // If slave cell is owner, the face remains the same (but with
289 // a new neighbour - the newly created cell). Otherwise, the
290 // face is flipped.
291 else if (sc[facei] == own[curfaceID])
292 {
293 // Orientation is good, just change neighbour
294 ref.setAction
295 (
296 polyModifyFace
297 (
298 faces[curfaceID], // modified face
299 curfaceID, // label of face being modified
300 own[curfaceID], // owner
301 addedCells[facei], // neighbour
302 false, // face flip
303 mesh.boundaryMesh().whichPatch(curfaceID),// patch for face
304 false, // remove from zone
305 faceZoneID_.index(), // zone for face
306 mfFlip[facei] // face flip in zone
307 )
308 );
309
310 if (debug > 1)
311 {
312 Pout<< "modify face, no flip " << curfaceID
313 << " own: " << own[curfaceID]
314 << " nei: " << addedCells[facei]
315 << endl;
316 }
317 }
318 else
319 {
320 // Change in orientation; flip face
321 ref.setAction
322 (
323 polyModifyFace
324 (
325 faces[curfaceID].reverseFace(), // modified face
326 curfaceID, // label of face being modified
327 nei[curfaceID], // owner
328 addedCells[facei], // neighbour
329 true, // face flip
330 mesh.boundaryMesh().whichPatch(curfaceID), // patch for face
331 false, // remove from zone
332 faceZoneID_.index(), // zone for face
333 !mfFlip[facei] // face flip in zone
334 )
335 );
336
337 if (debug > 1)
338 {
339 Pout<< "modify face, with flip " << curfaceID
340 << " own: " << own[curfaceID]
341 << " nei: " << addedCells[facei]
342 << endl;
343 }
344 }
345 }
346
347 // Create new faces from edges
348
349 const edgeList& zoneLocalEdges = masterFaceLayer.edges();
350
351 const labelListList& edgeFaces = masterFaceLayer.edgeFaces();
352
353 label nInternalEdges = masterFaceLayer.nInternalEdges();
354
355 const labelList& meshEdges =
356 mesh.faceZones()[faceZoneID_.index()].meshEdges();
357
358 // Do all internal edges
359 for (label curEdgeID = 0; curEdgeID < nInternalEdges; curEdgeID++)
360 {
361 face newFace(4);
362
363 newFace[0] = mp[zoneLocalEdges[curEdgeID].start()];
364 newFace[1] = mp[zoneLocalEdges[curEdgeID].end()];
365 newFace[2] = addedPoints[zoneLocalEdges[curEdgeID].end()];
366 newFace[3] = addedPoints[zoneLocalEdges[curEdgeID].start()];
367
368 ref.setAction
369 (
370 polyAddFace
371 (
372 newFace, // face
373 addedCells[edgeFaces[curEdgeID][0]], // owner
374 addedCells[edgeFaces[curEdgeID][1]], // neighbour
375 -1, // master point
376 meshEdges[curEdgeID], // master edge
377 -1, // master face
378 false, // flip flux
379 -1, // patch for face
380 -1, // zone for face
381 false // face zone flip
382 )
383 );
384
385 if (debug > 1)
386 {
387 Pout<< "Add internal face off edge: " << newFace
388 << " own: " << addedCells[edgeFaces[curEdgeID][0]]
389 << " nei: " << addedCells[edgeFaces[curEdgeID][1]]
390 << endl;
391 }
392 }
393
394 // Prepare creation of faces from boundary edges.
395 // Note:
396 // A tricky part of the algorithm is finding the patch into which the
397 // newly created face will be added. For this purpose, take the edge
398 // and grab all the faces coming from it. From the set of faces
399 // eliminate the internal faces and find the boundary face from the rest.
400 // If there is more than one boundary face (excluding the ones in
401 // the master zone), the patch with the lower label is selected.
402
403 const labelListList& meshEdgeFaces = mesh.edgeFaces();
404
405 const faceZoneMesh& zoneMesh = mesh.faceZones();
406
407 // Do all boundary edges
408
409 for
410 (
411 label curEdgeID = nInternalEdges;
412 curEdgeID < zoneLocalEdges.size();
413 curEdgeID++
414 )
415 {
416 face newFace(4);
417 newFace[0] = mp[zoneLocalEdges[curEdgeID].start()];
418 newFace[1] = mp[zoneLocalEdges[curEdgeID].end()];
419 newFace[2] = addedPoints[zoneLocalEdges[curEdgeID].end()];
420 newFace[3] = addedPoints[zoneLocalEdges[curEdgeID].start()];
421
422 // Determine the patch for insertion
423 const labelList& curFaces = meshEdgeFaces[meshEdges[curEdgeID]];
424
425 label patchID = -1;
426 label zoneID = -1;
427 label extrudeFaceID = -1;
428
429 forAll(curFaces, curFacei)
430 {
431 extrudeFaceID = curFaces[curFacei];
432
433 if (!mesh.isInternalFace(extrudeFaceID))
434 {
435 // Face not internal. Check if it is not in the zone
436 const label cfZone = zoneMesh.whichZone(extrudeFaceID);
437 if (cfZone != faceZoneID_.index())
438 {
439 // Found the face in a boundary patch which is not in zone
440 patchID = mesh.boundaryMesh().whichPatch(extrudeFaceID);
441 zoneID = cfZone;
442
443 break;
444 }
445 }
446 }
447
448 if (patchID < 0)
449 {
451 << "Cannot find patch for edge " << meshEdges[curEdgeID]
452 << ". Edge: " << mesh.edges()[meshEdges[curEdgeID]]
453 << abort(FatalError);
454 }
455
456 ref.setAction
457 (
458 polyAddFace
459 (
460 newFace, // face
461 addedCells[edgeFaces[curEdgeID][0]], // owner
462 -1, // neighbour
463 -1, // master point
464 meshEdges[curEdgeID], // master edge
465 -1, // master face
466 false, // flip flux
467 patchID, // patch for face
468 zoneID, // zone
469 false // zone face flip
470 )
471 );
472
473 if (debug > 1)
474 {
475 Pout<< "add boundary face: " << newFace
476 << " into patch " << patchID
477 << " own: " << addedCells[edgeFaces[curEdgeID][0]]
478 << endl;
479 }
480
481 // Handle duplicate boundary faces (for e.g. cyclicACMI)
482 if (patchID != -1)
483 {
484 for (const label cf : curFaces)
485 {
486 if (cf != extrudeFaceID && !mesh.isInternalFace(cf))
487 {
488 // Check if it is not in the zone and duplicate of
489 // extrudeFaceID
490 const label cfZone = zoneMesh.whichZone(cf);
491 if
492 (
493 cfZone != faceZoneID_.index()
494 && face::compare(faces[cf], faces[extrudeFaceID]) == 1
495 )
496 {
497 // Found the face in a boundary patch which is not
498 // in zone (so would not be extruded above)
500 zoneID = cfZone;
501
502 ref.setAction
503 (
504 polyAddFace
505 (
506 newFace, // face
507 addedCells[edgeFaces[curEdgeID][0]],// owner
508 -1, // neighbour
509 -1, // point
510 meshEdges[curEdgeID], // edge
511 -1, // face
512 false, // flip flux
513 patchID, // patch
514 zoneID, // zone
515 false // zone flip
516 )
517 );
518
519 if (debug > 1)
520 {
521 Pout<< "add duplicate boundary face: " << newFace
522 << " into patch " << patchID
523 << " own: "
524 << addedCells[edgeFaces[curEdgeID][0]] << endl;
525 }
526 }
527 }
528 }
529 }
530 }
531
532 // Modify the remaining faces of the master cells to reconnect to the new
533 // layer of faces.
534 // Algorithm: Go through all the cells of the master zone and make
535 // a map of faces to avoid duplicates. Do not insert the faces in
536 // the master patch (as they have already been dealt with). Make
537 // a master layer point renumbering map, which for every point in
538 // the master layer gives its new label. Loop through all faces in
539 // the map and attempt to renumber them using the master layer
540 // point renumbering map. Once the face is renumbered, compare it
541 // with the original face; if they are the same, the face has not
542 // changed; if not, modify the face but keep all of its old
543 // attributes (apart from the vertex numbers).
544
545 // Create the map of faces in the master cell layer
546 labelHashSet masterCellFaceMap(primitiveMesh::facesPerCell_*mc.size());
547
548 const cellList& cells = mesh.cells();
549
550 forAll(mc, celli)
551 {
552 const labelList& curFaces = cells[mc[celli]];
553
554 forAll(curFaces, facei)
555 {
556 // Check if the face belongs to the master zone; if not add it
557 if (zoneMesh.whichZone(curFaces[facei]) != faceZoneID_.index())
558 {
559 masterCellFaceMap.insert(curFaces[facei]);
560 }
561 }
562 }
563
564 // Create the master layer point map
565 Map<label> masterLayerPointMap(2*mp.size());
566
567 forAll(mp, pointi)
568 {
569 masterLayerPointMap.insert
570 (
571 mp[pointi],
572 addedPoints[pointi]
573 );
574 }
575
576 // Grab the list of faces of the master layer
577 const labelList masterCellFaces = masterCellFaceMap.toc();
578
579 forAll(masterCellFaces, facei)
580 {
581 // Attempt to renumber the face using the masterLayerPointMap.
582 // Missing point remain the same
583
584 const label curFaceID = masterCellFaces[facei];
585
586 const face& oldFace = faces[curFaceID];
587
588 face newFace(oldFace.size());
589
590 bool changed = false;
591
592 forAll(oldFace, pointi)
593 {
594 if (masterLayerPointMap.found(oldFace[pointi]))
595 {
596 changed = true;
597
598 newFace[pointi] = masterLayerPointMap.find(oldFace[pointi])();
599 }
600 else
601 {
602 newFace[pointi] = oldFace[pointi];
603 }
604 }
605
606 // If the face has changed, create a modification entry
607 if (changed)
608 {
609 // Get face zone and its flip
610 label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
611 bool modifiedFaceZoneFlip = false;
612
613 if (modifiedFaceZone >= 0)
614 {
615 modifiedFaceZoneFlip =
616 mesh.faceZones()[modifiedFaceZone].flipMap()
617 [
618 mesh.faceZones()[modifiedFaceZone].whichFace(curFaceID)
619 ];
620 }
621
622 if (mesh.isInternalFace(curFaceID))
623 {
624 ref.setAction
625 (
626 polyModifyFace
627 (
628 newFace, // modified face
629 curFaceID, // label of face being modified
630 own[curFaceID], // owner
631 nei[curFaceID], // neighbour
632 false, // face flip
633 -1, // patch for face
634 false, // remove from zone
635 modifiedFaceZone, // zone for face
636 modifiedFaceZoneFlip // face flip in zone
637 )
638 );
639
640 if (debug > 1)
641 {
642 Pout<< "modifying stick-out face. Internal Old face: "
643 << oldFace
644 << " new face: " << newFace
645 << " own: " << own[curFaceID]
646 << " nei: " << nei[curFaceID]
647 << endl;
648 }
649 }
650 else
651 {
652 ref.setAction
653 (
654 polyModifyFace
655 (
656 newFace, // modified face
657 curFaceID, // label of face being modified
658 own[curFaceID], // owner
659 -1, // neighbour
660 false, // face flip
661 mesh.boundaryMesh().whichPatch(curFaceID),
662 // patch for face
663 false, // remove from zone
664 modifiedFaceZone, // zone for face
665 modifiedFaceZoneFlip // face flip in zone
666 )
667 );
668
669 if (debug > 1)
670 {
671 Pout<< "modifying stick-out face. Boundary Old face: "
672 << oldFace
673 << " new face: " << newFace
674 << " own: " << own[curFaceID]
675 << " patch: "
676 << mesh.boundaryMesh().whichPatch(curFaceID)
677 << endl;
678 }
679 }
680 }
681 }
682
683 if (debug)
684 {
685 Pout<< "void layerAdditionRemoval::addCellLayer(polyTopoChange&) const "
686 << " for object " << name() << ": "
687 << "Finished adding cell layer" << endl;
688 }
689}
690
691
692// ************************************************************************* //
label index() const
The index of the first matching items, -1 if no matches.
Definition: DynamicID.H:123
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:288
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:281
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
const word & name() const
Return name of this modifier.
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:498
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const polyMesh & mesh() const
Return the mesh reference.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
static const unsigned facesPerCell_
Estimated number of faces per cell.
const labelListList & edgeFaces() const
const cellList & cells() const
A class for managing temporary objects.
Definition: tmp.H:65
rDeltaT ref()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
const cellShapeList & cells
const labelIOList & zoneID
const dimensionedScalar mp
Proton mass.
List< label > labelList
A List of labels.
Definition: List.H:66
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition: List.H:64
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
error FatalError
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333