lduPrimitiveMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "lduPrimitiveMesh.H"
30 #include "processorLduInterface.H"
31 #include "edgeHashes.H"
32 #include "labelPair.H"
33 #include "processorGAMGInterface.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(lduPrimitiveMesh, 0);
40 
41  //- Less operator for pairs of <processor><index>
42  class procLess
43  {
44  const labelPairList& lst_;
45 
46  public:
47 
48  procLess(const labelPairList& lst)
49  :
50  lst_(lst)
51  {}
52 
53  bool operator()(const label a, const label b)
54  {
55  return lst_[a].first() < lst_[b].first();
56  }
57  };
58 
59 }
60 
61 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
62 
64 (
65  const label size,
66  const labelUList& l,
67  const labelUList& u
68 )
69 {
70  forAll(l, facei)
71  {
72  if (u[facei] < l[facei])
73  {
75  << "Reversed face. Problem at face " << facei
76  << " l:" << l[facei] << " u:" << u[facei]
77  << abort(FatalError);
78  }
79  if (l[facei] < 0 || u[facei] < 0 || u[facei] >= size)
80  {
82  << "Illegal cell label. Problem at face " << facei
83  << " l:" << l[facei] << " u:" << u[facei]
84  << abort(FatalError);
85  }
86  }
87 
88  for (label facei=1; facei < l.size(); facei++)
89  {
90  if (l[facei-1] > l[facei])
91  {
93  << "Lower not in incremental cell order."
94  << " Problem at face " << facei
95  << " l:" << l[facei] << " u:" << u[facei]
96  << " previous l:" << l[facei-1]
97  << abort(FatalError);
98  }
99  else if (l[facei-1] == l[facei])
100  {
101  // Same cell.
102  if (u[facei-1] > u[facei])
103  {
105  << "Upper not in incremental cell order."
106  << " Problem at face " << facei
107  << " l:" << l[facei] << " u:" << u[facei]
108  << " previous u:" << u[facei-1]
109  << abort(FatalError);
110  }
111  }
112  }
113 }
114 
115 
116 Foam::label Foam::lduPrimitiveMesh::totalSize
117 (
119 )
120 {
121  label size = 0;
122 
123  for (const lduPrimitiveMesh& msh : meshes)
124  {
125  size += msh.lduAddr().size();
126  }
127  return size;
128 }
129 
130 
132 (
133  const label nCells,
134  const labelUList& lower,
135  const labelUList& upper
136 )
137 {
138  labelList nNbrs(nCells, Zero);
139 
140  // Count number of upper neighbours
141  forAll(lower, facei)
142  {
143  if (upper[facei] < lower[facei])
144  {
146  << "Problem at face:" << facei
147  << " lower:" << lower[facei]
148  << " upper:" << upper[facei]
149  << exit(FatalError);
150  }
151  nNbrs[lower[facei]]++;
152  }
153 
154  // Construct cell-upper cell addressing
155  labelList offsets(nCells+1);
156  offsets[0] = 0;
157  forAll(nNbrs, celli)
158  {
159  offsets[celli+1] = offsets[celli]+nNbrs[celli];
160  }
161 
162  nNbrs = offsets;
163 
164  labelList cellToFaces(offsets.last());
165  forAll(upper, facei)
166  {
167  label celli = lower[facei];
168  cellToFaces[nNbrs[celli]++] = facei;
169  }
170 
171  // Sort
172 
173  labelList oldToNew(lower.size());
174 
175  labelList order;
176  labelList nbr;
177 
178  label newFacei = 0;
179 
180  for (label celli = 0; celli < nCells; celli++)
181  {
182  label startOfCell = offsets[celli];
183  label nNbr = offsets[celli+1] - startOfCell;
184 
185  nbr.setSize(nNbr);
186  forAll(nbr, i)
187  {
188  nbr[i] = upper[cellToFaces[offsets[celli]+i]];
189  }
190  sortedOrder(nbr, order);
191 
192  for (const label index : order)
193  {
194  oldToNew[cellToFaces[startOfCell + index]] = newFacei++;
195  }
196  }
197 
198  return oldToNew;
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
203 
204 Foam::lduPrimitiveMesh::lduPrimitiveMesh
205 (
206  const label nCells,
207  labelList& l,
208  labelList& u,
209  const label comm,
210  bool reuse
211 )
212 :
213  lduAddressing(nCells),
214  lowerAddr_(l, reuse),
215  upperAddr_(u, reuse),
216  comm_(comm)
217 {}
218 
219 
221 (
222  lduInterfacePtrsList& interfaces,
223  const lduSchedule& ps
224 )
225 {
226  interfaces_ = interfaces;
227  patchSchedule_ = ps;
228 
229  // Create interfaces
230  primitiveInterfaces_.setSize(interfaces_.size());
231  forAll(interfaces_, i)
232  {
233  if (interfaces_.set(i))
234  {
235  primitiveInterfaces_.set(i, &interfaces_[i]);
236  }
237  }
238 }
239 
240 
241 Foam::lduPrimitiveMesh::lduPrimitiveMesh
242 (
243  const label nCells,
244  labelList& l,
245  labelList& u,
246  PtrList<const lduInterface>& primitiveInterfaces,
247  const lduSchedule& ps,
248  const label comm
249 )
250 :
251  lduAddressing(nCells),
252  lowerAddr_(l, true),
253  upperAddr_(u, true),
254  primitiveInterfaces_(0),
255  patchSchedule_(ps),
256  comm_(comm)
257 {
258  primitiveInterfaces_.transfer(primitiveInterfaces);
259 
260  // Create interfaces
261  interfaces_.setSize(primitiveInterfaces_.size());
262  forAll(primitiveInterfaces_, i)
263  {
264  if (primitiveInterfaces_.set(i))
265  {
266  interfaces_.set(i, &primitiveInterfaces_[i]);
267  }
268  }
269 }
270 
271 
272 Foam::lduPrimitiveMesh::lduPrimitiveMesh
273 (
274  const label comm,
275  const labelList& procAgglomMap,
276 
277  const labelList& procIDs,
278  const lduMesh& myMesh,
279  const PtrList<lduPrimitiveMesh>& otherMeshes,
280 
281  labelList& cellOffsets,
282  labelList& faceOffsets,
284  labelListList& boundaryMap,
285  labelListListList& boundaryFaceMap
286 )
287 :
288  lduAddressing(myMesh.lduAddr().size() + totalSize(otherMeshes)),
289  lowerAddr_(0),
290  upperAddr_(0),
291  interfaces_(0),
292  patchSchedule_(0),
293  comm_(comm)
294 {
295  const label currentComm = myMesh.comm();
296 
297  forAll(otherMeshes, i)
298  {
299  if (otherMeshes[i].comm() != currentComm)
300  {
302  << "Communicator " << otherMeshes[i].comm()
303  << " at index " << i
304  << " differs from that of predecessor "
305  << currentComm
306  << endl; //exit(FatalError);
307  }
308  }
309 
310  const label nMeshes = otherMeshes.size()+1;
311 
312  const label myAgglom = procAgglomMap[UPstream::myProcNo(currentComm)];
313 
315  {
316  Pout<< "I am " << UPstream::myProcNo(currentComm)
317  << " agglomerating into " << myAgglom
318  << " as are " << findIndices(procAgglomMap, myAgglom)
319  << endl;
320  }
321 
322 
323  forAll(procIDs, i)
324  {
325  if (procAgglomMap[procIDs[i]] != procAgglomMap[procIDs[0]])
326  {
328  << "Processor " << procIDs[i]
329  << " agglomerates to " << procAgglomMap[procIDs[i]]
330  << " whereas other processors " << procIDs
331  << " agglomerate to "
332  << labelUIndList(procAgglomMap, procIDs)
333  << exit(FatalError);
334  }
335  }
336 
337 
338  // Cells get added in order.
339  cellOffsets.setSize(nMeshes+1);
340  cellOffsets[0] = 0;
341  for (label procMeshI = 0; procMeshI < nMeshes; procMeshI++)
342  {
343  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
344 
345  cellOffsets[procMeshI+1] =
346  cellOffsets[procMeshI]
347  + procMesh.lduAddr().size();
348  }
349 
350 
351  // Faces initially get added in order, sorted later
352  labelList internalFaceOffsets(nMeshes+1);
353  internalFaceOffsets[0] = 0;
354  for (label procMeshI = 0; procMeshI < nMeshes; procMeshI++)
355  {
356  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
357 
358  internalFaceOffsets[procMeshI+1] =
359  internalFaceOffsets[procMeshI]
360  + procMesh.lduAddr().lowerAddr().size();
361  }
362 
363  // Count how faces get added. Interfaces inbetween get merged.
364 
365  // Merged interfaces: map from two coarse processors back to
366  // - procMeshes
367  // - interface in procMesh
368  // (estimate size from number of patches of mesh0)
369  EdgeMap<labelPairList> mergedMap(2*myMesh.interfaces().size());
370 
371  // Unmerged interfaces: map from two coarse processors back to
372  // - procMeshes
373  // - interface in procMesh
374  EdgeMap<labelPairList> unmergedMap(mergedMap.size());
375 
376  boundaryMap.setSize(nMeshes);
377  boundaryFaceMap.setSize(nMeshes);
378 
379 
380  label nOtherInterfaces = 0;
381  labelList nCoupledFaces(nMeshes, Zero);
382 
383  for (label procMeshI = 0; procMeshI < nMeshes; procMeshI++)
384  {
385  const lduInterfacePtrsList interfaces =
386  mesh(myMesh, otherMeshes, procMeshI).interfaces();
387 
388  // Initialise all boundaries as merged
389  boundaryMap[procMeshI].setSize(interfaces.size(), -1);
390  boundaryFaceMap[procMeshI].setSize(interfaces.size());
391 
392  // Get sorted order of processors
393  forAll(interfaces, intI)
394  {
395  if (interfaces.set(intI))
396  {
397  const lduInterface& ldui = interfaces[intI];
398 
399  if (isA<processorLduInterface>(ldui))
400  {
401  const processorLduInterface& pldui =
402  refCast<const processorLduInterface>(ldui);
403 
404  label agglom0 = procAgglomMap[pldui.myProcNo()];
405  label agglom1 = procAgglomMap[pldui.neighbProcNo()];
406 
407  const edge procEdge(agglom0, agglom1);
408 
409  if (agglom0 != myAgglom && agglom1 != myAgglom)
410  {
412  << "At mesh from processor " << procIDs[procMeshI]
413  << " have interface " << intI
414  << " with myProcNo:" << pldui.myProcNo()
415  << " with neighbProcNo:" << pldui.neighbProcNo()
416  << exit(FatalError);
417  }
418  else if (agglom0 == myAgglom && agglom1 == myAgglom)
419  {
420  // Merged interface
421  if (debug)
422  {
423  Pout<< "merged interface: myProcNo:"
424  << pldui.myProcNo()
425  << " nbr:" << pldui.neighbProcNo()
426  << " size:" << ldui.faceCells().size()
427  << endl;
428  }
429 
430  const label nbrProcMeshI =
431  procIDs.find(pldui.neighbProcNo());
432 
433  if (procMeshI < nbrProcMeshI)
434  {
435  // I am 'master' since get lowest numbered cells
436  nCoupledFaces[procMeshI] += ldui.faceCells().size();
437  }
438 
439  mergedMap(procEdge).append
440  (
441  labelPair(procMeshI, intI)
442  );
443  }
444  else
445  {
446  if (debug)
447  {
448  Pout<< "external interface: myProcNo:"
449  << pldui.myProcNo()
450  << " nbr:" << pldui.neighbProcNo()
451  << " size:" << ldui.faceCells().size()
452  << endl;
453  }
454 
455  unmergedMap(procEdge).append
456  (
457  labelPair(procMeshI, intI)
458  );
459  }
460  }
461  else
462  {
463  // Still external (non proc) interface
465  << "At mesh from processor " << procIDs[procMeshI]
466  << " have interface " << intI
467  << " of unhandled type " << interfaces[intI].type()
468  << exit(FatalError);
469 
470  ++nOtherInterfaces;
471  }
472  }
473  }
474  }
475 
476 
477 
478  if (debug)
479  {
480  Pout<< "Remaining interfaces:" << endl;
481  forAllConstIters(unmergedMap, iter)
482  {
483  Pout<< " agglom procEdge:" << iter.key() << endl;
484  const labelPairList& elems = iter.val();
485  forAll(elems, i)
486  {
487  label procMeshI = elems[i][0];
488  label interfacei = elems[i][1];
489  const lduInterfacePtrsList interfaces =
490  mesh(myMesh, otherMeshes, procMeshI).interfaces();
491 
492  const processorLduInterface& pldui =
493  refCast<const processorLduInterface>
494  (
495  interfaces[interfacei]
496  );
497 
498  Pout<< " proc:" << procIDs[procMeshI]
499  << " interfacei:" << interfacei
500  << " between:" << pldui.myProcNo()
501  << " and:" << pldui.neighbProcNo()
502  << endl;
503  }
504  Pout<< endl;
505  }
506  }
507  if (debug)
508  {
509  Pout<< "Merged interfaces:" << endl;
510  forAllConstIters(mergedMap, iter)
511  {
512  Pout<< " agglom procEdge:" << iter.key() << endl;
513  const labelPairList& elems = iter.val();
514 
515  forAll(elems, i)
516  {
517  label procMeshI = elems[i][0];
518  label interfacei = elems[i][1];
519  const lduInterfacePtrsList interfaces =
520  mesh(myMesh, otherMeshes, procMeshI).interfaces();
521  const processorLduInterface& pldui =
522  refCast<const processorLduInterface>
523  (
524  interfaces[interfacei]
525  );
526 
527  Pout<< " proc:" << procIDs[procMeshI]
528  << " interfacei:" << interfacei
529  << " between:" << pldui.myProcNo()
530  << " and:" << pldui.neighbProcNo()
531  << endl;
532  }
533  Pout<< endl;
534  }
535  }
536 
537 
538  // Adapt faceOffsets for internal interfaces
539  faceOffsets.setSize(nMeshes+1);
540  faceOffsets[0] = 0;
541  faceMap.setSize(nMeshes);
542  for (label procMeshI = 0; procMeshI < nMeshes; procMeshI++)
543  {
544  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
545  label nInternal = procMesh.lduAddr().lowerAddr().size();
546 
547  faceOffsets[procMeshI+1] =
548  faceOffsets[procMeshI]
549  + nInternal
550  + nCoupledFaces[procMeshI];
551 
552  labelList& map = faceMap[procMeshI];
553  map.setSize(nInternal);
554  forAll(map, i)
555  {
556  map[i] = faceOffsets[procMeshI] + i;
557  }
558  }
559 
560 
561  // Combine upper and lower
562  lowerAddr_.setSize(faceOffsets.last(), -1);
563  upperAddr_.setSize(lowerAddr_.size(), -1);
564 
565 
566  // Old internal faces and resolved coupled interfaces
567  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
568 
569  for (label procMeshI = 0; procMeshI < nMeshes; procMeshI++)
570  {
571  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
572 
573  const labelUList& l = procMesh.lduAddr().lowerAddr();
574  const labelUList& u = procMesh.lduAddr().upperAddr();
575 
576  // Add internal faces
577  label allFacei = faceOffsets[procMeshI];
578 
579  forAll(l, facei)
580  {
581  lowerAddr_[allFacei] = cellOffsets[procMeshI]+l[facei];
582  upperAddr_[allFacei] = cellOffsets[procMeshI]+u[facei];
583  allFacei++;
584  }
585 
586 
587  // Add merged interfaces
588  const lduInterfacePtrsList interfaces = procMesh.interfaces();
589 
590  forAll(interfaces, intI)
591  {
592  if (interfaces.set(intI))
593  {
594  const lduInterface& ldui = interfaces[intI];
595 
596  if (isA<processorLduInterface>(ldui))
597  {
598  const processorLduInterface& pldui =
599  refCast<const processorLduInterface>(ldui);
600 
601  // Look up corresponding interfaces
602  label myP = pldui.myProcNo();
603  label nbrP = pldui.neighbProcNo();
604  label nbrProcMeshI = procIDs.find(nbrP);
605 
606  if (procMeshI < nbrProcMeshI)
607  {
608  // I am 'master' since my cell numbers will be lower
609  // since cells get added in procMeshI order.
610 
611  const label agglom0 = procAgglomMap[myP];
612  const label agglom1 = procAgglomMap[nbrP];
613 
614  const auto fnd =
615  mergedMap.cfind(edge(agglom0, agglom1));
616 
617  if (fnd.found())
618  {
619  const labelPairList& elems = fnd();
620 
621  // Find nbrP in elems
622  label nbrIntI = -1;
623  forAll(elems, i)
624  {
625  label proci = elems[i][0];
626  label interfacei = elems[i][1];
627  const lduInterfacePtrsList interfaces =
628  mesh
629  (
630  myMesh,
631  otherMeshes,
632  proci
633  ).interfaces();
634  const processorLduInterface& pldui =
635  refCast<const processorLduInterface>
636  (
637  interfaces[interfacei]
638  );
639 
640  if
641  (
642  elems[i][0] == nbrProcMeshI
643  && pldui.neighbProcNo() == procIDs[procMeshI]
644  )
645  {
646  nbrIntI = elems[i][1];
647  break;
648  }
649  }
650 
651 
652  if (nbrIntI == -1)
653  {
655  << "elems:" << elems << abort(FatalError);
656  }
657 
658 
659  const lduInterfacePtrsList nbrInterfaces = mesh
660  (
661  myMesh,
662  otherMeshes,
663  nbrProcMeshI
664  ).interfaces();
665 
666 
667  const labelUList& faceCells =
668  ldui.faceCells();
669  const labelUList& nbrFaceCells =
670  nbrInterfaces[nbrIntI].faceCells();
671 
672  if (faceCells.size() != nbrFaceCells.size())
673  {
675  << "faceCells:" << faceCells
676  << " nbrFaceCells:" << nbrFaceCells
677  << abort(FatalError);
678  }
679 
680 
681  labelList& bfMap =
682  boundaryFaceMap[procMeshI][intI];
683  labelList& nbrBfMap =
684  boundaryFaceMap[nbrProcMeshI][nbrIntI];
685 
686  bfMap.setSize(faceCells.size());
687  nbrBfMap.setSize(faceCells.size());
688 
689  forAll(faceCells, pfI)
690  {
691  lowerAddr_[allFacei] =
692  cellOffsets[procMeshI]+faceCells[pfI];
693  bfMap[pfI] = allFacei;
694  upperAddr_[allFacei] =
695  cellOffsets[nbrProcMeshI]+nbrFaceCells[pfI];
696  nbrBfMap[pfI] = (-allFacei-1);
697  allFacei++;
698  }
699  }
700  }
701  }
702  }
703  }
704  }
705 
706 
707  // Sort upper-tri order
708  {
709  labelList oldToNew
710  (
711  upperTriOrder
712  (
713  cellOffsets.last(), //nCells
714  lowerAddr_,
715  upperAddr_
716  )
717  );
718 
719  forAll(faceMap, procMeshI)
720  {
721  labelList& map = faceMap[procMeshI];
722  forAll(map, i)
723  {
724  if (map[i] >= 0)
725  {
726  map[i] = oldToNew[map[i]];
727  }
728  else
729  {
730  label allFacei = -map[i]-1;
731  map[i] = -oldToNew[allFacei]-1;
732  }
733  }
734  }
735 
736 
737  inplaceReorder(oldToNew, lowerAddr_);
738  inplaceReorder(oldToNew, upperAddr_);
739 
740  forAll(boundaryFaceMap, proci)
741  {
742  const labelList& bMap = boundaryMap[proci];
743  forAll(bMap, intI)
744  {
745  if (bMap[intI] == -1)
746  {
747  // Merged interface
748  labelList& bfMap = boundaryFaceMap[proci][intI];
749 
750  forAll(bfMap, i)
751  {
752  if (bfMap[i] >= 0)
753  {
754  bfMap[i] = oldToNew[bfMap[i]];
755  }
756  else
757  {
758  label allFacei = -bfMap[i]-1;
759  bfMap[i] = (-oldToNew[allFacei]-1);
760  }
761  }
762  }
763  }
764  }
765  }
766 
767 
768  // Kept interfaces
769  // ~~~~~~~~~~~~~~~
770 
771  interfaces_.setSize(unmergedMap.size() + nOtherInterfaces);
772  primitiveInterfaces_.setSize(interfaces_.size());
773 
774  label allInterfacei = 0;
775 
776  forAllConstIters(unmergedMap, iter)
777  {
778  const labelPairList& elems = iter.val();
779 
780  // Sort processors in increasing order so both sides walk through in
781  // same order.
782  labelPairList procPairs(elems.size());
783  forAll(elems, i)
784  {
785  const labelPair& elem = elems[i];
786  label procMeshI = elem[0];
787  label interfacei = elem[1];
788  const lduInterfacePtrsList interfaces = mesh
789  (
790  myMesh,
791  otherMeshes,
792  procMeshI
793  ).interfaces();
794 
795  const processorLduInterface& pldui =
796  refCast<const processorLduInterface>
797  (
798  interfaces[interfacei]
799  );
800  label myProcNo = pldui.myProcNo();
801  label nbrProcNo = pldui.neighbProcNo();
802  procPairs[i] = labelPair
803  (
804  min(myProcNo, nbrProcNo),
805  max(myProcNo, nbrProcNo)
806  );
807  }
808 
809  labelList order(sortedOrder(procPairs));
810 
811  // Count
812  label n = 0;
813 
814  forAll(order, i)
815  {
816  const labelPair& elem = elems[order[i]];
817  label procMeshI = elem[0];
818  label interfacei = elem[1];
819  const lduInterfacePtrsList interfaces = mesh
820  (
821  myMesh,
822  otherMeshes,
823  procMeshI
824  ).interfaces();
825 
826  n += interfaces[interfacei].faceCells().size();
827  }
828 
829  // Size
830  labelField allFaceCells(n);
831  labelField allFaceRestrictAddressing(n);
832  n = 0;
833 
834  // Fill
835  forAll(order, i)
836  {
837  const labelPair& elem = elems[order[i]];
838  label procMeshI = elem[0];
839  label interfacei = elem[1];
840  const lduInterfacePtrsList interfaces = mesh
841  (
842  myMesh,
843  otherMeshes,
844  procMeshI
845  ).interfaces();
846 
847  boundaryMap[procMeshI][interfacei] = allInterfacei;
848  labelList& bfMap = boundaryFaceMap[procMeshI][interfacei];
849 
850  const labelUList& l = interfaces[interfacei].faceCells();
851  bfMap.setSize(l.size());
852 
853  forAll(l, facei)
854  {
855  allFaceCells[n] = cellOffsets[procMeshI]+l[facei];
856  allFaceRestrictAddressing[n] = n;
857  bfMap[facei] = n;
858  n++;
859  }
860  }
861 
862 
863  // Find out local and remote processor in new communicator
864 
865  label neighbProcNo = -1;
866 
867  // See what the two processors map onto
868 
869  if (iter.key()[0] == myAgglom)
870  {
871  if (iter.key()[1] == myAgglom)
872  {
874  << "problem procEdge:" << iter.key()
875  << exit(FatalError);
876  }
877 
878  neighbProcNo = iter.key()[1];
879  }
880  else
881  {
882  if (iter.key()[1] != myAgglom)
883  {
885  << "problem procEdge:" << iter.key()
886  << exit(FatalError);
887  }
888 
889  neighbProcNo = iter.key()[0];
890  }
891 
892  primitiveInterfaces_.set
893  (
894  allInterfacei,
896  (
897  allInterfacei,
898  interfaces_,
899  allFaceCells,
900  allFaceRestrictAddressing,
901  comm_,
902  myAgglom,
903  neighbProcNo,
904  tensorField(), // forwardT
905  Pstream::msgType() // tag
906  )
907  );
908  interfaces_.set(allInterfacei, &primitiveInterfaces_[allInterfacei]);
909 
910  if (debug)
911  {
912  Pout<< "Created " << interfaces_[allInterfacei].type()
913  << " interface at " << allInterfacei
914  << " comm:" << comm_
915  << " myProcNo:" << myAgglom
916  << " neighbProcNo:" << neighbProcNo
917  << " nFaces:" << allFaceCells.size()
918  << endl;
919  }
920 
921 
922  allInterfacei++;
923  }
924 
925 
926  patchSchedule_ = nonBlockingSchedule<processorGAMGInterface>(interfaces_);
927 
928  if (debug)
929  {
930  checkUpperTriangular(cellOffsets.last(), lowerAddr_, upperAddr_);
931  }
932 }
933 
934 
935 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
936 
938 (
939  const lduMesh& myMesh,
940  const PtrList<lduPrimitiveMesh>& otherMeshes,
941  const label meshI
942 )
943 {
944  return (meshI == 0 ? myMesh : otherMeshes[meshI-1]);
945 }
946 
947 
949 (
950  const label comm,
951  const lduMesh& mesh,
952  const labelList& procIDs,
953  PtrList<lduPrimitiveMesh>& otherMeshes
954 )
955 {
956  // Force calculation of schedule (since does parallel comms)
957  (void)mesh.lduAddr().patchSchedule();
958 
959  if (Pstream::myProcNo(comm) == procIDs[0])
960  {
961  otherMeshes.setSize(procIDs.size()-1);
962 
963  // Slave meshes
964  for (label i = 1; i < procIDs.size(); i++)
965  {
966  //Pout<< "on master :"
967  // << " receiving from slave " << procIDs[i] << endl;
968 
969  IPstream fromSlave
970  (
972  procIDs[i],
973  0, // bufSize
975  comm
976  );
977 
978  label nCells = readLabel(fromSlave);
979  labelList lowerAddr(fromSlave);
980  labelList upperAddr(fromSlave);
981  boolList validInterface(fromSlave);
982 
983 
984  // Construct mesh without interfaces
985  otherMeshes.set
986  (
987  i-1,
988  new lduPrimitiveMesh
989  (
990  nCells,
991  lowerAddr,
992  upperAddr,
993  comm,
994  true // reuse
995  )
996  );
997 
998  // Construct GAMGInterfaces
999  lduInterfacePtrsList newInterfaces(validInterface.size());
1000  forAll(validInterface, intI)
1001  {
1002  if (validInterface[intI])
1003  {
1004  word coupleType(fromSlave);
1005 
1006  newInterfaces.set
1007  (
1008  intI,
1010  (
1011  coupleType,
1012  intI,
1013  otherMeshes[i-1].rawInterfaces(),
1014  fromSlave
1015  ).ptr()
1016  );
1017  }
1018  }
1019 
1020  otherMeshes[i-1].addInterfaces
1021  (
1022  newInterfaces,
1023  nonBlockingSchedule<processorGAMGInterface>
1024  (
1025  newInterfaces
1026  )
1027  );
1028  }
1029  }
1030  else if (procIDs.found(Pstream::myProcNo(comm)))
1031  {
1032  // Send to master
1033 
1034  const lduAddressing& addressing = mesh.lduAddr();
1035  lduInterfacePtrsList interfaces(mesh.interfaces());
1036  boolList validInterface(interfaces.size());
1037  forAll(interfaces, intI)
1038  {
1039  validInterface[intI] = interfaces.set(intI);
1040  }
1041 
1042  OPstream toMaster
1043  (
1045  procIDs[0],
1046  0,
1047  Pstream::msgType(),
1048  comm
1049  );
1050 
1051  toMaster
1052  << addressing.size()
1053  << addressing.lowerAddr()
1054  << addressing.upperAddr()
1055  << validInterface;
1056 
1057  forAll(interfaces, intI)
1058  {
1059  if (interfaces.set(intI))
1060  {
1061  const GAMGInterface& interface = refCast<const GAMGInterface>
1062  (
1063  interfaces[intI]
1064  );
1065 
1066  toMaster << interface.type();
1067  interface.write(toMaster);
1068  }
1069  }
1070  }
1071 }
1072 
1073 
1074 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::lduAddressing
The class contains the addressing required by the lduMatrix: upper, lower and losort.
Definition: lduAddressing.H:114
Foam::UPtrList::size
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:97
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::lduPrimitiveMesh::gather
static void gather(const label comm, const lduMesh &mesh, const labelList &procIDs, PtrList< lduPrimitiveMesh > &otherMeshes)
Gather meshes from other processors onto procIDs[0].
Definition: lduPrimitiveMesh.C:949
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::lduPrimitiveMesh::addInterfaces
void addInterfaces(lduInterfacePtrsList &interfaces, const lduSchedule &ps)
Add interfaces to a mesh constructed without.
Definition: lduPrimitiveMesh.C:221
interface
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Foam::lduInterface
An abstract base class for implicitly-coupled interfaces e.g. processor and cyclic patches.
Definition: lduInterface.H:54
Foam::processorLduInterface
An abstract base class for processor coupled interfaces.
Definition: processorLduInterface.H:53
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::tensorField
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
Definition: primitiveFieldsFwd.H:57
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: PtrList.H:138
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::lduPrimitiveMesh::upperTriOrder
static labelList upperTriOrder(const label nCells, const labelUList &lower, const labelUList &upper)
Calculate upper-triangular order.
Definition: lduPrimitiveMesh.C:132
Foam::lduPrimitiveMesh::checkUpperTriangular
static void checkUpperTriangular(const label size, const labelUList &l, const labelUList &u)
Check if in upper-triangular ordering.
Definition: lduPrimitiveMesh.C:64
Foam::lduAddressing::size
label size() const
Return number of equations.
Definition: lduAddressing.H:171
Foam::lduAddressing::upperAddr
virtual const labelUList & upperAddr() const =0
Return upper addressing.
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
processorLduInterface.H
Foam::lduMesh::interfaces
virtual lduInterfacePtrsList interfaces() const =0
Return a list of pointers for each patch.
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::stringOps::lower
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1186
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< label >
Foam::GAMGInterface
Abstract base class for GAMG agglomerated interfaces.
Definition: GAMGInterface.H:53
Foam::fvMesh::lduAddr
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:685
Foam::procLess::operator()
bool operator()(const label a, const label b)
Definition: lduPrimitiveMesh.C:53
Foam::procLess::procLess
procLess(const labelPairList &lst)
Definition: lduPrimitiveMesh.C:48
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:105
Foam::UPtrList< const lduInterface >
meshes
PtrList< fvMesh > meshes(regionNames.size())
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::lduAddressing::patchSchedule
virtual const lduSchedule & patchSchedule() const =0
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::processorLduInterface::neighbProcNo
virtual int neighbProcNo() const =0
Return neighbour processor number (rank in communicator)
Foam::UPstream::commsTypes::scheduled
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
processorGAMGInterface.H
Foam::procLess
Less operator for pairs of <processor><index>
Definition: lduPrimitiveMesh.C:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::processorGAMGInterface
GAMG agglomerated processor interface.
Definition: processorGAMGInterface.H:52
Foam::lduPrimitiveMesh::mesh
static const lduMesh & mesh(const lduMesh &mesh0, const PtrList< lduPrimitiveMesh > &otherMeshes, const label meshI)
Select either mesh0 (meshI is 0) or otherMeshes[meshI-1].
Definition: lduPrimitiveMesh.C:938
Foam::lduPrimitiveMesh
Simplest contrete lduMesh which stores the addressing needed by lduMatrix.
Definition: lduPrimitiveMesh.H:52
Foam::EdgeMap
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:51
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:541
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::processorLduInterface::myProcNo
virtual int myProcNo() const =0
Return processor number (rank in communicator)
edgeHashes.H
Foam::UPtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: UPtrList.H:170
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:464
Foam::Pair< label >
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::lduAddressing::lowerAddr
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
Foam::readLabel
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:66
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::fvMesh::interfaces
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.H:308
Foam::UList< label >
Foam::GAMGInterface::New
static autoPtr< GAMGInterface > New(const label index, const lduInterfacePtrsList &coarseInterfaces, const lduInterface &fineInterface, const labelField &localRestrictAddressing, const labelField &neighbourRestrictAddressing, const label fineLevelIndex, const label coarseComm)
Return a pointer to a new interface created on freestore given.
Definition: GAMGInterfaceNew.C:37
Foam::lduMesh::comm
virtual label comm() const =0
Return communicator used for parallel communication.
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::findIndices
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
Foam::inplaceReorder
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
Definition: ListOpsTemplates.C:124
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
lduPrimitiveMesh.H
Foam::stringOps::upper
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1202
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::lduMesh::lduAddr
virtual const lduAddressing & lduAddr() const =0
Return ldu addressing.
Foam::lduMesh
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:62
labelPair.H
Foam::lduInterface::faceCells
virtual const labelUList & faceCells() const =0
Return faceCell addressing.
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58