tetMeshOptimisationParallel.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 "demandDrivenData.H"
29 #include "tetMeshOptimisation.H"
30 #include "partTetMesh.H"
31 #include "VRWGraph.H"
32 #include "helperFunctionsPar.H"
33 
34 // #define DEBUGSearch
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
44 {
45  //- make sure that processor nodes are selected on all processors
46  //- where they exist
47  const LongList<direction>& smoothVertex = tetMesh_.smoothVertex();
48  const DynList<label>& neiProcs = tetMesh_.neiProcs();
49  const labelLongList& globalPointLabel = tetMesh_.globalPointLabel();
50  const VRWGraph& pProcs = tetMesh_.pointAtProcs();
51  const Map<label>& globalToLocal = tetMesh_.globalToLocalPointAddressing();
52  const labelLongList& pAtParallelBoundaries =
54 
55  std::map<label, labelLongList> selectedNegativeNodes;
56  forAll(neiProcs, procI)
57  selectedNegativeNodes.insert
58  (
59  std::make_pair(neiProcs[procI], labelLongList())
60  );
61 
62  forAll(pAtParallelBoundaries, i)
63  {
64  const label nI = pAtParallelBoundaries[i];
65 
66  if( !negativeNode[nI] )
67  continue;
68  if( !(smoothVertex[nI] & partTetMesh::PARALLELBOUNDARY) )
69  continue;
70 
71  forAllRow(pProcs, nI, procI)
72  {
73  const label neiProc = pProcs(nI, procI);
74 
75  if( neiProc == Pstream::myProcNo() )
76  continue;
77 
78  selectedNegativeNodes[neiProc].append(globalPointLabel[nI]);
79  }
80  }
81 
82  labelLongList receivedNodes;
83  help::exchangeMap(selectedNegativeNodes, receivedNodes);
84 
85  forAll(receivedNodes, i)
86  negativeNode[globalToLocal[receivedNodes[i]]] = true;
87 }
88 
90 (
91  std::map<label, DynList<parPartTet> >& m,
92  boolList* negativeNodePtr
93 )
94 {
95  const DynList<label>& neiProcs = tetMesh_.neiProcs();
96  const labelLongList& globalPointLabel = tetMesh_.globalPointLabel();
97  const VRWGraph& pProcs = tetMesh_.pointAtProcs();
98  const Map<label>& globalToLocal = tetMesh_.globalToLocalPointAddressing();
99  const labelLongList& pAtParallelBoundaries =
100  tetMesh_.pointsAtProcessorBoundaries();
101 
102  const LongList<point>& points = tetMesh_.points();
103  const LongList<direction>& smoothVertex = tetMesh_.smoothVertex();
104  const LongList<partTet>& tets = tetMesh_.tets();
105 
106  //- create the map holding data which will be sent to other procs
107  std::map<label, LongList<parPartTet> > exchangeData;
108  forAll(neiProcs, procI)
109  exchangeData.insert
110  (
111  std::make_pair(neiProcs[procI], LongList<parPartTet>())
112  );
113 
114  //- create storage in the m map
115  m.clear();
116  forAll(pAtParallelBoundaries, i)
117  {
118  const label pI = pAtParallelBoundaries[i];
119  if( !(smoothVertex[pI] & partTetMesh::SMOOTH) )
120  continue;
121  if( negativeNodePtr && !(*negativeNodePtr)[pI] )
122  continue;
123 
124  m.insert(std::make_pair(pI, DynList<parPartTet>()));
125  }
126 
127  //- store local data into the maps
128  forAll(tets, tetI)
129  {
130  const partTet& tet = tets[tetI];
131 
132  DynList<label> sendToProcs;
133  for(label i=0;i<4;++i)
134  {
135  const label vI = tet[i];
136 
137  if( pProcs.sizeOfRow(vI) == 0 )
138  continue;
139  if( !(smoothVertex[vI] & partTetMesh::SMOOTH) )
140  continue;
141  if( negativeNodePtr && !(*negativeNodePtr)[vI] )
142  continue;
143 
144  //- add processor labels
145  forAllRow(pProcs, vI, procI)
146  sendToProcs.appendIfNotIn(pProcs(vI, procI));
147 
148  //- add data into the map of proc bnd points
149  DynList<parPartTet>& data = m[vI];
150  data.append
151  (
152  parPartTet
153  (
154  labelledPoint(globalPointLabel[tet[0]], points[tet[0]]),
155  labelledPoint(globalPointLabel[tet[1]], points[tet[1]]),
156  labelledPoint(globalPointLabel[tet[2]], points[tet[2]]),
157  labelledPoint(globalPointLabel[tet[3]], points[tet[3]])
158  )
159  );
160  }
161 
162  if( sendToProcs.size() != 0 )
163  {
164  //- add data into the map
165  forAll(sendToProcs, procI)
166  {
167  const label neiProc = sendToProcs[procI];
168  if( neiProc == Pstream::myProcNo() )
169  continue;
170 
171  exchangeData[neiProc].append
172  (
173  parPartTet
174  (
175  labelledPoint(globalPointLabel[tet[0]], points[tet[0]]),
176  labelledPoint(globalPointLabel[tet[1]], points[tet[1]]),
177  labelledPoint(globalPointLabel[tet[2]], points[tet[2]]),
178  labelledPoint(globalPointLabel[tet[3]], points[tet[3]])
179  )
180  );
181  }
182  }
183  }
184 
185  //- exchange data with other processors
186  LongList<parPartTet> receivedData;
187  help::exchangeMap(exchangeData, receivedData);
188 
189  forAll(receivedData, i)
190  {
191  const parPartTet& tet = receivedData[i];
192 
193  for(label i=0;i<4;++i)
194  {
195  const label gpI = tet[i].pointLabel();
196 
197  if(
198  globalToLocal.found(gpI) &&
199  (smoothVertex[globalToLocal[gpI]] &partTetMesh::SMOOTH)
200  )
201  {
202  DynList<parPartTet>& data = m[globalToLocal[gpI]];
203  data.append(tet);
204  }
205  }
206  }
207 }
208 
210 {
212  const labelLongList& bufferLayerPoints = tetMesh_.bufferLayerPoints();
213  const VRWGraph& pProcs = tetMesh_.pointAtProcs();
214  const labelLongList& globalPointLabel = tetMesh_.globalPointLabel();
215  const Map<label>& globalToLocal = tetMesh_.globalToLocalPointAddressing();
216  const DynList<label>& neiProcs = tetMesh_.neiProcs();
217 
218  //- create the map
219  std::map<label, LongList<labelledPoint> > exchangeData;
220  forAll(neiProcs, i)
221  exchangeData.insert
222  (
223  std::make_pair(neiProcs[i], LongList<labelledPoint>())
224  );
225 
226  //- add points into the map
227  forAll(bufferLayerPoints, pI)
228  {
229  const label pointI = bufferLayerPoints[pI];
230 
231  forAllRow(pProcs, pointI, i)
232  {
233  const label neiProc = pProcs(pointI, i);
234 
235  if( neiProc == Pstream::myProcNo() )
236  continue;
237 
238  exchangeData[neiProc].append
239  (
240  labelledPoint(globalPointLabel[pointI], points[pointI])
241  );
242  }
243  }
244 
245  LongList<labelledPoint> receivedData;
246  help::exchangeMap(exchangeData, receivedData);
247 
248  forAll(receivedData, i)
249  {
250  const labelledPoint& lp = receivedData[i];
251 
253  (
254  globalToLocal[lp.pointLabel()],
255  lp.coordinates()
256  );
257  }
258 }
259 
261 (
262  const boolList* negativeNodePtr
263 )
264 {
265  const LongList<point>& points = tetMesh_.points();
266  const DynList<label>& neiProcs = tetMesh_.neiProcs();
267  const VRWGraph& pProcs = tetMesh_.pointAtProcs();
268  const Map<label>& globalToLocal = tetMesh_.globalToLocalPointAddressing();
269  const labelLongList& globalPointLabel = tetMesh_.globalPointLabel();
270  const LongList<direction>& smoothVertex = tetMesh_.smoothVertex();
271  const labelLongList& pAtParallelBoundaries =
272  tetMesh_.pointsAtProcessorBoundaries();
273 
274  //- create the map
275  std::map<label, LongList<labelledPoint> > exchangeData;
276  forAll(neiProcs, procI)
277  exchangeData.insert
278  (
279  std::make_pair(neiProcs[procI], LongList<labelledPoint>())
280  );
281 
282  //- fill in the data
283  std::map<label, labelledPoint> parallelBndPoints;
284  forAll(pAtParallelBoundaries, i)
285  {
286  const label pI = pAtParallelBoundaries[i];
287 
288  if( !(smoothVertex[pI] & partTetMesh::PARALLELBOUNDARY) )
289  continue;
290 
291  parallelBndPoints[pI] = labelledPoint(1, points[pI]);
292 
293  if( negativeNodePtr && !(*negativeNodePtr)[pI] )
294  continue;
295 
296  forAllRow(pProcs, pI, procI)
297  {
298  const label neiProc = pProcs(pI, procI);
299  if( neiProc == Pstream::myProcNo() )
300  continue;
301 
302  exchangeData[neiProc].append
303  (
304  labelledPoint(globalPointLabel[pI], points[pI])
305  );
306  }
307  }
308 
309  //- send points to other processors
310  LongList<labelledPoint> receivedData;
311  help::exchangeMap(exchangeData, receivedData);
312 
313  //- gather the data
314  forAll(receivedData, i)
315  {
316  const labelledPoint& lp = receivedData[i];
317 
318  std::map<label, labelledPoint>::iterator iter =
319  parallelBndPoints.find(globalToLocal[lp.pointLabel()]);
320 
321  if( iter == parallelBndPoints.end() )
322  continue;
323 
324  ++iter->second.pointLabel();
325  iter->second.coordinates() += lp.coordinates();
326  }
327 
328  //- move the point to the averaged position
329  for
330  (
331  std::map<label, labelledPoint>::iterator it=parallelBndPoints.begin();
332  it!=parallelBndPoints.end();
333  ++it
334  )
335  {
336  const label pI = it->first;
337 
338  const point newP = it->second.coordinates() / it->second.pointLabel();
339  tetMesh_.updateVertex(pI, newP);
340  }
341 }
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 } // End namespace Foam
346 
347 // ************************************************************************* //
Foam::labelledPoint::coordinates
const point & coordinates() const
return point coordinates
Definition: labelledPoint.H:93
Foam::parPartTet
Definition: parPartTet.H:49
helperFunctionsPar.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
VRWGraph.H
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
tetMeshOptimisation.H
Foam::tetMeshOptimisation::exchangeData
void exchangeData(std::map< label, DynList< parPartTet > > &m, boolList *negativeNodePtr=NULL)
exchange data with other processors
Definition: tetMeshOptimisationParallel.C:90
Foam::partTet
Definition: partTet.H:53
Foam::Map< label >
Foam::LongList< direction >
Foam::tetMeshOptimisation::tetMesh_
partTetMesh & tetMesh_
reference to the tet mesh
Definition: tetMeshOptimisation.H:62
Foam::partTetMesh::pointAtProcs
const VRWGraph & pointAtProcs() const
Definition: partTetMesh.H:218
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::partTetMesh::globalToLocalPointAddressing
const Map< label > & globalToLocalPointAddressing() const
Definition: partTetMesh.H:226
Foam::DynList::appendIfNotIn
void appendIfNotIn(const T &e)
Definition: DynListI.H:324
partTetMesh.H
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::partTetMesh::SMOOTH
@ SMOOTH
Definition: partTetMesh.H:161
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::DynList< label >
Foam::partTetMesh::PARALLELBOUNDARY
@ PARALLELBOUNDARY
Definition: partTetMesh.H:164
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::partTetMesh::bufferLayerPoints
const labelLongList & bufferLayerPoints() const
Definition: partTetMesh.H:250
Foam::tetMeshOptimisation::updateBufferLayerPoints
void updateBufferLayerPoints()
update buffer layer points
Definition: tetMeshOptimisationParallel.C:209
Foam::partTetMesh::pointsAtProcessorBoundaries
const labelLongList & pointsAtProcessorBoundaries() const
Definition: partTetMesh.H:242
Foam::labelledPoint
Definition: labelledPoint.H:50
Foam::Vector< scalar >
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::help::exchangeMap
void exchangeMap(const std::map< label, ListType > &m, LongList< T > &data, const Pstream::commsTypes commsType)
Definition: helperFunctionsPar.C:129
Foam::tetMeshOptimisation::unifyNegativePoints
void unifyNegativePoints(boolList &negativeNode) const
Definition: tetMeshOptimisationParallel.C:43
Foam::DynList::size
label size() const
Definition: DynListI.H:235
Foam::partTetMesh::neiProcs
const DynList< label > & neiProcs() const
Definition: partTetMesh.H:234
Foam::partTetMesh::points
const LongList< point > & points() const
access to points, tets and other data
Definition: partTetMesh.H:174
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::labelLongList
LongList< label > labelLongList
Definition: labelLongList.H:46
Foam::partTetMesh::smoothVertex
const LongList< direction > & smoothVertex() const
Definition: partTetMesh.H:189
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52
Foam::tetMeshOptimisation::unifyCoordinatesParallel
void unifyCoordinatesParallel(const boolList *negativeNodePtr=NULL)
Definition: tetMeshOptimisationParallel.C:261
Foam::partTetMesh::updateVertex
void updateVertex(const label pointI, const point &newP)
move the vertex to a new position
Definition: partTetMesh.C:381
Foam::labelledPoint::pointLabel
label pointLabel() const
return point label
Definition: labelledPoint.H:82
Foam::partTetMesh::globalPointLabel
const labelLongList & globalPointLabel() const
Definition: partTetMesh.H:210
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304