removeCellLayer.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) 2018-2019 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 \*---------------------------------------------------------------------------*/
28 
29 #include "layerAdditionRemoval.H"
30 #include "polyMesh.H"
31 #include "primitiveMesh.H"
32 #include "polyTopoChange.H"
33 #include "oppositeFace.H"
34 #include "polyTopoChanger.H"
35 #include "polyRemoveCell.H"
36 #include "polyRemoveFace.H"
37 #include "polyRemovePoint.H"
38 #include "polyModifyFace.H"
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 bool Foam::layerAdditionRemoval::validCollapse() const
43 {
44  // Check for valid layer collapse
45  // - no boundary-to-boundary collapse
46 
47  if (debug)
48  {
49  Pout<< "Checking layer collapse for object " << name() << endl;
50  }
51 
52  // Grab the face collapse mapping
53  const polyMesh& mesh = topoChanger().mesh();
54 
55  const labelList& ftc = facesPairing();
56  const labelList& mf = mesh.faceZones()[faceZoneID_.index()];
57 
58  label nBoundaryHits = 0;
59 
60  forAll(mf, facei)
61  {
62  if
63  (
64  !mesh.isInternalFace(mf[facei])
65  && !mesh.isInternalFace(ftc[facei])
66  )
67  {
68  nBoundaryHits++;
69  }
70  }
71 
72 
73  if (debug)
74  {
75  Pout<< "Finished checking layer collapse for object "
76  << name() <<". Number of boundary-on-boundary hits: "
77  << nBoundaryHits << endl;
78  }
79 
80  if (returnReduce(nBoundaryHits, sumOp<label>()) > 0)
81  {
82  return false;
83  }
84  else
85  {
86  return true;
87  }
88 }
89 
90 
91 void Foam::layerAdditionRemoval::removeCellLayer
92 (
93  polyTopoChange& ref
94 ) const
95 {
96  // Algorithm for layer removal. Second phase: topological change
97  // 2) Go through all the faces of the master cell layer and remove
98  // the ones that are not in the master face zone.
99  //
100  // 3) Grab all the faces coming out of points that are collapsed
101  // and select the ones that are not marked for removal. Go
102  // through the remaining faces and replace the point to remove by
103  // the equivalent point in the master face zone.
104  if (debug)
105  {
106  Pout<< "Removing the cell layer for object " << name() << endl;
107  }
108 
109  const polyMesh& mesh = topoChanger().mesh();
110 
111  const labelList& own = mesh.faceOwner();
112  const labelList& nei = mesh.faceNeighbour();
113 
114  const labelList& ptc = pointsPairing();
115  const labelList& ftc = facesPairing();
116 
117  // Remove all the cells from the master layer
118  const labelList& mc =
119  mesh.faceZones()[faceZoneID_.index()].masterCells();
120 
121  forAll(mc, facei)
122  {
123  label slaveSideCell = own[ftc[facei]];
124 
125  if (mesh.isInternalFace(ftc[facei]) && slaveSideCell == mc[facei])
126  {
127  // Owner cell of the face is being removed.
128  // Grab the neighbour instead
129  slaveSideCell = nei[ftc[facei]];
130  }
131 
132  ref.setAction(polyRemoveCell(mc[facei], slaveSideCell));
133  }
134 
135  // Remove all the faces from the master layer cells which are not in
136  // the master face layer
137  labelHashSet facesToRemoveMap(mc.size()*primitiveMesh::facesPerCell_);
138 
139  const cellList& cells = mesh.cells();
140 
141  forAll(mc, celli)
142  {
143  const cell& curCell = cells[mc[celli]];
144 
145  forAll(curCell, facei)
146  {
147  // Check if the face is in the master zone. If not, remove it
148  if
149  (
150  mesh.faceZones().whichZone(curCell[facei])
151  != faceZoneID_.index()
152  )
153  {
154  facesToRemoveMap.insert(curCell[facei]);
155  }
156  }
157  }
158 
159  for (const label facei : facesToRemoveMap)
160  {
161  ref.setAction(polyRemoveFace(facei));
162  }
163 
164  // Remove all points that will be collapsed
165  forAll(ptc, pointi)
166  {
167  ref.setAction(polyRemovePoint(ptc[pointi]));
168  }
169 
170  // Grab all faces coming off points to be deleted. If the face
171  // has not been deleted, replace the removed point with the
172  // equivalent from the master layer.
173 
174  // Make a map of all point to be removed, giving the master point label
175  // for each of them
176 
177  Map<label> removedPointMap(2*ptc.size());
178 
179  const labelList& meshPoints =
180  mesh.faceZones()[faceZoneID_.index()]().meshPoints();
181 
182  forAll(ptc, pointi)
183  {
184  removedPointMap.insert(ptc[pointi], meshPoints[pointi]);
185  }
186 
187  const labelListList& pf = mesh.pointFaces();
188 
189  const faceList& faces = mesh.faces();
190 
191  // Make a list of faces to be modified using the map to avoid duplicates
192  labelHashSet facesToModify(ptc.size()*primitiveMesh::facesPerPoint_);
193 
194  forAll(ptc, pointi)
195  {
196  const labelList& curFaces = pf[ptc[pointi]];
197 
198  forAll(curFaces, facei)
199  {
200  if (!facesToRemoveMap.found(curFaces[facei]))
201  {
202  facesToModify.insert(curFaces[facei]);
203  }
204  }
205  }
206 
207  labelList ftm = facesToModify.toc();
208 
209  if (debug > 1)
210  {
211  Pout<< "faces to modify: " << ftm << endl;
212  }
213 
214  forAll(ftm, facei)
215  {
216  // For every face to modify, copy the face and re-map the vertices.
217  // It is known all the faces will be changed since they hang off
218  // re-mapped vertices
219  label curFaceID = ftm[facei];
220 
221  face newFace(faces[curFaceID]);
222 
223  forAll(newFace, pointi)
224  {
225  const auto rpmIter = removedPointMap.cfind(newFace[pointi]);
226 
227  if (rpmIter.found())
228  {
229  // Point mapped. Replace it
230  newFace[pointi] = rpmIter();
231  }
232  }
233 
234  if (debug > 1)
235  {
236  Pout<< "face label: " << curFaceID
237  << " old face: " << faces[curFaceID]
238  << " new face: " << newFace << endl;
239  }
240 
241  // Get face zone and its flip
242  label modifiedFaceZone = mesh.faceZones().whichZone(curFaceID);
243  bool modifiedFaceZoneFlip = false;
244 
245  if (modifiedFaceZone >= 0)
246  {
247  const faceZone& fz = mesh.faceZones()[modifiedFaceZone];
248  modifiedFaceZoneFlip = fz.flipMap()[fz.whichFace(curFaceID)];
249  }
250 
251  label newNeighbour = -1;
252 
253  if (curFaceID < mesh.nInternalFaces())
254  {
255  newNeighbour = nei[curFaceID];
256  }
257 
258  // Modify the face
259  ref.setAction
260  (
261  polyModifyFace
262  (
263  newFace, // modified face
264  curFaceID, // label of face being modified
265  own[curFaceID], // owner
266  newNeighbour, // neighbour
267  false, // face flip
268  mesh.boundaryMesh().whichPatch(curFaceID),// patch for face
269  false, // remove from zone
270  modifiedFaceZone, // zone for face
271  modifiedFaceZoneFlip // face flip in zone
272  )
273  );
274  }
275 
276  // Modify the faces in the master layer to point past the removed cells
277 
278  const labelList& mf = mesh.faceZones()[faceZoneID_.index()];
279  const boolList& mfFlip = mesh.faceZones()[faceZoneID_.index()].flipMap();
280 
281  forAll(mf, facei)
282  {
283  // Grab the owner and neighbour of the faces to be collapsed and get rid
284  // of the cell to be removed
285  label masterSideCell = own[mf[facei]];
286 
287  if (masterSideCell == mc[facei])
288  {
289  if (mesh.isInternalFace(mf[facei]))
290  {
291  // Owner cell of the face is being removed.
292  // Grab the neighbour instead
293  masterSideCell = nei[mf[facei]];
294  }
295  else
296  {
297  masterSideCell = -1;
298  }
299  }
300 
301  label slaveSideCell = own[ftc[facei]];
302 
303  if (slaveSideCell == mc[facei])
304  {
305  if (mesh.isInternalFace(ftc[facei]))
306  {
307  // Owner cell of the face is being removed.
308  // Grab the neighbour instead
309  slaveSideCell = nei[ftc[facei]];
310  }
311  else
312  {
313  slaveSideCell = -1;
314  }
315  }
316 
317  // Find out if the face needs to be flipped
318  label newOwner = -1;
319  label newNeighbour = -1;
320  bool flipFace = false;
321  label newPatchID = -1;
322  label newZoneID = -1;
323 
324  // Is any of the faces a boundary face? If so, grab the patch
325  // A boundary-to-boundary collapse is checked for in validCollapse()
326  // and cannot happen here.
327 
328  if (!mesh.isInternalFace(mf[facei]))
329  {
330  // Master is the boundary face: it gets a new owner but no flip
331  newOwner = slaveSideCell;
332  newNeighbour = -1;
333  flipFace = false;
334  newPatchID = mesh.boundaryMesh().whichPatch(mf[facei]);
335  newZoneID = mesh.faceZones().whichZone(mf[facei]);
336  }
337  else if (!mesh.isInternalFace(ftc[facei]))
338  {
339  // Slave is the boundary face: grab its patch
340  newOwner = slaveSideCell;
341  newNeighbour = -1;
342 
343  // Find out if the face flip is necessary
344  if (own[mf[facei]] == slaveSideCell)
345  {
346  flipFace = false;
347  }
348  else
349  {
350  flipFace = true;
351  }
352 
353  newPatchID = mesh.boundaryMesh().whichPatch(ftc[facei]);
354 
355  // The zone of the master face is preserved
356  newZoneID = mesh.faceZones().whichZone(mf[facei]);
357  }
358  else
359  {
360  // Both faces are internal: flip is decided based on which of the
361  // new cells around it has a lower label
362  newOwner = min(masterSideCell, slaveSideCell);
363  newNeighbour = max(masterSideCell, slaveSideCell);
364 
365  if (newOwner == own[mf[facei]] || newNeighbour == nei[mf[facei]])
366  {
367  flipFace = false;
368  }
369  else
370  {
371  flipFace = true;
372  }
373 
374  newPatchID = -1;
375 
376  // The zone of the master face is preserved
377  newZoneID = mesh.faceZones().whichZone(mf[facei]);
378  }
379 
380  // Modify the face and flip if necessary
381  face newFace = faces[mf[facei]];
382  bool zoneFlip = mfFlip[facei];
383 
384  if (flipFace)
385  {
386  newFace.flip();
387  zoneFlip = !zoneFlip;
388  }
389 
390  if (debug > 1)
391  {
392  Pout<< "Modifying face " << mf[facei]
393  << " newFace: " << newFace << nl
394  << " newOwner: " << newOwner
395  << " newNeighbour: " << newNeighbour
396  << " flipFace: " << flipFace
397  << " newPatchID: " << newPatchID
398  << " newZoneID: " << newZoneID << nl
399  << " oldOwn: " << own[mf[facei]];
400  if (newPatchID == -1)
401  {
402  Pout<< " oldNei: " << nei[mf[facei]];
403  }
404  Pout<< endl;
405  }
406 
407  ref.setAction
408  (
409  polyModifyFace
410  (
411  newFace, // modified face
412  mf[facei], // label of face being modified
413  newOwner, // owner
414  newNeighbour, // neighbour
415  flipFace, // flip
416  newPatchID, // patch for face
417  false, // remove from zone
418  newZoneID, // new zone
419  zoneFlip // face flip in zone
420  )
421  );
422  }
423 }
424 
425 
426 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
oppositeFace.H
polyRemovePoint.H
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::primitiveMesh::facesPerPoint_
static const unsigned facesPerPoint_
Estimated number of faces per point.
Definition: primitiveMesh.H:410
Foam::polyMeshModifier::topoChanger
const polyTopoChanger & topoChanger() const
Return reference to morph engine.
Definition: polyMeshModifier.C:63
polyRemoveFace.H
polyTopoChanger.H
polyTopoChange.H
Foam::primitiveMesh::pointFaces
const labelListList & pointFaces() const
Definition: primitiveMeshPointFaces.C:34
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
primitiveMesh.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
ref
rDeltaT ref()
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
polyMesh.H
layerAdditionRemoval.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::primitiveMesh::facesPerCell_
static const unsigned facesPerCell_
Estimated number of faces per cell.
Definition: primitiveMesh.H:404
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::cellList
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:289
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
polyRemoveCell.H
polyModifyFace.H
Foam::polyMeshModifier::name
const word & name() const
Return name of this modifier.
Definition: polyMeshModifier.H:150
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::polyTopoChanger::mesh
const polyMesh & mesh() const
Return the mesh reference.
Definition: polyTopoChanger.H:111
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::DynamicID::index
label index() const
The index of the first matching items, -1 if no matches.
Definition: DynamicID.H:123