polyMeshGenAddressingParallelAddressing.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 \*---------------------------------------------------------------------------*/
25 
26 #include "polyMeshGenAddressing.H"
27 #include "demandDrivenData.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
37 {
38  if( !Pstream::parRun() )
40  (
41  "void polyMeshGenAddressing::calcGlobalPointLabels() const"
42  ) << "Cannot calculate global point addressing "
43  << "for a serial run!" << exit(FatalError);
44 
49  globalPointLabel = -1;
50 
51  if( !pProcsPtr_ )
52  pProcsPtr_ = new VRWGraph();
53  VRWGraph& pProcs = *pProcsPtr_;
54  pProcs.setSize(mesh_.points().size());
55 
59 
60  if( !pointNeiProcsPtr_ )
62  DynList<label>& pNeiProcs = *pointNeiProcsPtr_;
63 
64  const faceListPMG& faces = mesh_.faces();
65  const PtrList<processorBoundaryPatch>& procBoundaries =
67 
68  //- find which processors contain a given bnd point
70  patchPoints.setSize(procBoundaries.size());
71  forAll(procBoundaries, patchI)
72  {
73  const label start = procBoundaries[patchI].patchStart();
74  const label end = start + procBoundaries[patchI].patchSize();
75 
76  std::map<label, std::pair<label, label> >& patchPointsMap =
77  patchPoints[patchI];
78 
79  for(label faceI=start;faceI<end;++faceI)
80  {
81  const face& f = faces[faceI];
82 
83  forAll(f, pI)
84  {
85  std::map<label, std::pair<label, label> >::iterator it =
86  patchPointsMap.find(f[pI]);
87 
88  if( it != patchPointsMap.end() )
89  continue;
90 
91  const std::pair<label, label> pp
92  (
93  faceI-start,
94  (f.size()-pI)%f.size()
95  );
96  patchPointsMap.insert(std::make_pair(f[pI], pp));
97 
98  pProcs.appendIfNotIn(f[pI], procBoundaries[patchI].myProcNo());
99  pProcs.appendIfNotIn(f[pI], procBoundaries[patchI].neiProcNo());
100  }
101  }
102  }
103 
104  //- exchange data with other processor and update bProcs
105  //- continue this process as long as there is some change
106  bool finished;
107  do
108  {
109  finished = true;
110 
111  forAll(procBoundaries, patchI)
112  {
113  const std::map<label, std::pair<label, label> >& patchPointsMap =
114  patchPoints[patchI];
115  std::map<label, std::pair<label, label> >::const_iterator it;
116 
117  labelLongList dataToSend;
118  for(it=patchPointsMap.begin();it!=patchPointsMap.end();++it)
119  {
120  //- data is sent as follows
121  //- 1. face position in patch
122  //- 2. local point position in face
123  //- 3. number of patches for point
124  //- 4. patch labels
125  dataToSend.append(it->second.first);
126  dataToSend.append(it->second.second);
127  dataToSend.append(pProcs.sizeOfRow(it->first));
128  forAllRow(pProcs, it->first, i)
129  dataToSend.append(pProcs(it->first, i));
130  }
131 
132  OPstream toOtherProc
133  (
135  procBoundaries[patchI].neiProcNo(),
136  dataToSend.byteSize()
137  );
138  toOtherProc << dataToSend;
139  }
140 
141  forAll(procBoundaries, patchI)
142  {
143  IPstream fromOtherProc
144  (
146  procBoundaries[patchI].neiProcNo()
147  );
148  labelList receivedData;
149  fromOtherProc >> receivedData;
150 
151  const label start = procBoundaries[patchI].patchStart();
152 
153  label counter(0);
154  while( counter < receivedData.size() )
155  {
156  const face& f = faces[start+receivedData[counter++]];
157  const label pI = receivedData[counter++];
158  const label nProcs = receivedData[counter++];
159  for(label i=0;i<nProcs;++i)
160  {
161  const label neiProc = receivedData[counter++];
162  if( !pProcs.contains(f[pI], neiProc) )
163  {
164  pProcs.append(f[pI], neiProc);
165  finished = false;
166  }
167  }
168  }
169  }
170 
171  reduce(finished, minOp<bool>());
172  } while( !finished );
173 
174  //- start calculating global point labels
175  label nLocalPoints(0);
176  boolList localPoints(mesh_.points().size(), true);
177  forAll(pProcs, pI)
178  if( pProcs.sizeOfRow(pI) != 0 )
179  {
180  label pMin = pProcs(pI, 0);
181  forAllRow(pProcs, pI, procI)
182  pMin = Foam::min(pMin, pProcs(pI, procI));
183 
184  if( pMin == Pstream::myProcNo() )
185  {
186  ++nLocalPoints;
187  }
188  else
189  {
190  localPoints[pI] = false;
191  }
192  }
193  else
194  {
195  ++nLocalPoints;
196  }
197 
198  //- give local points their labels
199  label startPoint(0);
200  labelList nPointsAtProc(Pstream::nProcs());
201  nPointsAtProc[Pstream::myProcNo()] = nLocalPoints;
202  Pstream::gatherList(nPointsAtProc);
203  Pstream::scatterList(nPointsAtProc);
204  for(label i=0;i<Pstream::myProcNo();++i)
205  startPoint += nPointsAtProc[i];
206 
207  forAll(localPoints, pI)
208  if( localPoints[pI] )
209  globalPointLabel[pI] = startPoint++;
210 
211  //- send labels to non-local points
212  do
213  {
214  finished = true;
215 
216  forAll(procBoundaries, patchI)
217  {
218  const std::map<label, std::pair<label, label> >& patchPointsMap =
219  patchPoints[patchI];
220  std::map<label, std::pair<label, label> >::const_iterator it;
221 
222  labelLongList dataToSend;
223  for(it=patchPointsMap.begin();it!=patchPointsMap.end();++it)
224  {
225  if( globalPointLabel[it->first] != -1 )
226  {
227  //- data is sent as follows
228  //- 1. face position in patch
229  //- 2. local point position in face
230  //- 3. global point label
231  dataToSend.append(it->second.first);
232  dataToSend.append(it->second.second);
233  dataToSend.append(globalPointLabel[it->first]);
234  }
235  }
236 
237  OPstream toOtherProc
238  (
240  procBoundaries[patchI].neiProcNo(),
241  dataToSend.byteSize()
242  );
243  toOtherProc << dataToSend;
244  }
245 
246  forAll(procBoundaries, patchI)
247  {
248  IPstream fromOtherProc
249  (
251  procBoundaries[patchI].neiProcNo()
252  );
253  labelList receivedData;
254  fromOtherProc >> receivedData;
255 
256  const label start = procBoundaries[patchI].patchStart();
257 
258  label counter(0);
259  while( counter < receivedData.size() )
260  {
261  const face& f = faces[start+receivedData[counter++]];
262  const label pI = receivedData[counter++];
263  const label globalLabel = receivedData[counter++];
264 
265  if( globalPointLabel[f[pI]] == -1 )
266  {
267  globalPointLabel[f[pI]] = globalLabel;
268  finished = false;
269  }
270  else if( globalPointLabel[f[pI]] != globalLabel )
271  {
273  (
274  "void polyMeshGenAddressing::"
275  "calcGlobalPointLabels() const"
276  ) << "Point labels on processors do not match!"
277  << abort(FatalError);
278  }
279  }
280  }
281 
282  reduce(finished, minOp<bool>());
283  } while( !finished );
284 
285  //- create global to local addressing and pointNeiProcs addressing
286  forAll(globalPointLabel, pointI)
287  if( pProcs.sizeOfRow(pointI) != 0 )
288  {
289  globalToLocal.insert(globalPointLabel[pointI], pointI);
290 
291  forAllRow(pProcs, pointI, i)
292  pNeiProcs.appendIfNotIn(pProcs(pointI, i));
293  }
294 }
295 
297 {
298  if( !globalFaceLabelPtr_ )
302  globalFaceLabel = -1;
303 
304  if( !Pstream::parRun() )
305  return;
306 
307  const label nIntFaces = mesh_.nInternalFaces();
308 
309  const PtrList<processorBoundaryPatch>& procBoundaries =
311 
312  label nBoundaryFaces = 0;
313  forAll(procBoundaries, patchI)
314  {
315  if( procBoundaries[patchI].owner() )
316  nBoundaryFaces += procBoundaries[patchI].patchSize();
317  }
318 
319  label startFace(0);
320  labelList numberOfInternalFaces(Pstream::nProcs());
321  numberOfInternalFaces[Pstream::myProcNo()] = nIntFaces + nBoundaryFaces;
322 
323  Pstream::gatherList(numberOfInternalFaces);
324  Pstream::scatterList(numberOfInternalFaces);
325  for(label i=0;i<Pstream::myProcNo();++i)
326  startFace += numberOfInternalFaces[i];
327 
328  //- calculate labels for internal faces
329  for(label faceI=0;faceI<nIntFaces;++faceI)
330  globalFaceLabel[faceI] = startFace++;
331 
332  //- calculate labels for processor boundaries
333  forAll(procBoundaries, patchI)
334  {
335  if( procBoundaries[patchI].owner() )
336  {
337  const label start = procBoundaries[patchI].patchStart();
338  labelList dataToSend(procBoundaries[patchI].patchSize());
339  forAll(dataToSend, bfI)
340  {
341  globalFaceLabel[start+bfI] = startFace;
342  dataToSend[bfI] = startFace++;
343  }
344 
345  OPstream toOtherProc
346  (
348  procBoundaries[patchI].neiProcNo(),
349  dataToSend.byteSize()
350  );
351  toOtherProc << dataToSend;
352  }
353  }
354 
355  forAll(procBoundaries, patchI)
356  {
357  if( !procBoundaries[patchI].owner() )
358  {
359  IPstream fromOtherProc
360  (
362  procBoundaries[patchI].neiProcNo()
363  );
364 
365  labelList receivedLabels;
366  fromOtherProc >> receivedLabels;
367 
368  const label start = procBoundaries[patchI].patchStart();
369  forAll(receivedLabels, bfI)
370  {
371  globalFaceLabel[start+bfI] = receivedLabels[bfI];
372  }
373  }
374  }
375 
376  //- create global labels for boundary faces
377  startFace = 0;
378  forAll(numberOfInternalFaces, procI)
379  startFace += numberOfInternalFaces[procI];
380 
381  const PtrList<boundaryPatch>& boundaries = mesh_.boundaries();
382  labelListList numberOfBoundaryFaces(Pstream::nProcs());
383  numberOfBoundaryFaces[Pstream::myProcNo()].setSize(boundaries.size());
384  forAll(boundaries, patchI)
385  {
386  numberOfBoundaryFaces[Pstream::myProcNo()][patchI] =
387  boundaries[patchI].patchSize();
388  }
389  Pstream::gatherList(numberOfBoundaryFaces);
390  Pstream::scatterList(numberOfBoundaryFaces);
391 
392  forAll(boundaries, patchI)
393  {
394  const label start = boundaries[patchI].patchStart();
395  const label nBoundaryFaces = boundaries[patchI].patchSize();
396 
397  for(label procI=0;procI<Pstream::myProcNo();++procI)
398  startFace += numberOfBoundaryFaces[procI][patchI];
399 
400  for(label bfI=0;bfI<nBoundaryFaces;++bfI)
401  globalFaceLabel[start+bfI] = startFace++;
402 
403  for(label procI=Pstream::myProcNo()+1;procI<Pstream::nProcs();++procI)
404  startFace += numberOfBoundaryFaces[procI][patchI];
405  }
406 }
407 
409 {
410  if( !globalCellLabelPtr_ )
414  globalCellLabel = -1;
415 
416  if( !Pstream::parRun() )
417  return;
418 
419  label startLabel(0);
420  labelList nCellsAtProc(Pstream::nProcs());
421  nCellsAtProc[Pstream::myProcNo()] = globalCellLabel.size();
422  Pstream::gatherList(nCellsAtProc);
423  Pstream::scatterList(nCellsAtProc);
424  for(label i=0;i<Pstream::myProcNo();++i)
425  startLabel += nCellsAtProc[i];
426 
427  forAll(globalCellLabel, cellI)
428  globalCellLabel[cellI] = startLabel++;
429 }
430 
432 {
433  if( !Pstream::parRun() )
435  (
436  "void polyMeshGenAddressing::calcGlobalEdgeLabels() const"
437  ) << "Cannot calculate global edge addressing "
438  << "for a serial run?!?!" << exit(FatalError);
439 
440  if( !globalEdgeLabelPtr_ )
443  globalEdgeLabel.setSize(this->edges().size());
444  globalEdgeLabel = -1;
445 
446  if( !eProcsPtr_ )
447  eProcsPtr_ = new VRWGraph();
448  VRWGraph& eProcs = *eProcsPtr_;
449  eProcs.setSize(this->edges().size());
450 
454 
455  if( !edgeNeiProcsPtr_ )
457  DynList<label>& eNeiProcs = *edgeNeiProcsPtr_;
458 
459  const faceListPMG& faces = mesh_.faces();
460  const VRWGraph& faceEdges = this->faceEdges();
461  const PtrList<processorBoundaryPatch>& procBoundaries =
463 
464  //- find which processor contain a given bnd edge
465  forAll(procBoundaries, patchI)
466  {
467  const label start = procBoundaries[patchI].patchStart();
468  const label end = start + procBoundaries[patchI].patchSize();
469  for(label faceI=start;faceI<end;++faceI)
470  {
471  forAllRow(faceEdges, faceI, eI)
472  {
473  eProcs.appendIfNotIn
474  (
475  faceEdges(faceI, eI),
476  procBoundaries[patchI].myProcNo()
477  );
478  eProcs.appendIfNotIn
479  (
480  faceEdges(faceI, eI),
481  procBoundaries[patchI].neiProcNo()
482  );
483  }
484  }
485  }
486 
487  //- exchange data with other processor and update eProcs
488  //- continue this process as long as there is some change
489  bool finished;
490  do
491  {
492  finished = true;
493 
494  forAll(procBoundaries, patchI)
495  {
496  const label start = procBoundaries[patchI].patchStart();
497  const label end = start + procBoundaries[patchI].patchSize();
498 
499 /* label nToSend(0);
500  for(label faceI=start;faceI<end;++faceI)
501  {
502  forAllRow(faceEdges, faceI, eI)
503  {
504  nToSend += 3;
505  nToSend += eProcs.sizeOfRow(faceEdges(faceI, eI));
506  }
507  }
508 */
509  labelLongList dataToSend;
510  //nToSend = 0;
511  for(label faceI=start;faceI<end;++faceI)
512  {
513  forAllRow(faceEdges, faceI, eI)
514  {
515  const label edgeI = faceEdges(faceI, eI);
516  const face& f = faces[faceI];
517  //- data is sent as follows
518  //- 1. face position in patch
519  //- 2. local edge position in face
520  //- 3. number of processors for edge
521  //- 4. processor labels
522  dataToSend.append(faceI-start);
523  dataToSend.append((f.size()-eI-1)%f.size());
524  dataToSend.append(eProcs.sizeOfRow(edgeI));
525  forAllRow(eProcs, edgeI, i)
526  dataToSend.append(eProcs(edgeI, i));
527  }
528  }
529 
530  OPstream toOtherProc
531  (
533  procBoundaries[patchI].neiProcNo(),
534  dataToSend.byteSize()
535  );
536  toOtherProc << dataToSend;
537  }
538 
539  forAll(procBoundaries, patchI)
540  {
541  IPstream fromOtherProc
542  (
544  procBoundaries[patchI].neiProcNo()
545  );
546  labelList receivedData;
547  fromOtherProc >> receivedData;
548 
549  const label start = procBoundaries[patchI].patchStart();
550 
551  label counter(0);
552  while( counter < receivedData.size() )
553  {
554  const label faceI = start+receivedData[counter++];
555  const label eI = receivedData[counter++];
556 
557  const label edgeI = faceEdges(faceI, eI);
558 
559  const label nProcs = receivedData[counter++];
560  for(label i=0;i<nProcs;++i)
561  {
562  const label neiProc = receivedData[counter++];
563  if( !eProcs.contains(edgeI, neiProc) )
564  {
565  eProcs.append(edgeI, neiProc);
566  finished = false;
567  }
568  }
569  }
570  }
571 
572  reduce(finished, minOp<bool>());
573  } while( !finished );
574 
575  //- start calculating global edge labels
576  label nLocalEdges(0);
577  boolList localEdges(eProcs.size(), true);
578  forAll(eProcs, edgeI)
579  if( eProcs.sizeOfRow(edgeI) != 0 )
580  {
581  label pMin = eProcs(edgeI, 0);
582  forAllRow(eProcs, edgeI, procI)
583  if( eProcs(edgeI, procI) < pMin )
584  pMin = eProcs(edgeI, procI);
585 
586  if( pMin == Pstream::myProcNo() )
587  {
588  ++nLocalEdges;
589  }
590  else
591  {
592  localEdges[edgeI] = false;
593  }
594  }
595  else
596  {
597  ++nLocalEdges;
598  }
599 
600  //- give local points their labels
601  label startEdge(0);
602  labelList nEdgesAtProc(Pstream::nProcs());
603  nEdgesAtProc[Pstream::myProcNo()] = nLocalEdges;
604  Pstream::gatherList(nEdgesAtProc);
605  Pstream::scatterList(nEdgesAtProc);
606  for(label i=0;i<Pstream::myProcNo();++i)
607  startEdge += nEdgesAtProc[i];
608 
609  forAll(localEdges, pI)
610  if( localEdges[pI] )
611  globalEdgeLabel[pI] = startEdge++;
612 
613  //- send labels to non-local edges
614  do
615  {
616  finished = true;
617 
618  forAll(procBoundaries, patchI)
619  {
620  const label start = procBoundaries[patchI].patchStart();
621  const label end = start + procBoundaries[patchI].patchSize();
622 
623  labelLongList dataToSend;
624  for(label faceI=start;faceI<end;++faceI)
625  {
626  forAllRow(faceEdges, faceI, eI)
627  {
628  const label edgeI = faceEdges(faceI, eI);
629  if( globalEdgeLabel[edgeI] != -1 )
630  {
631  const face& f = faces[faceI];
632  //- data is sent as follows
633  //- 1. face position in patch
634  //- 2. local edge position in face
635  //- 3. number global label for the edge
636  dataToSend.append(faceI-start);
637  dataToSend.append((f.size()-eI-1)%f.size());
638  dataToSend.append(globalEdgeLabel[edgeI]);
639  }
640  }
641  }
642 
643  OPstream toOtherProc
644  (
646  procBoundaries[patchI].neiProcNo(),
647  dataToSend.byteSize()
648  );
649  toOtherProc << dataToSend;
650  }
651 
652  forAll(procBoundaries, patchI)
653  {
654  IPstream fromOtherProc
655  (
657  procBoundaries[patchI].neiProcNo()
658  );
659  labelList receivedData;
660  fromOtherProc >> receivedData;
661 
662  const label start = procBoundaries[patchI].patchStart();
663 
664  label counter(0);
665  while( counter < receivedData.size() )
666  {
667  const label faceI = start+receivedData[counter++];
668  const label eI = receivedData[counter++];
669 
670  const label edgeI = faceEdges(faceI, eI);
671 
672  const label globalLabel = receivedData[counter++];
673 
674  if( globalEdgeLabel[edgeI] == -1 )
675  {
676  globalEdgeLabel[edgeI] = globalLabel;
677  finished = false;
678  }
679  else if( globalEdgeLabel[edgeI] != globalLabel )
680  {
682  (
683  "void polyMeshGenAddressing::"
684  "calcGlobalEdgeLabels() const"
685  ) << "Edge labels on processors do not match!"
686  << abort(FatalError);
687  }
688  }
689  }
690 
691  reduce(finished, minOp<bool>());
692  } while( !finished );
693 
694  //- create globalToLocal addressing
695  forAll(eProcs, edgeI)
696  if( eProcs.sizeOfRow(edgeI) != 0 )
697  {
698  globalToLocal.insert(globalEdgeLabel[edgeI], edgeI);
699 
700  forAllRow(eProcs, edgeI, i)
701  eNeiProcs.appendIfNotIn(eProcs(edgeI, i));
702  }
703 }
704 
705 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
706 
708 {
710  {
711  # ifdef USE_OMP
712  if( omp_in_parallel() )
714  (
715  "const labelLongList& polyMeshGenAddressing::"
716  "globalPointLabel() const"
717  ) << "Calculating addressing inside a parallel region."
718  << " This is not thread safe" << exit(FatalError);
719  # endif
720 
722  }
723 
724  return *globalPointLabelPtr_;
725 }
726 
728 {
729  if( !globalFaceLabelPtr_ )
730  {
731  # ifdef USE_OMP
732  if( omp_in_parallel() )
734  (
735  "const labelLongList& polyMeshGenAddressing"
736  "::globalFaceLabel() const"
737  ) << "Calculating addressing inside a parallel region."
738  << " This is not thread safe" << exit(FatalError);
739  # endif
740 
742  }
743 
744  return *globalFaceLabelPtr_;
745 }
746 
748 {
749  if( !globalCellLabelPtr_ )
751 
752  return *globalCellLabelPtr_;
753 }
754 
756 {
757  if( !globalEdgeLabelPtr_ )
758  {
759  # ifdef USE_OMP
760  if( omp_in_parallel() )
762  (
763  "const labelLongList& polyMeshGenAddressing"
764  "::globalEdgeLabel() const"
765  ) << "Calculating addressing inside a parallel region."
766  << " This is not thread safe" << exit(FatalError);
767  # endif
768 
770  }
771 
772  return *globalEdgeLabelPtr_;
773 }
774 
776 {
777  if( !globalPointLabelPtr_ )
778  {
779  # ifdef USE_OMP
780  if( omp_in_parallel() )
782  (
783  "const VRWGraph& polyMeshGenAddressing::pointAtProcs() const"
784  ) << "Calculating addressing inside a parallel region."
785  << " This is not thread safe" << exit(FatalError);
786  # endif
787 
789  }
790 
791  return *pProcsPtr_;
792 }
793 
795 {
796  if( !pointNeiProcsPtr_ )
797  {
798  # ifdef USE_OMP
799  if( omp_in_parallel() )
801  (
802  "const DynList<label>& polyMeshGenAddressing"
803  "::pointNeiProcs() const"
804  ) << "Calculating addressing inside a parallel region."
805  << " This is not thread safe" << exit(FatalError);
806  # endif
807 
809  }
810 
811  return *pointNeiProcsPtr_;
812 }
813 
815 {
817  {
818  # ifdef USE_OMP
819  if( omp_in_parallel() )
821  (
822  "const Map<label>& polyMeshGenAddressing"
823  "::globalToLocalPointAddressing() const"
824  ) << "Calculating addressing inside a parallel region."
825  << " This is not thread safe" << exit(FatalError);
826  # endif
827 
829  }
830 
832 }
833 
835 {
836  if( !globalEdgeLabelPtr_ )
837  {
838  # ifdef USE_OMP
839  if( omp_in_parallel() )
841  (
842  "const VRWGraph& polyMeshGenAddressing::edgeAtProcs() const"
843  ) << "Calculating addressing inside a parallel region."
844  << " This is not thread safe" << exit(FatalError);
845  # endif
846 
848  }
849 
850  return *eProcsPtr_;
851 }
852 
854 {
855  if( !edgeNeiProcsPtr_ )
856  {
857  # ifdef USE_OMP
858  if( omp_in_parallel() )
860  (
861  "const DynList<label>& polyMeshGenAddressing::edgeNeiProcs() const"
862  ) << "Calculating addressing inside a parallel region."
863  << " This is not thread safe" << exit(FatalError);
864  # endif
865 
867  }
868 
869  return *edgeNeiProcsPtr_;
870 }
871 
873 {
875  {
876  # ifdef USE_OMP
877  if( omp_in_parallel() )
879  (
880  "const Map<label>& polyMeshGenAddressing"
881  "::globalToLocalEdgeAddressing() const"
882  ) << "Calculating addressing inside a parallel region."
883  << " This is not thread safe" << exit(FatalError);
884  # endif
885 
887  }
888 
890 }
891 
892 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
893 
894 } // End namespace Foam
895 
896 // ************************************************************************* //
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::polyMeshGenAddressing::globalFaceLabel
const labelLongList & globalFaceLabel() const
Definition: polyMeshGenAddressingParallelAddressing.C:727
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
Foam::polyMeshGenAddressing::mesh_
const polyMeshGenCells & mesh_
reference to the mesh
Definition: polyMeshGenAddressing.H:74
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:205
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::minOp
Definition: ops.H:173
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
pMin
dimensionedScalar pMin("pMin", dimPressure, fluid)
Foam::polyMeshGenAddressing::calcGlobalCellLabels
void calcGlobalCellLabels() const
calculate global cell labels
Definition: polyMeshGenAddressingParallelAddressing.C:408
Foam::polyMeshGenAddressing::globalEdgeLabel
const labelLongList & globalEdgeLabel() const
Definition: polyMeshGenAddressingParallelAddressing.C:755
Foam::polyMeshGenAddressing::globalToLocalEdgeAddressing
const Map< label > & globalToLocalEdgeAddressing() const
Definition: polyMeshGenAddressingParallelAddressing.C:872
Foam::polyMeshGenAddressing::faceEdges
const VRWGraph & faceEdges() const
Definition: polyMeshGenAddressingFaceEdges.C:112
Foam::Map< label >
Foam::polyMeshGenPoints::points
const pointFieldPMG & points() const
access to points
Definition: polyMeshGenPointsI.H:44
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::polyMeshGenAddressing::globalFaceLabelPtr_
labelLongList * globalFaceLabelPtr_
global face labels
Definition: polyMeshGenAddressing.H:135
Foam::polyMeshGenAddressing::pointNeiProcsPtr_
DynList< label > * pointNeiProcsPtr_
neighbour processors sharing a point with this processor
Definition: polyMeshGenAddressing.H:152
Foam::VRWGraph::appendIfNotIn
void appendIfNotIn(const label rowI, const label)
Append an element to the given row if it does not exist there.
Definition: VRWGraphI.H:346
Foam::VRWGraph::contains
bool contains(const label rowI, const label e) const
check if the element is in the given row (takes linear time)
Definition: VRWGraphI.H:511
Foam::polyMeshGenFaces::faces
const faceListPMG & faces() const
access to faces
Definition: polyMeshGenFacesI.H:43
Foam::polyMeshGenAddressing::pointAtProcs
const VRWGraph & pointAtProcs() const
Definition: polyMeshGenAddressingParallelAddressing.C:775
Foam::polyMeshGenFaces::nInternalFaces
label nInternalFaces() const
return number of internal faces
Definition: polyMeshGenFacesI.H:48
Foam::LongList::setSize
void setSize(const label)
Reset size of List.
Definition: LongListI.H:223
Foam::polyMeshGenAddressing::globalPointLabelPtr_
labelLongList * globalPointLabelPtr_
global point labels
Definition: polyMeshGenAddressing.H:132
Foam::LongList< label >
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::polyMeshGenAddressing::globalToLocalEdgeAddressingPtr_
Map< label > * globalToLocalEdgeAddressingPtr_
global to local edge addressing
Definition: polyMeshGenAddressing.H:158
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::polyMeshGenAddressing::pointNeiProcs
const DynList< label > & pointNeiProcs() const
Definition: polyMeshGenAddressingParallelAddressing.C:794
Foam::VRWGraph::size
label size() const
Returns the number of rows.
Definition: VRWGraphI.H:122
Foam::polyMeshGenAddressing::edgeNeiProcsPtr_
DynList< label > * edgeNeiProcsPtr_
neighbour processors sharing an edge with the current processor
Definition: polyMeshGenAddressing.H:161
Foam::DynList::appendIfNotIn
void appendIfNotIn(const T &e)
Definition: DynListI.H:324
Foam::polyMeshGenAddressing::calcGlobalEdgeLabels
void calcGlobalEdgeLabels() const
calculate global edge labels
Definition: polyMeshGenAddressingParallelAddressing.C:431
Foam::VRWGraph::setSize
void setSize(const label)
Reset the number of rows.
Definition: VRWGraphI.H:132
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
Foam::pointFieldPMG::size
label size() const
return the number of used elements
Definition: pointFieldPMGI.H:71
Foam::polyMeshGenAddressing::edgeNeiProcs
const DynList< label > & edgeNeiProcs() const
Definition: polyMeshGenAddressingParallelAddressing.C:853
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::FatalError
error FatalError
Foam::polyMeshGenAddressing::pProcsPtr_
VRWGraph * pProcsPtr_
processors containing a vertex
Definition: polyMeshGenAddressing.H:146
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyMeshGenAddressing::eProcsPtr_
VRWGraph * eProcsPtr_
processors containing an edge
Definition: polyMeshGenAddressing.H:155
Foam::polyMeshGenAddressing::edges
const edgeList & edges() const
Return mesh edges.
Definition: polyMeshGenAddressingEdges.C:157
Foam::DynList< label >
Foam::polyMeshGenAddressing::calcGlobalPointLabels
void calcGlobalPointLabels() const
calculate global point labels
Definition: polyMeshGenAddressingParallelAddressing.C:36
Foam::polyMeshGenAddressing::calcGlobalFaceLabels
void calcGlobalFaceLabels() const
calculate global face labels
Definition: polyMeshGenAddressingParallelAddressing.C:296
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::polyMeshGenAddressing::globalToLocalPointAddressingPtr_
Map< label > * globalToLocalPointAddressingPtr_
global to local point adressing
Definition: polyMeshGenAddressing.H:149
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
Foam::LongList::byteSize
label byteSize() const
Return the binary size in number of characters of the UList.
Definition: LongListI.H:209
f
labelList f(nPoints)
Foam::faceListPMG::size
label size() const
return the number of used elements
Definition: faceListPMGI.H:73
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::polyMeshGenAddressing::globalCellLabel
const labelLongList & globalCellLabel() const
Definition: polyMeshGenAddressingParallelAddressing.C:747
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::polyMeshGenAddressing::edgeAtProcs
const VRWGraph & edgeAtProcs() const
Definition: polyMeshGenAddressingParallelAddressing.C:834
Foam::VRWGraph::append
void append(const label rowI, const label)
Append an element to the given row.
Definition: VRWGraphI.H:303
Foam::polyMeshGenFaces::procBoundaries
const PtrList< processorBoundaryPatch > & procBoundaries() const
inter-processor boundaries
Definition: polyMeshGenFacesI.H:106
Foam::polyMeshGenAddressing::globalEdgeLabelPtr_
labelLongList * globalEdgeLabelPtr_
global edge labels
Definition: polyMeshGenAddressing.H:141
Foam::polyMeshGenAddressing::globalPointLabel
const labelLongList & globalPointLabel() const
Definition: polyMeshGenAddressingParallelAddressing.C:707
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:50
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::faceListPMG
Definition: faceListPMG.H:50
Foam::VRWGraph
Definition: VRWGraph.H:101
Foam::labelLongList
LongList< label > labelLongList
Definition: labelLongList.H:46
Foam::polyMeshGenAddressing::globalCellLabelPtr_
labelLongList * globalCellLabelPtr_
global cell labels
Definition: polyMeshGenAddressing.H:138
polyMeshGenAddressing.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cellListPMG::size
label size() const
return the number of used elements
Definition: cellListPMGI.H:56