polyMeshGenModifierAddBufferCells.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 "polyMeshGenModifier.H"
29 #include "polyMeshGenAddressing.H"
30 #include "demandDrivenData.H"
31 #include "labelledPoint.H"
32 #include "helperFunctions.H"
33 
34 // #define DEBUGSearch
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
44 {
45  if( !Pstream::parRun() )
46  return;
47 
48  Info << "Adding buffer layers" << endl;
49 
50  const labelList& owner = mesh_.owner();
52  faceListPMG& faces = facesAccess();
53  const cellListPMG& cells = mesh_.cells();
54  const PtrList<processorBoundaryPatch>& procBoundaries = mesh_.procBoundaries();
55  const polyMeshGenAddressing& addressing = mesh_.addressingData();
56  const labelLongList& globalPointLabel = addressing.globalPointLabel();
57  const Map<label>& globalToLocal = addressing.globalToLocalPointAddressing();
58 
59  //- receive vertices
60  forAll(procBoundaries, patchI)
61  {
62  labelHashSet pointsToSend;
63 
64  label faceI = procBoundaries[patchI].patchStart();
65  const label end = faceI + procBoundaries[patchI].patchSize();
66  for(;faceI<end;++faceI)
67  {
68  const cell& c = cells[owner[faceI]];
69  forAll(c, fI)
70  {
71  const face& f = faces[c[fI]];
72 
73  forAll(f, pI)
74  pointsToSend.insert(f[pI]);
75  }
76  }
77 
78  faceI = 0;
79  List<labelledPoint> ptsToSend(pointsToSend.size());
80  forAllConstIter(labelHashSet, pointsToSend, it)
81  ptsToSend[faceI++] =
83  (
84  globalPointLabel[it.key()],
85  points[it.key()]
86  );
87 
88  OPstream toOtherProc
89  (
91  procBoundaries[patchI].neiProcNo()
92  );
93 
94  toOtherProc << ptsToSend;
95  }
96 
97  Map<label> globalToLocalReceived;
98  forAll(procBoundaries, patchI)
99  {
100  List<labelledPoint> receivedPoints;
101  IPstream fromOtherProc
102  (
104  procBoundaries[patchI].neiProcNo()
105  );
106 
107  fromOtherProc >> receivedPoints;
108 
109  forAll(receivedPoints, i)
110  {
111  if( globalToLocal.found(receivedPoints[i].pointLabel()) )
112  continue;
113  if( globalToLocalReceived.found(receivedPoints[i].pointLabel()) )
114  continue;
115 
116  globalToLocalReceived.insert
117  (
118  receivedPoints[i].pointLabel(),
119  points.size()
120  );
121  points.append(receivedPoints[i].coordinates());
122  }
123  }
124 
125  //- send cells to other processors
126  forAll(procBoundaries, patchI)
127  {
128  labelHashSet cellsToSend;
129 
130  label faceI = procBoundaries[patchI].patchStart();
131  const label end = faceI + procBoundaries[patchI].patchSize();
132  for(;faceI<end;++faceI)
133  cellsToSend.insert(owner[faceI]);
134 
135  labelLongList flattenedCells;
136  forAllConstIter(labelHashSet, cellsToSend, it)
137  {
138  const cell& c = cells[it.key()];
139 
140  //- the number of faces in the cell
141  flattenedCells.append(c.size());
142 
143  //- add faces
144  forAll(c, fI)
145  {
146  const face& f = faces[c[fI]];
147 
148  //- the number of vertices in the face
149  flattenedCells.append(f.size());
150  forAll(f, pI)
151  flattenedCells.append(globalPointLabel[f[pI]]);
152  }
153  }
154 
155  OPstream toOtherProc
156  (
158  procBoundaries[patchI].neiProcNo()
159  );
160 
161  toOtherProc << flattenedCells;
162  }
163 
164  forAll(procBoundaries, patchI)
165  {
166  word subsetName = "processor_";
167  subsetName += help::scalarToText(procBoundaries[patchI].neiProcNo());
168  const label subsetID = mesh_.addCellSubset(subsetName);
169 
170  labelList receivedCells;
171 
172  IPstream fromOtherProc
173  (
175  procBoundaries[patchI].neiProcNo()
176  );
177 
178  fromOtherProc >> receivedCells;
179 
180  label counter(0);
181  while( counter < receivedCells.size() )
182  {
183  faceList c(receivedCells[counter++]);
184  forAll(c, fI)
185  {
186  c[fI].setSize(receivedCells[counter++]);
187  forAll(c[fI], pI)
188  {
189  const label gpI = receivedCells[counter++];
190 
191  if( globalToLocal.found(gpI) )
192  {
193  c[fI][pI] = globalToLocal[gpI];
194  }
195  else
196  {
197  c[fI][pI] = globalToLocalReceived[gpI];
198  }
199  }
200  }
201 
202  mesh_.addCellToSubset(subsetID, cells.size());
203  addCell(c);
204  }
205  }
206 
207  mesh_.clearOut();
208 
209  Info << "Finished adding buffer layers" << endl;
210 }
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace Foam
215 
216 // ************************************************************************* //
Foam::polyMeshGenFaces::owner
const labelList & owner() const
owner and neighbour cells for faces
Definition: polyMeshGenFacesI.H:67
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::word
A class for handling words, derived from string.
Definition: word.H:59
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::polyMeshGenAddressing::globalToLocalPointAddressing
const Map< label > & globalToLocalPointAddressing() const
Definition: polyMeshGenAddressingParallelAddressing.C:814
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::help::scalarToText
word scalarToText(const scalar s)
convert the scalar value into text
Definition: helperFunctionsStringConversion.C:61
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
polyMeshGenModifier.H
Foam::Map< label >
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::polyMeshGenCells::clearOut
void clearOut() const
clear all pointer data
Definition: polyMeshGenCells.C:254
Foam::polyMeshGenModifier::mesh_
polyMeshGen & mesh_
reference to the mesh
Definition: polyMeshGenModifier.H:56
Foam::cellListPMG
Definition: cellListPMG.H:49
Foam::polyMeshGenAddressing
Definition: polyMeshGenAddressing.H:69
Foam::HashSet< label, Hash< label > >
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::LongList< label >
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::polyMeshGenCells::addCellToSubset
void addCellToSubset(const label, const label)
Definition: polyMeshGenCellsI.H:45
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::Info
messageStream Info
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::polyMeshGenCells::addCellSubset
label addCellSubset(const word &)
Definition: polyMeshGenCells.C:351
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
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
labelledPoint.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
helperFunctions.H
Foam::labelledPoint
Definition: labelledPoint.H:50
f
labelList f(nPoints)
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
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::polyMeshGenModifier::addCell
void addCell(const faceList &cellFaces)
Definition: polyMeshGenModifierAddCells.C:172
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::polyMeshGenAddressing::globalPointLabel
const labelLongList & globalPointLabel() const
Definition: polyMeshGenAddressingParallelAddressing.C:707
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::polyMeshGenModifier::facesAccess
faceListPMG & facesAccess()
access to mesh faces
Definition: polyMeshGenModifier.H:113
Foam::polyMeshGenModifier::addBufferCells
void addBufferCells()
Definition: polyMeshGenModifierAddBufferCells.C:43
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::pointFieldPMG
Definition: pointFieldPMG.H:50
polyMeshGenAddressing.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56