boundaryLayersFacesAndCells.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 "boundaryLayers.H"
29 #include "meshSurfaceEngine.H"
30 #include "helperFunctions.H"
31 #include "demandDrivenData.H"
32 #include "VRWGraphList.H"
33 
34 #include <map>
35 
36 //#define DEBUGLayer
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
46 {
47  Info << "Starting creating layer cells" << endl;
48 
49  const meshSurfaceEngine& mse = surfaceEngine();
50  const faceList::subList& bFaces = mse.boundaryFaces();
51  const VRWGraph& faceEdges = mse.faceEdges();
52  const VRWGraph& edgeFaces = mse.edgeFaces();
53  const labelList& boundaryFacePatches = mse.boundaryFacePatches();
54  const labelList& faceOwners = mse.faceOwners();
55 
56  //- this is used for parallel runs
57  const Map<label>* otherProcPatchPtr(NULL);
58 
59  if( Pstream::parRun() )
60  {
61  createNewFacesParallel(treatPatches);
62 
63  otherProcPatchPtr = &mse.otherEdgeFacePatch();
64  }
65 
66  //- create lists for new boundary faces
67  VRWGraph newBoundaryFaces;
68  labelLongList newBoundaryOwners;
69  labelLongList newBoundaryPatches;
70 
71  //- create storage for new cells
72  VRWGraphList cellsToAdd;
73 
74  //- create layer cells and store boundary faces
75  const label nOldCells = mesh_.cells().size();
76  forAll(bFaces, bfI)
77  if( treatPatches[boundaryFacePatches[bfI]] )
78  {
79  const face& f = bFaces[bfI];
80 
81  faceList cellFaces(f.size() + 2);
82 
83  label fI(0);
84 
85  //- store boundary face
86  cellFaces[fI++] = f.reverseFace();
87 
88  //- create parallel face
89  face newF(f.size());
90  forAll(f, pI)
91  newF[pI] = newLabelForVertex_[f[pI]];
92  cellFaces[fI++] = newF;
93 
94  newBoundaryFaces.appendList(newF);
95  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
96  newBoundaryPatches.append(boundaryFacePatches[bfI]);
97 
98  //- create quad faces
99  newF.setSize(4);
100  forAll(f, pI)
101  {
102  newF[0] = f[pI];
103  newF[1] = f.nextLabel(pI);
104  newF[2] = newLabelForVertex_[f.nextLabel(pI)];
105  newF[3] = newLabelForVertex_[f[pI]];
106 
107  cellFaces[fI++] = newF;
108 
109  //- check if the face is at the boundary
110  //- of the treated partitions
111  const label edgeI = faceEdges(bfI, pI);
112  if( edgeFaces.sizeOfRow(edgeI) == 2 )
113  {
114  label neiFace = edgeFaces(edgeI, 0);
115  if( neiFace == bfI )
116  neiFace = edgeFaces(edgeI, 1);
117 
118  if( !treatPatches[boundaryFacePatches[neiFace]] )
119  {
120  newBoundaryFaces.appendList(newF);
121  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
122  newBoundaryPatches.append(boundaryFacePatches[neiFace]);
123  }
124  }
125  else if( edgeFaces.sizeOfRow(edgeI) == 1 )
126  {
127  const Map<label>& otherProcPatch = *otherProcPatchPtr;
128  if( !treatPatches[otherProcPatch[edgeI]] )
129  {
130  newBoundaryFaces.appendList(newF);
131  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
132  newBoundaryPatches.append(otherProcPatch[edgeI]);
133  }
134  }
135  }
136 
137  cellsToAdd.appendGraph(cellFaces);
138  }
139  else
140  {
141  # ifdef DEBUGLayer
142  Info << "Storing original boundary face "
143  << bfI << " into patch " << boundaryFacePatches[bfI] << endl;
144  # endif
145 
146  newBoundaryFaces.appendList(bFaces[bfI]);
147  newBoundaryOwners.append(faceOwners[bfI]);
148  newBoundaryPatches.append(boundaryFacePatches[bfI]);
149  }
150 
151  //- create mesh modifier
152  polyMeshGenModifier meshModifier(mesh_);
153 
154  meshModifier.addCells(cellsToAdd);
155  cellsToAdd.clear();
156  meshModifier.reorderBoundaryFaces();
157  meshModifier.replaceBoundary
158  (
159  patchNames_,
160  newBoundaryFaces,
161  newBoundaryOwners,
162  newBoundaryPatches
163  );
164 
165  PtrList<boundaryPatch>& boundaries = meshModifier.boundariesAccess();
166  forAll(boundaries, patchI)
167  boundaries[patchI].patchType() = patchTypes_[patchI];
168 
169  //- delete meshSurfaceEngine
170  this->clearOut();
171 
172  # ifdef DEBUGLayer
173  mesh_.addressingData().checkMesh(true);
174  # endif
175 
176  Info << "Finished creating layer cells" << endl;
177 }
178 
180 (
181  const boolList& treatPatches
182 )
183 {
184  const meshSurfaceEngine& mse = surfaceEngine();
185  const faceList::subList& bFaces = mse.boundaryFaces();
186  const VRWGraph& faceEdges = mse.faceEdges();
187  const VRWGraph& edgeFaces = mse.edgeFaces();
188  const labelList& boundaryFacePatches = mse.boundaryFacePatches();
189  const labelList& globalEdgeLabel = mse.globalBoundaryEdgeLabel();
190  const Map<label>& globalToLocal = mse.globalToLocalBndEdgeAddressing();
191 
192  const Map<label>& otherProcPatch = mse.otherEdgeFacePatch();
193  const Map<label>& otherFaceProc = mse.otherEdgeFaceAtProc();
194 
195  //- the next stage is the generation of processor faces
196  //- this step can be done without any communication, but only if the faces
197  //- are added in the same order on both processors
198  //- this will be achieved by sorting edges according to their global labes
199  //- another difficulty here is that new processor patches may occur
200  //- during this procedure
201  Map<label> otherProcToProcPatch;
202  forAll(mesh_.procBoundaries(), patchI)
203  {
204  const processorBoundaryPatch& wp = mesh_.procBoundaries()[patchI];
205  otherProcToProcPatch.insert(wp.neiProcNo(), patchI);
206  }
207 
208  label nTreatedEdges(0);
209  boolList treatEdge(edgeFaces.size(), false);
210  for
211  (
212  Map<label>::const_iterator iter=globalToLocal.begin();
213  iter!=globalToLocal.end();
214  ++iter
215  )
216  {
217  const label beI = iter();
218 
219  if( edgeFaces.sizeOfRow(beI) != 1 )
220  continue;
221 
222  if(
223  treatPatches[boundaryFacePatches[edgeFaces(beI, 0)]] &&
224  treatPatches[otherProcPatch[beI]]
225  )
226  {
227  ++nTreatedEdges;
228  treatEdge[beI] = true;
229  }
230  }
231 
232  //- create a list of treated edges and sort the list
233  labelList treatedEdgeLabels(nTreatedEdges);
234  nTreatedEdges = 0;
235  forAll(treatEdge, beI)
236  if( treatEdge[beI] )
237  {
238  treatedEdgeLabels[nTreatedEdges++] = globalEdgeLabel[beI];
239  }
240  treatedEdgeLabels.setSize(nTreatedEdges);
241 
242  sort(treatedEdgeLabels);
243 
244  //- create additional processor patches if needed
245  forAll(treatedEdgeLabels, eI)
246  {
247  const label beI = globalToLocal[treatedEdgeLabels[eI]];
248 
249  if( !otherProcToProcPatch.found(otherFaceProc[beI]) )
250  {
251  otherProcToProcPatch.insert
252  (
253  otherFaceProc[beI],
254  polyMeshGenModifier(mesh_).addProcessorPatch
255  (
256  otherFaceProc[beI]
257  )
258  );
259  }
260  }
261 
262  //- create new processor faces
263  VRWGraph newProcFaces;
264  labelLongList faceProcPatch;
265  FixedList<label, 4> newF;
266  forAll(treatedEdgeLabels, geI)
267  {
268  const label beI = globalToLocal[treatedEdgeLabels[geI]];
269 
270  if( edgeFaces.sizeOfRow(beI) == 0 )
271  continue;
272 
273  const label bfI = edgeFaces(beI, 0);
274  const label pos = faceEdges.containsAtPosition(bfI, beI);
275  const edge e = bFaces[bfI].faceEdge(pos);
276 
277  if( otherFaceProc[beI] > Pstream::myProcNo() )
278  {
279  newF[0] = e.start();
280  newF[1] = e.end();
281  if( patchKey_.size() != 0 )
282  {
283  newF[2] =
284  findNewNodeLabel(e.end(), patchKey_[otherProcPatch[beI]]);
285  newF[3] =
286  findNewNodeLabel(e.start(), patchKey_[otherProcPatch[beI]]);
287  }
288  else
289  {
290  newF[2] = newLabelForVertex_[e.end()];
291  newF[3] = newLabelForVertex_[e.start()];
292  }
293  }
294  else
295  {
296  newF[0] = e.end();
297  if( patchKey_.size() != 0 )
298  {
299  newF[1] =
300  findNewNodeLabel
301  (
302  e.end(),
303  patchKey_[boundaryFacePatches[bfI]]
304  );
305  newF[2] =
306  findNewNodeLabel
307  (
308  e.start(),
309  patchKey_[boundaryFacePatches[bfI]]
310  );
311  }
312  else
313  {
314  newF[1] = newLabelForVertex_[e.end()];
315  newF[2] = newLabelForVertex_[e.start()];
316  }
317  newF[3] = e.start();
318  }
319 
320  newProcFaces.appendList(newF);
321  faceProcPatch.append(otherProcToProcPatch[otherFaceProc[beI]]);
322  }
323 
324  //- add faces into the mesh
325  polyMeshGenModifier(mesh_).addProcessorFaces(newProcFaces, faceProcPatch);
326 }
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 } // End namespace Foam
331 
332 // ************************************************************************* //
Foam::polyMeshGenCells::addressingData
const polyMeshGenAddressing & addressingData() const
addressing which may be needed
Definition: polyMeshGenCells.C:327
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::boundaryLayers::patchNames_
wordList patchNames_
patch names
Definition: boundaryLayers.H:82
Foam::boundaryLayers::surfaceEngine
const meshSurfaceEngine & surfaceEngine() const
Return const reference to meshSurfaceEngine.
Definition: boundaryLayers.C:54
Foam::meshSurfaceEngine::globalToLocalBndEdgeAddressing
const Map< label > & globalToLocalBndEdgeAddressing() const
global boundary edge label to local label. Only for processor edges
Definition: meshSurfaceEngineI.H:510
Foam::meshSurfaceEngine::globalBoundaryEdgeLabel
const labelList & globalBoundaryEdgeLabel() const
global boundary edge label
Definition: meshSurfaceEngineI.H:489
Foam::VRWGraphList::size
label size() const
Returns the number of graphs.
Definition: VRWGraphListI.H:96
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::boundaryLayers::mesh_
polyMeshGen & mesh_
Reference to the mesh.
Definition: boundaryLayers.H:63
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::boundaryLayers::patchTypes_
wordList patchTypes_
patch types
Definition: boundaryLayers.H:85
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::polyMeshGenModifier::addProcessorFaces
void addProcessorFaces(const VRWGraph &procFaces, const labelLongList &facePatches)
add additional faces into processor patches
Definition: polyMeshGenModifierAddProcessorFaces.C:41
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::meshSurfaceEngine::boundaryFacePatches
const labelList & boundaryFacePatches() const
patch label for each boundary face
Definition: meshSurfaceEngineI.H:123
Foam::meshSurfaceEngine::faceEdges
const VRWGraph & faceEdges() const
Definition: meshSurfaceEngineI.H:353
Foam::Map< label >
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::boundaryLayers::newLabelForVertex_
labelLongList newLabelForVertex_
label of a new node (helper)
Definition: boundaryLayers.H:95
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::VRWGraphList
Definition: VRWGraphList.H:51
Foam::VRWGraphList::appendGraph
void appendGraph(const GraphType &l)
Append a graph at the end of the graphList.
Definition: VRWGraphListI.H:123
Foam::polyMeshGenModifier::boundariesAccess
PtrList< boundaryPatch > & boundariesAccess()
access to boundary data
Definition: polyMeshGenModifier.H:131
Foam::LongList< label >
Foam::meshSurfaceEngine::boundaryFaces
const faceList::subList & boundaryFaces() const
Definition: meshSurfaceEngineI.H:103
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::Info
messageStream Info
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::meshSurfaceEngine::otherEdgeFaceAtProc
const Map< label > & otherEdgeFaceAtProc() const
Definition: meshSurfaceEngineI.H:568
Foam::VRWGraph::size
label size() const
Returns the number of rows.
Definition: VRWGraphI.H:122
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::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::meshSurfaceEngine::edgeFaces
const VRWGraph & edgeFaces() const
Definition: meshSurfaceEngineI.H:334
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
meshSurfaceEngine.H
Foam::processorBoundaryPatch
Definition: processorBoundaryPatch.H:47
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::boundaryLayers::clearOut
void clearOut()
delete meshSurfaceEngine
Definition: boundaryLayers.H:225
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::VRWGraph::containsAtPosition
label containsAtPosition(const label rowI, const label e) const
Definition: VRWGraphI.H:529
Foam::processorBoundaryPatch::neiProcNo
label neiProcNo() const
return the neighbour processor
Definition: processorBoundaryPatch.H:118
Foam::VRWGraphList::clear
void clear()
Clear the graph.
Definition: VRWGraphListI.H:115
helperFunctions.H
f
labelList f(nPoints)
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::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::FixedList< label, 4 >
Foam::boundaryLayers::createNewFacesParallel
void createNewFacesParallel(const boolList &treatPatches)
Definition: boundaryLayersFacesAndCells.C:180
Foam::polyMeshGenModifier::addCells
void addCells(const LongList< faceList > &cellFaces)
add cells (vertices must be added)
Definition: polyMeshGenModifierAddCells.C:37
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::sort
void sort(UList< T > &)
Definition: UList.C:107
Foam::boundaryLayers::createNewFacesAndCells
void createNewFacesAndCells(const boolList &treatPatches)
create a layer of cells
Definition: boundaryLayersFacesAndCells.C:45
Foam::meshSurfaceEngine::otherEdgeFacePatch
const Map< label > & otherEdgeFacePatch() const
Definition: meshSurfaceEngineI.H:588
Foam::VRWGraph
Definition: VRWGraph.H:101
VRWGraphList.H
boundaryLayers.H
Foam::meshSurfaceEngine
Definition: meshSurfaceEngine.H:54
Foam::cellListPMG::size
label size() const
return the number of used elements
Definition: cellListPMGI.H:56
Foam::meshSurfaceEngine::faceOwners
const labelList & faceOwners() const
Definition: meshSurfaceEngineI.H:143
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190