checkNonMappableCellConnections.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh 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  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "polyMeshGenModifier.H"
30 #include "helperFunctions.H"
31 #include "meshSurfaceEngine.H"
32 
33 # ifdef USE_OMP
34 #include <omp.h>
35 # endif
36 
37 //#define DEBUGCheck
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 // * * * * * * * * * * Private member functions * * * * * * * * * * * * * * * //
45 
47 {
48  const faceListPMG& faces = mesh_.faces();
49  const cellListPMG& cells = mesh_.cells();
50  const labelList& owner = mesh_.owner();
51 
52  cellType_.setSize(cells.size());
54 
55  //- find boundary cells
56  const PtrList<boundaryPatch>& boundaries = mesh_.boundaries();
57  forAll(boundaries, patchI)
58  {
59  const label start = boundaries[patchI].patchStart();
60  const label end = start + boundaries[patchI].patchSize();
61 
62  for(label faceI=start;faceI<end;++faceI)
63  cellType_[owner[faceI]] = BNDCELL;
64  }
65 
66  //- find boundary cells with all vertices at the boundary
68  const labelList& bp = mse.bp();
69 
70  # ifdef USE_OMP
71  # pragma omp parallel for schedule(dynamic, 1000)
72  # endif
73  for(label cellI=cells.size()-1;cellI>=0;--cellI)
74  {
75  if( cellType_[cellI] & INTERNALCELL )
76  continue;
77 
78  const cell& c = cells[cellI];
79 
80  //- mark boundary cells with all vertices at the boundary
81  const labelList cellPoints = c.labels(faces);
82  bool allBoundary(true);
83  forAll(cellPoints, cpI)
84  {
85  if( bp[cellPoints[cpI]] < 0 )
86  {
87  allBoundary = false;
88  break;
89  }
90  }
91 
92  if( allBoundary )
93  {
94  cellType_[cellI] |= ALLBNDVERTEXCELL;
95  }
96  else
97  {
98  continue;
99  }
100 
101  //- check if the internal faces are connected into a single group
102  //- over their edges
103  DynList<label> internalFaces;
104  forAll(c, fI)
105  {
106  if( c[fI] < mesh_.nInternalFaces() )
107  {
108  internalFaces.append(c[fI]);
109  }
110  else if( mesh_.faceIsInProcPatch(c[fI]) != -1 )
111  {
112  internalFaces.append(c[fI]);
113  }
114  }
115 
116  Map<label> faceGroup(internalFaces.size());
117  label nGroup(0);
118  forAll(internalFaces, i)
119  {
120  if( faceGroup.found(internalFaces[i]) )
121  continue;
122 
123  DynList<label> front;
124  front.append(internalFaces[i]);
125  faceGroup.insert(internalFaces[i], nGroup);
126 
127  while( front.size() )
128  {
129  const label fLabel = front.removeLastElement();
130 
131  forAll(internalFaces, j)
132  {
133  const label nei = internalFaces[j];
134  if( faceGroup.found(nei) )
135  continue;
136 
137  if( help::shareAnEdge(faces[fLabel], faces[nei]) )
138  {
139  front.append(nei);
140  faceGroup.insert(nei, nGroup);
141  }
142  }
143  }
144 
145  ++nGroup;
146  }
147 
148  if( nGroup > 1 )
149  cellType_[cellI] |= INTERNALFACEGROUP;
150  }
151 }
152 
153 // * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * * * //
154 // Constructors
155 
157 (
159 )
160 :
161  mesh_(mesh),
162  cellType_()
163 {
164 }
165 
166 // * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * * * //
167 
169 {
170 }
171 
172 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 
175 {
176  badCells.clear();
177 
178  //- classify cell types
179  findCellTypes();
180 
181  //- select ALLBNDVERTEXCELL and INTERNALFACEGROUP cells
182  //- with at least one INTERNALCELL neighbour
183  //- these cells do not need to stay in the mesh
184  const cellListPMG& cells = mesh_.cells();
185  const labelList& owner = mesh_.owner();
186  const labelList& neighbour = mesh_.neighbour();
187  const label nInternalFaces = mesh_.nInternalFaces();
188  const PtrList<processorBoundaryPatch>& procBoundaries = mesh_.procBoundaries();
189 
190  labelListList otherProcType;
191  if( Pstream::parRun() )
192  {
193  //- exchange cell types at processor boundaries
194  otherProcType.setSize(procBoundaries.size());
195 
196  //- send data to other processors
197  forAll(procBoundaries, patchI)
198  {
199  label start = procBoundaries[patchI].patchStart();
200  labelList patchCellType(procBoundaries[patchI].patchSize());
201 
202  forAll(patchCellType, faceI)
203  patchCellType[faceI] = cellType_[owner[start++]];
204 
205  OPstream toOtherProc
206  (
208  procBoundaries[patchI].neiProcNo(),
209  patchCellType.byteSize()
210  );
211 
212  toOtherProc << patchCellType;
213  }
214 
215  //- receive data from other processors
216  forAll(procBoundaries, patchI)
217  {
218  IPstream fromOtherProc
219  (
221  procBoundaries[patchI].neiProcNo()
222  );
223 
224  labelList& otherTypes = otherProcType[patchI];
225  fromOtherProc >> otherTypes;
226  }
227  }
228 
229  # ifdef USE_OMP
230  # pragma omp parallel for schedule(dynamic, 40)
231  # endif
232  for(label cellI=cellType_.size()-1;cellI>=0;--cellI)
233  {
234  if( cellType_[cellI] & INTERNALFACEGROUP )
235  {
236  # ifdef USE_OMP
237  # pragma omp critical
238  # endif
239  badCells.insert(cellI);
240  }
241  else if( cellType_[cellI] & (ALLBNDVERTEXCELL+INTERNALFACEGROUP) )
242  {
243  //- mark cells which have only one internal neighbour
244  const cell& c = cells[cellI];
245 
246  bool hasInternalNeighbour(false);
247  label nNeiCells(0);
248 
249  forAll(c, fI)
250  {
251  const label faceI = c[fI];
252 
253  if( faceI < nInternalFaces )
254  {
255  ++nNeiCells;
256 
257  label nei = neighbour[c[fI]];
258  if( nei == cellI )
259  nei = owner[c[fI]];
260 
261  if( cellType_[nei] & INTERNALCELL )
262  {
263  hasInternalNeighbour = true;
264  break;
265  }
266  }
267  else if( mesh_.faceIsInProcPatch(faceI) != -1 )
268  {
269  ++nNeiCells;
270 
271  const label patchI = mesh_.faceIsInProcPatch(faceI);
272  const label j = faceI - procBoundaries[patchI].patchStart();
273 
274  if( otherProcType[patchI][j] & INTERNALCELL )
275  {
276  hasInternalNeighbour = true;
277  break;
278  }
279  }
280  }
281 
282  if( hasInternalNeighbour || (nNeiCells == 1) )
283  {
284  # ifdef USE_OMP
285  # pragma omp critical
286  # endif
287  badCells.insert(cellI);
288  }
289  }
290  }
291 }
292 
294 {
295  labelHashSet badCells;
296 
297  label nRemoved;
298  bool changed(false);
299 
300  do
301  {
302  findCells(badCells);
303 
304  nRemoved = badCells.size();
305  reduce(nRemoved, sumOp<label>());
306 
307  Info << "Found " << nRemoved << " non-mappable cells" << endl;
308 
309  if( nRemoved != 0 )
310  {
311  boolList removeCell(mesh_.cells().size(), false);
312  forAllConstIter(labelHashSet, badCells, it)
313  removeCell[it.key()] = true;
314 
316 
317  changed = true;
318  }
319  } while( nRemoved );
320 
321  return changed;
322 }
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 } // End namespace Foam
327 
328 // ************************************************************************* //
Foam::polyMeshGenFaces::neighbour
const labelList & neighbour() const
Definition: polyMeshGenFacesI.H:86
Foam::polyMeshGenFaces::owner
const labelList & owner() const
owner and neighbour cells for faces
Definition: polyMeshGenFacesI.H:67
Foam::meshSurfaceEngine::bp
const labelList & bp() const
Definition: meshSurfaceEngineI.H:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:50
Foam::checkNonMappableCellConnections::~checkNonMappableCellConnections
~checkNonMappableCellConnections()
Definition: checkNonMappableCellConnections.C:168
Foam::checkNonMappableCellConnections::removeCells
bool removeCells()
find and remove problematic cells
Definition: checkNonMappableCellConnections.C:293
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
polyMeshGenModifier.H
Foam::Map< label >
Foam::polyMeshGen
Definition: polyMeshGen.H:46
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::HashSet< label, Hash< label > >
Foam::DynList::removeLastElement
T removeLastElement()
Return and remove the last element.
Definition: DynListI.H:384
Foam::polyMeshGenFaces::faces
const faceListPMG & faces() const
access to faces
Definition: polyMeshGenFacesI.H:43
Foam::polyMeshGenFaces::nInternalFaces
label nInternalFaces() const
return number of internal faces
Definition: polyMeshGenFacesI.H:48
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
checkNonMappableCellConnections.H
Foam::checkNonMappableCellConnections::INTERNALCELL
@ INTERNALCELL
Definition: checkNonMappableCellConnections.H:67
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
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::polyMeshGenFaces::faceIsInProcPatch
label faceIsInProcPatch(const label faceLabel) const
return processor patch label for the given face label
Definition: polyMeshGenFaces.C:168
Foam::Info
messageStream Info
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::checkNonMappableCellConnections::findCellTypes
void findCellTypes()
classify cells
Definition: checkNonMappableCellConnections.C:46
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
meshSurfaceEngine.H
Foam::checkNonMappableCellConnections::ALLBNDVERTEXCELL
@ ALLBNDVERTEXCELL
Definition: checkNonMappableCellConnections.H:69
Foam::DynList< label >
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::sumOp
Definition: ops.H:162
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:473
helperFunctions.H
Foam::polyMeshGenFaces::boundaries
const PtrList< boundaryPatch > & boundaries() const
ordinary boundaries
Definition: polyMeshGenFacesI.H:111
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::checkNonMappableCellConnections::cellType_
labelList cellType_
type of cell
Definition: checkNonMappableCellConnections.H:60
Foam::checkNonMappableCellConnections::findCells
void findCells(labelHashSet &badCells)
find problematic cells
Definition: checkNonMappableCellConnections.C:174
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::polyMeshGenFaces::procBoundaries
const PtrList< processorBoundaryPatch > & procBoundaries() const
inter-processor boundaries
Definition: polyMeshGenFacesI.H:106
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::DynList::size
label size() const
Definition: DynListI.H:235
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::checkNonMappableCellConnections::mesh_
polyMeshGen & mesh_
Reference to polyMeshGen.
Definition: checkNonMappableCellConnections.H:57
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::checkNonMappableCellConnections::BNDCELL
@ BNDCELL
Definition: checkNonMappableCellConnections.H:68
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::meshSurfaceEngine
Definition: meshSurfaceEngine.H:54
Foam::cellListPMG::size
label size() const
return the number of used elements
Definition: cellListPMGI.H:56
Foam::polyMeshGenModifier::removeCells
void removeCells(const boolList &removeCell, const bool removeProcFaces=true)
remove cells
Definition: polyMeshGenModifierRemoveCells.C:41
Foam::help::shareAnEdge
bool shareAnEdge(const faceType1 &f1, const faceType2 &f2)
check if two faces share an edge
Definition: helperFunctionsTopologyManipulationI.H:324
Foam::checkNonMappableCellConnections::checkNonMappableCellConnections
checkNonMappableCellConnections(const checkNonMappableCellConnections &)
Disallow default bitwise copy construct.
Foam::checkNonMappableCellConnections::INTERNALFACEGROUP
@ INTERNALFACEGROUP
Definition: checkNonMappableCellConnections.H:70
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304