cartesianMeshExtractorPolyMesh.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 
28 #include "cartesianMeshExtractor.H"
29 #include "demandDrivenData.H"
30 #include "meshOctree.H"
31 #include "labelledPair.H"
32 
33 #include "meshSurfaceEngine.H"
34 #include "helperFunctions.H"
35 #include "helperFunctionsPar.H"
36 #include "decomposeFaces.H"
37 #include "decomposeCells.H"
38 
39 #include <map>
40 #include <sstream>
41 
42 //#define DEBUGMesh
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
52 {
53  Info << "Creating polyMesh from octree" << endl;
54 
55  const meshOctree& octree = octreeCheck_.octree();
56 
57  //- give labels to cubes which will be used as mesh cells
58  const List<direction>& cType = octreeCheck_.boxType();
59 
60  labelList& leafCellLabel = *leafCellLabelPtr_;
61  label nCells(0);
62  forAll(cType, leafI)
63  {
64  if(
65  Pstream::parRun() &&
66  (octree.returnLeaf(leafI).procNo() != Pstream::myProcNo())
67  )
68  continue;
69 
70  if( cType[leafI] & meshOctreeAddressing::MESHCELL )
71  {
72  leafCellLabel[leafI] = nCells++;
73  }
74  }
75 
76  //- access to mesh data
77  polyMeshGenModifier meshModifier(mesh_);
78  faceListPMG& faces = meshModifier.facesAccess();
79  cellListPMG& cells = meshModifier.cellsAccess();
80 
81  //- start creating octree mesh
82  cells.setSize(nCells);
83  List<direction> nFacesInCell(nCells, direction(0));
84  label nFaces(0);
85 
86  const VRWGraph& octreeFaces = octreeCheck_.octreeFaces();
88  const labelLongList& neighbour = octreeCheck_.octreeFaceNeighbour();
89 
90  //- map storing box label and a direction for each processor face
91  //- The map stores data in the same order on both sides of processor
92  //- boundaries. This is a consequence of Morton ordering of
93  //- leaf boxes in the octree.
94  std::map<label, labelLongList> procFaces;
95 
96  forAll(octreeFaces, faceI)
97  {
98  const label own = owner[faceI];
99  const label nei = neighbour[faceI];
100 
101  const label ownLabel = leafCellLabel[own];
102  label neiLabel(-1);
103  if( nei != -1 )
104  neiLabel = leafCellLabel[nei];
105 
106  if( (ownLabel != -1) && (neiLabel != -1) )
107  {
108  ++nFaces;
109  ++nFacesInCell[ownLabel];
110  ++nFacesInCell[neiLabel];
111  }
112  else if( ownLabel != -1 )
113  {
114  ++nFaces;
115  ++nFacesInCell[ownLabel];
116 
117  if( (nei != -1) && (cType[nei] & meshOctreeAddressing::MESHCELL) )
118  {
119  const label procNo = octree.returnLeaf(nei).procNo();
120 
121  procFaces[procNo].append(faceI);
122  }
123  }
124  else if( neiLabel != -1 )
125  {
126  ++nFaces;
127  ++nFacesInCell[neiLabel];
128 
129  if( (own != -1) && (cType[own] & meshOctreeAddressing::MESHCELL) )
130  {
131  const label procNo = octree.returnLeaf(own).procNo();
132 
133  procFaces[procNo].append(faceI);
134  }
135  }
136  }
137 
138  //- case is a serial run
139  faces.setSize(nFaces);
140  forAll(cells, cI)
141  cells[cI].setSize(nFacesInCell[cI]);
142  nFacesInCell = 0;
143 
144  //- calculate faces in processor patches
145  if( Pstream::parRun() )
146  {
147  PtrList<processorBoundaryPatch>& procBoundaries =
148  meshModifier.procBoundariesAccess();
149 
150  //- set the number of procBoundaries
151  procBoundaries.setSize(procFaces.size());
152  std::ostringstream ss;
153  ss << Pstream::myProcNo();
154  const word name("processor"+ss.str()+"to");
155  label nProcBoundaries(nFaces), patchI(0);
156 
157  //- allocate memory for processor patches
158  std::map<label, labelLongList>::const_iterator iter;
159  for(iter=procFaces.begin();iter!=procFaces.end();++iter)
160  {
161  const label procI = iter->first;
162 
163  std::ostringstream ssNei;
164  ssNei << procI;
165  procBoundaries.set
166  (
167  patchI,
169  (
170  name+ssNei.str(),
171  "processor",
172  iter->second.size(),
173  0,
175  procI
176  )
177  );
178 
179  nProcBoundaries -= iter->second.size();
180  ++patchI;
181  }
182 
183  //- create processor faces
184  //- they need to be created here because of the correct ordering
185  patchI = 0;
186  for(iter=procFaces.begin();iter!=procFaces.end();++iter)
187  {
188  procBoundaries[patchI].patchStart() = nProcBoundaries;
189 
190  const labelLongList& patchFaces = iter->second;
191 
192  forAll(patchFaces, pfI)
193  {
194  const label fLabel = patchFaces[pfI];
195  const label own = owner[fLabel];
196  const label nei = neighbour[fLabel];
197 
198  const label curCell = leafCellLabel[own];
199  label neiCell(-1);
200  if( nei != -1 )
201  neiCell = leafCellLabel[nei];
202 
203  //- create a processor face
204  if( neiCell == -1 )
205  {
206  //- add a face
207  faces[nProcBoundaries].setSize(octreeFaces[fLabel].size());
208  forAllRow(octreeFaces, fLabel, pI)
209  faces[nProcBoundaries][pI] = octreeFaces(fLabel, pI);
210  cells[curCell][nFacesInCell[curCell]++] = nProcBoundaries++;
211  }
212  else if( curCell == -1 )
213  {
214  //- add a reversed face
215  faces[nProcBoundaries].setSize(octreeFaces[fLabel].size());
216  label i(0);
217  faces[nProcBoundaries][i++] = octreeFaces(fLabel, 0);
218  for(label pI=octreeFaces.sizeOfRow(fLabel)-1;pI>0;--pI)
219  faces[nProcBoundaries][i++] = octreeFaces(fLabel, pI);
220  cells[neiCell][nFacesInCell[neiCell]++] = nProcBoundaries++;
221  }
222  else
223  {
225  (
226  "void cartesianMeshExtractor::createPolyMesh()"
227  ) << "Face " << octreeFaces[fLabel] << " causes problems!"
228  << abort(FatalError);
229  }
230  }
231 
232  if( procBoundaries[patchI].patchSize() !=
233  (nProcBoundaries - procBoundaries[patchI].patchStart())
234  )
236  (
237  "cartesianMeshExtractor::createPolyMesh()"
238  ) << "Invalid patch size!" << Pstream::myProcNo()
239  << abort(FatalError);
240 
241  ++patchI;
242  }
243  }
244 
245  nFaces = 0;
246 
247  forAll(octreeFaces, faceI)
248  {
249  const label own = owner[faceI];
250  const label nei = neighbour[faceI];
251 
252  const label ownLabel = leafCellLabel[own];
253  label neiLabel(-1);
254  if( nei != -1 )
255  neiLabel = leafCellLabel[nei];
256 
257  if( (ownLabel != -1) && (neiLabel != -1) )
258  {
259  //- internal face
260  faces[nFaces].setSize(octreeFaces.sizeOfRow(faceI));
261  forAllRow(octreeFaces, faceI, pI)
262  faces[nFaces][pI] = octreeFaces(faceI, pI);
263 
264  cells[ownLabel][nFacesInCell[ownLabel]++] = nFaces;
265  cells[neiLabel][nFacesInCell[neiLabel]++] = nFaces;
266  ++nFaces;
267  }
268  else if( ownLabel != -1 )
269  {
270  if( (nei != -1) && (cType[nei] & meshOctreeAddressing::MESHCELL) )
271  {
272  //- face at a parallel boundary
273  continue;
274  }
275 
276  //- boundary face
277  faces[nFaces].setSize(octreeFaces.sizeOfRow(faceI));
278  forAllRow(octreeFaces, faceI, pI)
279  faces[nFaces][pI] = octreeFaces(faceI, pI);
280 
281  cells[ownLabel][nFacesInCell[ownLabel]++] = nFaces;
282  ++nFaces;
283  }
284  else if( neiLabel != -1 )
285  {
286  if( (own != -1) && (cType[own] & meshOctreeAddressing::MESHCELL) )
287  {
288  //- face at a parallel boundary
289  continue;
290  }
291 
292  //- boundary face
293  faces[nFaces].setSize(octreeFaces.sizeOfRow(faceI));
294  faces[nFaces][0] = octreeFaces(faceI, 0);
295  for(label pI=octreeFaces.sizeOfRow(faceI)-1;pI>0;--pI)
296  faces[nFaces][octreeFaces.sizeOfRow(faceI)-pI] =
297  octreeFaces(faceI, pI);
298 
299  cells[neiLabel][nFacesInCell[neiLabel]++] = nFaces;
300  ++nFaces;
301  }
302  }
303 
304  # ifdef DEBUGMesh
305  label nProcBoundaries(0);
306  forAll(procBoundaries, patchI)
307  nProcBoundaries += procBoundaries[patchI].patchSize();
308 
309  if( faces.size() != (nProcBoundaries + nFaces) )
310  {
311  Serr << "Number of faces " << faces.size() << endl;
312  Serr << "Number of processor boundaries " << nProcBoundaries << endl;
313  Serr << "Number of domain faces " << nFaces << endl;
315  (
316  "void cartesianMeshExtractor::createPolyMesh()"
317  ) << Pstream::myProcNo() << "This mesh is invalid!"
318  << abort(FatalError);
319  }
320 
321  vectorField closedness(cells.size(), vector::zero);
322  const labelList& owner = mesh_.owner();
323  const labelList& neighbour = mesh_.neighbour();
324  forAll(owner, faceI)
325  if( owner[faceI] == -1 )
326  {
327  Info << "faces " << faces << endl;
329  (
330  "void cartesianMeshExtractor::createPolyMesh"
331  "("
332  "pointFieldPMG& points,"
333  "faceListPMG& faces,"
334  "cellListPMG& cells"
335  ")"
336  ) << "Face " << faceI
337  << " has no owner and neighbour!!" << abort(FatalError);
338  }
339 
340  forAll(faces, faceI)
341  {
342  const vector area = faces[faceI].normal(mesh_.points());
343  closedness[owner[faceI]] += area;
344  if( neighbour[faceI] != -1 )
345  closedness[neighbour[faceI]] -= area;
346  }
347 
348  forAll(closedness, cellI)
349  if( mag(closedness[cellI]) > 1e-10 )
350  Info << "Cell " << cellI << " is not closed by "
351  << closedness[cellI] << endl;
352 
353  # endif
354 
355  meshModifier.reorderBoundaryFaces();
356 
357  if( octree.isQuadtree() )
358  {
359  //- generate empty patches
360  //- search for faces with a dominant z coordinate and store them
361  //- into an empty patch
363  const vectorField& fNormals = mse.faceNormals();
364  const faceList::subList& bFaces = mse.boundaryFaces();
365  const labelList& fOwner = mse.faceOwners();
366  const vectorField& fCentres = mse.faceCentres();
367 
368  const boundBox& bb = octree.rootBox();
369  const scalar tZ = 0.05 * (bb.max().z() - bb.min().z());
370 
371  wordList patchNames(3);
372  patchNames[0] = "defaultFaces";
373  patchNames[1] = "unusedFacesBottom";
374  patchNames[2] = "unusedFacesTop";
375 
376  VRWGraph boundaryFaces;
377  labelLongList newFaceOwner;
378  labelLongList newFacePatch;
379 
380  forAll(fNormals, bfI)
381  {
382  //- store the face and its owner
383  boundaryFaces.appendList(bFaces[bfI]);
384  newFaceOwner.append(fOwner[bfI]);
385 
386  const vector& fNormal = fNormals[bfI];
387 
388  if( Foam::mag(fNormal.z()) > Foam::mag(fNormal.x() + fNormal.y()) )
389  {
390  if( Foam::mag(fCentres[bfI].z() - bb.min().z()) < tZ )
391  {
392  newFacePatch.append(1);
393  }
394  else if( Foam::mag(fCentres[bfI].z() - bb.max().z()) < tZ )
395  {
396  newFacePatch.append(2);
397  }
398  else
399  {
401  (
402  "void cartesianMeshExtractor::createPolyMesh()"
403  ) << "Cannot distribute the face!!" << exit(FatalError);
404  }
405  }
406  else
407  {
408  newFacePatch.append(0);
409  }
410  }
411 
412  //- replace the boundary with faces in correct patches
413  meshModifier.replaceBoundary
414  (
415  patchNames,
416  boundaryFaces,
417  newFaceOwner,
418  newFacePatch
419  );
420 
421  meshModifier.boundariesAccess()[1].patchType() = "empty";
422  meshModifier.boundariesAccess()[2].patchType() = "empty";
423  }
424 
425  Info << "Finished creating polyMesh" << endl;
426 }
427 
428 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
429 
430 } // End namespace Foam
431 
432 // ************************************************************************* //
Foam::meshSurfaceEngine::faceCentres
const vectorField & faceCentres() const
Definition: meshSurfaceEngineI.H:277
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::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
helperFunctionsPar.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::cartesianMeshExtractor::octreeCheck_
meshOctreeAddressing octreeCheck_
reference to the octree addressing
Definition: cartesianMeshExtractor.H:55
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
decomposeCells.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::meshOctree::isQuadtree
bool isQuadtree() const
is octree a quadtree or an octree
Definition: meshOctreeI.H:36
decomposeFaces.H
Foam::polyMeshGenModifier::replaceBoundary
void replaceBoundary(const wordList &patchNames, const VRWGraph &boundaryFaces, const labelLongList &faceOwners, const labelLongList &facePatches)
replace the boundary with new boundary faces
Definition: polyMeshGenModifierReplaceBoundary.C:42
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMeshGenPoints::points
const pointFieldPMG & points() const
access to points
Definition: polyMeshGenPointsI.H:44
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::cellListPMG
Definition: cellListPMG.H:49
patchFaces
labelList patchFaces(const polyBoundaryMesh &patches, const wordList &names)
Definition: extrudeMesh.C:148
cartesianMeshExtractor.H
Foam::meshOctreeAddressing::octreeFaceOwner
const labelLongList & octreeFaceOwner() const
return owners of octree faces
Definition: meshOctreeAddressingI.H:110
meshOctree.H
Foam::cartesianMeshExtractor::createPolyMesh
void createPolyMesh()
create mesh data
Definition: cartesianMeshExtractorPolyMesh.C:51
Foam::polyMeshGenModifier::boundariesAccess
PtrList< boundaryPatch > & boundariesAccess()
access to boundary data
Definition: polyMeshGenModifier.H:131
Foam::LongList< label >
Foam::PtrList::set
bool set(const label) const
Is element set.
labelledPair.H
Foam::polyMeshGenModifier::procBoundariesAccess
PtrList< processorBoundaryPatch > & procBoundariesAccess()
access to processor boundary data
Definition: polyMeshGenModifier.H:125
Foam::meshSurfaceEngine::boundaryFaces
const faceList::subList & boundaryFaces() const
Definition: meshSurfaceEngineI.H:103
Foam::polyMeshGenModifier::cellsAccess
cellListPMG & cellsAccess()
access to cells
Definition: polyMeshGenModifier.H:119
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::polyMeshGenModifier::reorderBoundaryFaces
void reorderBoundaryFaces()
Definition: polyMeshGenModifierReorderBoundaryFaces.C:44
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::Serr
OSstream Serr(cerr, "Serr")
Definition: IOstreams.H:52
Foam::faceListPMG::setSize
void setSize(const label nElmts)
set the number of used elements
Definition: faceListPMGI.H:78
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
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
patchNames
wordList patchNames(nPatches)
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::cartesianMeshExtractor::leafCellLabelPtr_
labelList * leafCellLabelPtr_
cell label for a given leaf
Definition: cartesianMeshExtractor.H:64
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
Foam::meshOctreeAddressing::MESHCELL
@ MESHCELL
Definition: meshOctreeAddressing.H:234
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
meshSurfaceEngine.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::processorBoundaryPatch
Definition: processorBoundaryPatch.H:47
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::meshOctreeAddressing::octreeFaceNeighbour
const labelLongList & octreeFaceNeighbour() const
return neighbours of octree faces
Definition: meshOctreeAddressingI.H:118
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::meshOctreeAddressing::boxType
const List< direction > & boxType() const
return which octree boxes are used for mesh creation
Definition: meshOctreeAddressingI.H:68
Foam::meshOctreeAddressing::octree
const meshOctree & octree() const
return const reference to meshOctree
Definition: meshOctreeAddressingI.H:89
Foam::meshOctree::rootBox
const boundBox & rootBox() const
return rootBox
Definition: meshOctreeI.H:135
helperFunctions.H
Foam::faceListPMG::size
label size() const
return the number of used elements
Definition: faceListPMGI.H:73
Foam::VRWGraph::appendList
void appendList(const ListType &l)
Append a list as a row at the end of the graph.
Definition: VRWGraphI.H:286
Foam::Vector< scalar >
Foam::List< direction >
Foam::meshOctree
Definition: meshOctree.H:55
Foam::PtrList::setSize
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
Foam::meshOctreeCubeBasic::procNo
short procNo() const
return processor number
Definition: meshOctreeCubeBasicI.H:80
Foam::cartesianMeshExtractor::mesh_
polyMeshGen & mesh_
reference to the mesh
Definition: cartesianMeshExtractor.H:58
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::polyMeshGenModifier::facesAccess
faceListPMG & facesAccess()
access to mesh faces
Definition: polyMeshGenModifier.H:113
Foam::meshOctreeAddressing::octreeFaces
const VRWGraph & octreeFaces() const
return octree faces, created for MESHCELL boxes
Definition: meshOctreeAddressingI.H:102
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::meshSurfaceEngine::faceNormals
const vectorField & faceNormals() const
Definition: meshSurfaceEngineI.H:258
Foam::meshOctree::returnLeaf
const meshOctreeCubeBasic & returnLeaf(const label) const
Definition: meshOctreeI.H:60
Foam::meshSurfaceEngine
Definition: meshSurfaceEngine.H:54
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::meshSurfaceEngine::faceOwners
const labelList & faceOwners() const
Definition: meshSurfaceEngineI.H:143