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-2023 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"
32 #include "edgeHashes.H"
33 #include "labelPair.H"
34 #include "processorGAMGInterface.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 (
48  const label size,
49  const labelUList& l,
50  const labelUList& u
51 )
52 {
53  forAll(l, facei)
54  {
55  if (u[facei] < l[facei])
56  {
58  << "Reversed face. Problem at face " << facei
59  << " l:" << l[facei] << " u:" << u[facei]
60  << abort(FatalError);
61  }
62  if (l[facei] < 0 || u[facei] < 0 || u[facei] >= size)
63  {
65  << "Illegal cell label. Problem at face " << facei
66  << " l:" << l[facei] << " u:" << u[facei]
67  << abort(FatalError);
68  }
69  }
70 
71  for (label facei=1; facei < l.size(); ++facei)
72  {
73  if (l[facei-1] > l[facei])
74  {
76  << "Lower not in incremental cell order."
77  << " Problem at face " << facei
78  << " l:" << l[facei] << " u:" << u[facei]
79  << " previous l:" << l[facei-1]
80  << abort(FatalError);
81  }
82  else if (l[facei-1] == l[facei])
83  {
84  // Same cell.
85  if (u[facei-1] > u[facei])
86  {
88  << "Upper not in incremental cell order."
89  << " Problem at face " << facei
90  << " l:" << l[facei] << " u:" << u[facei]
91  << " previous u:" << u[facei-1]
92  << abort(FatalError);
93  }
94  }
95  }
96 }
97 
98 
100 (
101  const lduMesh& mesh,
102  const globalIndex& globalNumbering
103 )
104 {
105  const lduAddressing& addr = mesh.lduAddr();
106  lduInterfacePtrsList interfaces = mesh.interfaces();
107 
108  const labelList globalIndices
109  (
110  identity
111  (
112  addr.size(),
113  globalNumbering.localStart(UPstream::myProcNo(mesh.comm()))
114  )
115  );
116 
117  // Get the interface cells
118  PtrList<labelList> nbrGlobalCells(interfaces.size());
119  {
120  const label startOfRequests = UPstream::nRequests();
121 
122  // Initialise transfer of restrict addressing on the interface
123  forAll(interfaces, inti)
124  {
125  if (interfaces.set(inti))
126  {
127  interfaces[inti].initInternalFieldTransfer
128  (
130  globalIndices
131  );
132  }
133  }
134 
135  UPstream::waitRequests(startOfRequests);
136 
137  forAll(interfaces, inti)
138  {
139  if (interfaces.set(inti))
140  {
141  nbrGlobalCells.set
142  (
143  inti,
144  new labelList
145  (
146  interfaces[inti].internalFieldTransfer
147  (
149  globalIndices
150  )
151  )
152  );
153  }
154  }
155  }
156 
157 
158  // Scan the neighbour list to find out how many times the cell
159  // appears as a neighbour of the face. Done this way to avoid guessing
160  // and resizing list
161  labelList nNbrs(addr.size(), Zero);
162 
163  const labelUList& nbr = addr.upperAddr();
164  const labelUList& own = addr.lowerAddr();
165 
166  {
167  forAll(nbr, facei)
168  {
169  nNbrs[nbr[facei]]++;
170  nNbrs[own[facei]]++;
171  }
172 
173  forAll(interfaces, inti)
174  {
175  if (interfaces.set(inti))
176  {
177  for (const label celli : interfaces[inti].faceCells())
178  {
179  nNbrs[celli]++;
180  }
181  }
182  }
183  }
184 
185 
186  // Create cell-cells addressing
187  labelListList cellCells(addr.size());
188 
189  forAll(cellCells, celli)
190  {
191  cellCells[celli].setSize(nNbrs[celli], -1);
192  }
193 
194  // Reset the list of number of neighbours to zero
195  nNbrs = 0;
196 
197  // Scatter the neighbour faces
198  forAll(nbr, facei)
199  {
200  const label c0 = own[facei];
201  const label c1 = nbr[facei];
202 
203  cellCells[c0][nNbrs[c0]++] = globalIndices[c1];
204  cellCells[c1][nNbrs[c1]++] = globalIndices[c0];
205  }
206  forAll(interfaces, inti)
207  {
208  if (interfaces.set(inti))
209  {
210  const labelUList& faceCells = interfaces[inti].faceCells();
211 
212  forAll(faceCells, facei)
213  {
214  const label c0 = faceCells[facei];
215  cellCells[c0][nNbrs[c0]++] = nbrGlobalCells[inti][facei];
216  }
217  }
218  }
219 
220  return cellCells;
221 }
222 
223 
224 Foam::label Foam::lduPrimitiveMesh::totalSize
225 (
226  const PtrList<lduPrimitiveMesh>& meshes
227 )
228 {
229  label size = 0;
230 
231  for (const lduPrimitiveMesh& msh : meshes)
232  {
233  size += msh.lduAddr().size();
234  }
235  return size;
236 }
237 
238 
240 (
241  const label nCells,
242  const labelUList& lower,
243  const labelUList& upper
244 )
245 {
246  labelList nNbrs(nCells, Zero);
247 
248  // Count number of upper neighbours
249  forAll(lower, facei)
250  {
251  if (upper[facei] < lower[facei])
252  {
254  << "Problem at face:" << facei
255  << " lower:" << lower[facei]
256  << " upper:" << upper[facei]
257  << exit(FatalError);
258  }
259  nNbrs[lower[facei]]++;
260  }
261 
262  // Construct cell-upper cell addressing
263  labelList offsets(nCells+1);
264  offsets[0] = 0;
265  forAll(nNbrs, celli)
266  {
267  offsets[celli+1] = offsets[celli]+nNbrs[celli];
268  }
269 
270  nNbrs = offsets;
271 
272  labelList cellToFaces(offsets.last());
273  forAll(upper, facei)
274  {
275  const label celli = lower[facei];
276  cellToFaces[nNbrs[celli]++] = facei;
277  }
278 
279  // Sort
280 
281  labelList oldToNew(lower.size());
282 
283  DynamicList<label> order;
284  DynamicList<label> nbr;
285 
286  label newFacei = 0;
287 
288  for (label celli = 0; celli < nCells; ++celli)
289  {
290  const label startOfCell = offsets[celli];
291  const label nNbr = offsets[celli+1] - startOfCell;
292 
293  nbr.resize(nNbr);
294  order.resize(nNbr);
295 
296  forAll(nbr, i)
297  {
298  nbr[i] = upper[cellToFaces[offsets[celli]+i]];
299  }
300  sortedOrder(nbr, order);
301 
302  for (const label index : order)
303  {
304  oldToNew[cellToFaces[startOfCell + index]] = newFacei++;
305  }
306  }
307 
308  return oldToNew;
309 }
310 
311 
312 Foam::label Foam::lduPrimitiveMesh::findConnectedInterface
313 (
314  const lduMesh& myMesh,
315  const PtrList<lduPrimitiveMesh>& otherMeshes,
316  const labelPairList& procAndInterfaces,
317 
318  const label nbrProci,
319  const label myRank
320 ) const
321 {
322  // Find mesh, interfacei in procAndInterfaces
323  label nbrInti = -1;
324  for (const auto& procAndInterface : procAndInterfaces)
325  {
326  const label proci = procAndInterface[0];
327 
328  if (proci == nbrProci)
329  {
330  const label interfacei = procAndInterface[1];
331  const lduInterfacePtrsList interfaces =
332  mesh
333  (
334  myMesh,
335  otherMeshes,
336  proci
337  ).interfaces();
338  const processorLduInterface& pldui =
339  refCast<const processorLduInterface>
340  (
341  interfaces[interfacei]
342  );
343 
344  if (pldui.neighbProcNo() == myRank)
345  {
346  nbrInti = procAndInterface[1];
347  break;
348  }
349  }
350  }
351 
352 
353  if (nbrInti == -1)
354  {
356  << "procAndInterfaces:" << procAndInterfaces << abort(FatalError);
357  }
358  return nbrInti;
359 }
360 
361 
362 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
363 
364 Foam::lduPrimitiveMesh::lduPrimitiveMesh
365 (
366  const label nCells,
367  labelList& l,
368  labelList& u,
369  const label comm,
370  bool reuse
371 )
372 :
373  lduAddressing(nCells),
374  lowerAddr_(l, reuse),
375  upperAddr_(u, reuse),
376  comm_(comm)
377 {
378  if (debug && lowerAddr_.size())
379  {
380  if (max(lowerAddr_) >= nCells || min(lowerAddr_) < 0)
381  {
382  FatalErrorInFunction << "Illegal lower addressing."
383  << " nCells:" << nCells
384  << " max(lower):" << max(lowerAddr_)
385  << " min(lower):" << min(lowerAddr_)
386  << exit(FatalError);
387  }
388  }
389  if (debug && upperAddr_.size())
390  {
391  if (max(upperAddr_) >= nCells || min(upperAddr_) < 0)
392  {
393  FatalErrorInFunction << "Illegal upper addressing."
394  << " nCells:" << nCells
395  << " max(upper):" << max(upperAddr_)
396  << " min(upper):" << min(upperAddr_)
398  }
399  }
400 }
401 
402 
404 (
405  lduInterfacePtrsList& interfaces,
406  const lduSchedule& ps
407 )
408 {
409  interfaces_ = interfaces;
410  patchSchedule_ = ps;
411 
412  // Create interfaces
413  primitiveInterfaces_.setSize(interfaces_.size());
414  forAll(interfaces_, i)
415  {
416  if (interfaces_.set(i))
417  {
418  primitiveInterfaces_.set(i, &interfaces_[i]);
419  }
420  }
421 }
422 
423 
424 Foam::lduPrimitiveMesh::lduPrimitiveMesh
425 (
426  const label nCells
427 )
428 :
429  lduAddressing(nCells),
430  comm_(0)
431 {}
432 
433 
434 Foam::lduPrimitiveMesh::lduPrimitiveMesh
435 (
436  const label nCells,
437  labelList& l,
438  labelList& u,
439  PtrList<const lduInterface>& primitiveInterfaces,
440  const lduSchedule& ps,
441  const label comm
442 )
443 :
444  lduAddressing(nCells),
445  lowerAddr_(l, true),
446  upperAddr_(u, true),
447  primitiveInterfaces_(),
448  patchSchedule_(ps),
449  comm_(comm)
450 {
451  if (debug && lowerAddr_.size())
452  {
453  if (max(lowerAddr_) >= nCells || min(lowerAddr_) < 0)
454  {
455  FatalErrorInFunction << "Illegal lower addressing."
456  << " nCells:" << nCells
457  << " max(lower):" << max(lowerAddr_)
458  << " min(lower):" << min(lowerAddr_)
459  << exit(FatalError);
460  }
461  }
462  if (debug && upperAddr_.size())
463  {
464  if (max(upperAddr_) >= nCells || min(upperAddr_) < 0)
465  {
466  FatalErrorInFunction << "Illegal upper addressing."
467  << " nCells:" << nCells
468  << " max(upper):" << max(upperAddr_)
469  << " min(upper):" << min(upperAddr_)
470  << exit(FatalError);
471  }
472  }
473 
474 
475  primitiveInterfaces_.transfer(primitiveInterfaces);
476 
477  // Create interfaces
478  interfaces_.setSize(primitiveInterfaces_.size());
479  forAll(primitiveInterfaces_, i)
480  {
481  if (primitiveInterfaces_.set(i))
482  {
483  interfaces_.set(i, &primitiveInterfaces_[i]);
484  }
485  }
486 }
487 
488 
489 Foam::lduPrimitiveMesh::lduPrimitiveMesh
490 (
491  const label comm,
492  const labelList& procAgglomMap,
493 
494  const labelList& procIDs,
495  const lduMesh& myMesh,
496  const PtrList<lduPrimitiveMesh>& otherMeshes,
497 
498  labelList& cellOffsets,
499  labelList& faceOffsets,
501  labelListList& boundaryMap,
502  labelListListList& boundaryFaceMap
503 )
504 :
505  lduAddressing(myMesh.lduAddr().size() + totalSize(otherMeshes)),
506  comm_(comm)
507 {
508  const label currentComm = myMesh.comm();
509 
510  forAll(otherMeshes, i)
511  {
512  if (otherMeshes[i].comm() != currentComm)
513  {
515  << "Communicator " << otherMeshes[i].comm()
516  << " at index " << i
517  << " differs from that of predecessor "
518  << currentComm
519  << endl; //exit(FatalError);
520  }
521  }
522 
523  const label nMeshes = otherMeshes.size()+1;
524 
525  const label myAgglom = procAgglomMap[UPstream::myProcNo(currentComm)];
527  {
528  Pout<< "I am " << UPstream::myProcNo(currentComm)
529  << " agglomerating into " << myAgglom
530  << " as are " << findIndices(procAgglomMap, myAgglom)
531  << endl;
532  }
533 
534  const lduInterfacePtrsList myInterfaces = myMesh.interfaces();
535 
536  forAll(procIDs, i)
537  {
538  if (procAgglomMap[procIDs[i]] != procAgglomMap[procIDs[0]])
539  {
541  << "Processor " << procIDs[i]
542  << " agglomerates to " << procAgglomMap[procIDs[i]]
543  << " whereas other processors " << procIDs
544  << " agglomerate to "
545  << labelUIndList(procAgglomMap, procIDs)
546  << exit(FatalError);
547  }
548  }
549 
550 
551  // Cells get added in order.
552  cellOffsets.setSize(nMeshes+1);
553  cellOffsets[0] = 0;
554  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
555  {
556  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
557 
558  cellOffsets[procMeshI+1] =
559  cellOffsets[procMeshI]
560  + procMesh.lduAddr().size();
561  }
562 
563 
564  // Faces initially get added in order, sorted later
565  labelList internalFaceOffsets(nMeshes+1);
566  internalFaceOffsets[0] = 0;
567  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
568  {
569  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
570 
571  internalFaceOffsets[procMeshI+1] =
572  internalFaceOffsets[procMeshI]
573  + procMesh.lduAddr().lowerAddr().size();
574  }
575 
576  // Count how faces get added. Proc interfaces inbetween get merged.
577 
578  // Merged proc interfaces: map from two coarse processors back to
579  // - procMeshes
580  // - interface in procMesh
581  // (estimate size from number of patches of mesh0)
582  EdgeMap<labelPairList> mergedMap(2*myInterfaces.size());
583 
584  // Unmerged proc interfaces: map from two coarse processors back to
585  // - procMeshes
586  // - interface in procMesh
587  EdgeMap<labelPairList> unmergedMap(mergedMap.size());
588 
589  // (unmerged) global interfaces. These are present on all processors
590  // in the same order (and keep the order in the merged mesh)
591  List<DynamicList<label>> procToGlobal(nMeshes);
592 
593 
594  boundaryMap.setSize(nMeshes);
595  boundaryFaceMap.setSize(nMeshes);
596 
597 
598  labelList nCoupledFaces(nMeshes, Zero);
599 
600  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
601  {
603  mesh(myMesh, otherMeshes, procMeshI).interfaces();
604 
605  // Initialise all boundaries as merged
606  boundaryMap[procMeshI].setSize(interfaces.size(), -1);
607  boundaryFaceMap[procMeshI].setSize(interfaces.size());
608 
609  // Get sorted order of processors
610  forAll(interfaces, intI)
611  {
612  if (interfaces.set(intI))
613  {
614  const lduInterface& ldui = interfaces[intI];
615 
616  if (isA<processorLduInterface>(ldui))
617  {
618  if (isA<processorCyclicGAMGInterface>(ldui))
619  {
621  << "At mesh from processor " << procIDs[procMeshI]
622  << " have interface " << intI
623  << " of unhandled type " << ldui.type()
624  << ". Adapt decomposition to avoid these"
625  << exit(FatalError);
626  }
627 
628  const processorLduInterface& pldui =
629  refCast<const processorLduInterface>(ldui);
630 
631  label agglom0 = procAgglomMap[pldui.myProcNo()];
632  label agglom1 = procAgglomMap[pldui.neighbProcNo()];
633 
634  const edge procEdge(agglom0, agglom1);
635 
636  if (agglom0 != myAgglom && agglom1 != myAgglom)
637  {
639  << "At mesh from processor " << procIDs[procMeshI]
640  << " have interface " << intI
641  << " with myProcNo:" << pldui.myProcNo()
642  << " with neighbProcNo:" << pldui.neighbProcNo()
643  << exit(FatalError);
644  }
645  else if (agglom0 == myAgglom && agglom1 == myAgglom)
646  {
647  // Merged interface
648  if (debug)
649  {
650  Pout<< "merged proc interface: myProcNo:"
651  << pldui.myProcNo()
652  << " nbr:" << pldui.neighbProcNo()
653  << " size:" << ldui.faceCells().size()
654  << endl;
655  }
656 
657  const label nbrProcMeshI =
658  procIDs.find(pldui.neighbProcNo());
659 
660  if (procMeshI < nbrProcMeshI)
661  {
662  // I am 'master' since get lowest numbered cells
663  nCoupledFaces[procMeshI] += ldui.faceCells().size();
664  }
665 
666  mergedMap(procEdge).append
667  (
668  labelPair(procMeshI, intI)
669  );
670  }
671  else
672  {
673  if (debug)
674  {
675  Pout<< "external proc interface: myProcNo:"
676  << pldui.myProcNo()
677  << " nbr:" << pldui.neighbProcNo()
678  << " size:" << ldui.faceCells().size()
679  << endl;
680  }
681 
682  unmergedMap(procEdge).append
683  (
684  labelPair(procMeshI, intI)
685  );
686  }
687  }
688  else
689  {
690  // Still external (non proc) interface
691  procToGlobal[procMeshI].append(intI);
692  }
693  }
694  }
695  }
696 
697 
698  // Check that all processors have any global patches in the same order
699  if (debug)
700  {
701  const auto& global0 = procToGlobal[0];
702  for (label procMeshI = 1; procMeshI < nMeshes; ++procMeshI)
703  {
704  const auto& global = procToGlobal[procMeshI];
705  if (global != global0)
706  {
708  << "At mesh from processor " << procIDs[procMeshI]
709  << " have global interfaces " << global
710  << " which differ from those on processor "
711  << procIDs[procMeshI]
712  << " : " << global0
713  << exit(FatalError);
714  }
715  }
716  }
717 
718 
719  if (debug)
720  {
721  Pout<< "Global interfaces:" << endl;
722  const auto& global0 = procToGlobal[0];
723  for (const label intI : global0)
724  {
725  Pout<< " interfacei:" << intI
726  << " type:" << myInterfaces[intI].type()
727  << endl;
728 
729  }
730  Pout<< endl;
731 
732  Pout<< "Remaining interfaces:" << endl;
733  for (const auto& iter : unmergedMap.csorted())
734  {
735  Pout<< " agglom procEdge:" << iter.key() << endl;
736  const labelPairList& elems = iter.val();
737  forAll(elems, i)
738  {
739  label procMeshI = elems[i][0];
740  label interfacei = elems[i][1];
742  mesh(myMesh, otherMeshes, procMeshI).interfaces();
743 
744  const processorLduInterface& pldui =
745  refCast<const processorLduInterface>
746  (
747  interfaces[interfacei]
748  );
749 
750  Pout<< " proc:" << procIDs[procMeshI]
751  << " interfacei:" << interfacei
752  << " between:" << pldui.myProcNo()
753  << " and:" << pldui.neighbProcNo()
754  << " localsize:"
755  << interfaces[interfacei].faceCells().size()
756  << endl;
757  }
758  Pout<< endl;
759  }
760  Pout<< endl;
761 
762  Pout<< "Merged interfaces:" << endl;
763  //forAllConstIters(mergedMap, iter)
764  for (const auto& iter : mergedMap.csorted())
765  {
766  Pout<< " agglom procEdge:" << iter.key() << endl;
767  const labelPairList& elems = iter.val();
768 
769  forAll(elems, i)
770  {
771  label procMeshI = elems[i][0];
772  label interfacei = elems[i][1];
774  mesh(myMesh, otherMeshes, procMeshI).interfaces();
775  const processorLduInterface& pldui =
776  refCast<const processorLduInterface>
777  (
778  interfaces[interfacei]
779  );
780 
781  Pout<< " proc:" << procIDs[procMeshI]
782  << " interfacei:" << interfacei
783  << " between:" << pldui.myProcNo()
784  << " and:" << pldui.neighbProcNo()
785  << " localsize:"
786  << interfaces[interfacei].faceCells().size()
787  << endl;
788  }
789  Pout<< endl;
790  }
791  Pout<< endl;
792  }
793 
794 
795  // Adapt faceOffsets for internal interfaces
796  faceOffsets.setSize(nMeshes+1);
797  faceOffsets[0] = 0;
798  faceMap.setSize(nMeshes);
799  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
800  {
801  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
802  label nInternal = procMesh.lduAddr().lowerAddr().size();
803 
804  faceOffsets[procMeshI+1] =
805  faceOffsets[procMeshI]
806  + nInternal
807  + nCoupledFaces[procMeshI];
808 
809  labelList& map = faceMap[procMeshI];
810  map.setSize(nInternal);
811  forAll(map, i)
812  {
813  map[i] = faceOffsets[procMeshI] + i;
814  }
815  }
816 
817 
818  // Combine upper and lower
819  lowerAddr_.setSize(faceOffsets.last(), -1);
820  upperAddr_.setSize(lowerAddr_.size(), -1);
821 
822 
823  // Old internal faces and resolved coupled interfaces
824  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
825 
826  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
827  {
828  const lduMesh& procMesh = mesh(myMesh, otherMeshes, procMeshI);
829 
830  const labelUList& l = procMesh.lduAddr().lowerAddr();
831  const labelUList& u = procMesh.lduAddr().upperAddr();
832 
833  // Add internal faces
834  label allFacei = faceOffsets[procMeshI];
835 
836  forAll(l, facei)
837  {
838  lowerAddr_[allFacei] = cellOffsets[procMeshI]+l[facei];
839  upperAddr_[allFacei] = cellOffsets[procMeshI]+u[facei];
840  allFacei++;
841  }
842 
843 
844  // Add merged interfaces
845  const lduInterfacePtrsList interfaces = procMesh.interfaces();
846 
847  forAll(interfaces, intI)
848  {
849  if (interfaces.set(intI))
850  {
851  const lduInterface& ldui = interfaces[intI];
852 
853  if (isA<processorLduInterface>(ldui))
854  {
855  const processorLduInterface& pldui =
856  refCast<const processorLduInterface>(ldui);
857 
858  // Look up corresponding interfaces
859  label myP = pldui.myProcNo();
860  label nbrP = pldui.neighbProcNo();
861  label nbrProcMeshI = procIDs.find(nbrP);
862 
863  if (procMeshI < nbrProcMeshI)
864  {
865  // I am 'master' since my cell numbers will be lower
866  // since cells get added in procMeshI order.
867 
868  const label agglom0 = procAgglomMap[myP];
869  const label agglom1 = procAgglomMap[nbrP];
870 
871  const auto fnd =
872  mergedMap.cfind(edge(agglom0, agglom1));
873 
874  if (fnd.good())
875  {
876  const labelPairList& elems = fnd();
877 
878  const label nbrIntI = findConnectedInterface
879  (
880  myMesh,
881  otherMeshes,
882  elems,
883  nbrProcMeshI,
884  procIDs[procMeshI]
885  );
886 
887  const lduInterfacePtrsList nbrInterfaces = mesh
888  (
889  myMesh,
890  otherMeshes,
891  nbrProcMeshI
892  ).interfaces();
893 
894 
895  const labelUList& faceCells =
896  ldui.faceCells();
897  const labelUList& nbrFaceCells =
898  nbrInterfaces[nbrIntI].faceCells();
899 
900  if (faceCells.size() != nbrFaceCells.size())
901  {
903  << "faceCells:" << faceCells
904  << " nbrFaceCells:" << nbrFaceCells
905  << abort(FatalError);
906  }
907 
908 
909  labelList& bfMap =
910  boundaryFaceMap[procMeshI][intI];
911  labelList& nbrBfMap =
912  boundaryFaceMap[nbrProcMeshI][nbrIntI];
913 
914  bfMap.setSize(faceCells.size());
915  nbrBfMap.setSize(faceCells.size());
916 
917  forAll(faceCells, pfI)
918  {
919  lowerAddr_[allFacei] =
920  cellOffsets[procMeshI]+faceCells[pfI];
921  bfMap[pfI] = allFacei;
922  upperAddr_[allFacei] =
923  cellOffsets[nbrProcMeshI]+nbrFaceCells[pfI];
924  nbrBfMap[pfI] = (-allFacei-1);
925  allFacei++;
926  }
927  }
928  }
929  }
930  }
931  }
932  }
933 
934 
935  // Sort upper-tri order
936  {
937  labelList oldToNew
938  (
940  (
941  cellOffsets.last(), //nCells
942  lowerAddr_,
943  upperAddr_
944  )
945  );
946 
947  forAll(faceMap, procMeshI)
948  {
949  labelList& map = faceMap[procMeshI];
950  forAll(map, i)
951  {
952  if (map[i] >= 0)
953  {
954  map[i] = oldToNew[map[i]];
955  }
956  else
957  {
958  label allFacei = -map[i]-1;
959  map[i] = -oldToNew[allFacei]-1;
960  }
961  }
962  }
963 
964 
965  inplaceReorder(oldToNew, lowerAddr_);
966  inplaceReorder(oldToNew, upperAddr_);
967 
968  forAll(boundaryFaceMap, proci)
969  {
970  const labelList& bMap = boundaryMap[proci];
971  forAll(bMap, intI)
972  {
973  if (bMap[intI] == -1)
974  {
975  // Merged interface
976  labelList& bfMap = boundaryFaceMap[proci][intI];
977 
978  forAll(bfMap, i)
979  {
980  if (bfMap[i] >= 0)
981  {
982  bfMap[i] = oldToNew[bfMap[i]];
983  }
984  else
985  {
986  label allFacei = -bfMap[i]-1;
987  bfMap[i] = (-oldToNew[allFacei]-1);
988  }
989  }
990  }
991  }
992  }
993  }
994 
995 
996  // Kept interfaces
997  // ~~~~~~~~~~~~~~~
998 
999  interfaces_.setSize(unmergedMap.size() + procToGlobal[0].size());
1000  primitiveInterfaces_.setSize(interfaces_.size());
1001 
1002  label allInterfacei = 0;
1003 
1004 
1005  // Global interfaces
1006  // ~~~~~~~~~~~~~~~~~
1007  // (e.g. cyclicAMI)
1008 
1009  {
1010  const auto& global0 = procToGlobal[0];
1011 
1012  for (const label interfacei : global0)
1013  {
1014  //Pout<< " interfacei:" << interfacei
1015  // << " type:" << interfaces[interfacei].type()
1016  // << endl;
1017 
1018  // Just add all individual face-cells in processor order
1019 
1020  // Count
1021  label n = 0;
1022  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
1023  {
1024  const auto& procMesh = mesh(myMesh, otherMeshes, procMeshI);
1025  n += procMesh.interfaces()[interfacei].faceCells().size();
1026  }
1027 
1028  // Size
1029  labelField allFaceCells(n);
1030  labelField allFaceRestrictAddressing(n);
1031  labelList faceOffsets(nMeshes+1, 0);
1032  lduInterfacePtrsList allProcInterfaces(nMeshes);
1033  n = 0;
1034 
1035  // Fill
1036  for (label procMeshI = 0; procMeshI < nMeshes; ++procMeshI)
1037  {
1038  faceOffsets[procMeshI] = n;
1039 
1040  const auto& procMesh = mesh(myMesh, otherMeshes, procMeshI);
1041  const lduInterfacePtrsList interfaces = procMesh.interfaces();
1042 
1043  allProcInterfaces.set(procMeshI, &interfaces[interfacei]);
1044  boundaryMap[procMeshI][interfacei] = allInterfacei;
1045  labelList& bfMap = boundaryFaceMap[procMeshI][interfacei];
1046 
1047  const labelUList& l = interfaces[interfacei].faceCells();
1048  bfMap.setSize(l.size());
1049 
1050  forAll(l, facei)
1051  {
1052  allFaceCells[n] = cellOffsets[procMeshI]+l[facei];
1053  allFaceRestrictAddressing[n] = n;
1054  bfMap[facei] = n;
1055  n++;
1056  }
1057  }
1058 
1059  // For convenience populate last element with size
1060  faceOffsets.last() = n;
1061 
1062 
1063  const auto& myFineInterface =
1064  refCast<const GAMGInterface>(myInterfaces[interfacei]);
1065 
1066  autoPtr<GAMGInterface> ppPtr
1067  (
1068  myFineInterface.clone
1069  (
1070  allInterfacei,
1071  interfaces_,
1072  global0,
1073  allFaceCells,
1074  allFaceRestrictAddressing,
1075  faceOffsets,
1076  allProcInterfaces,
1077  comm_,
1078  myAgglom,
1079  procAgglomMap
1080  )
1081  );
1082 
1083  primitiveInterfaces_.set
1084  (
1085  allInterfacei,
1086  ppPtr.ptr()
1087  );
1088 
1089  interfaces_.set(allInterfacei, &primitiveInterfaces_[allInterfacei]);
1090 
1091  if (debug)
1092  {
1093  Pout<< "Created " << interfaces_[allInterfacei].type()
1094  << " interface at " << allInterfacei
1095  << " comm:" << comm_
1096  << " myProcNo:" << myAgglom
1097  << " nFaces:" << allFaceCells.size()
1098  << endl;
1099  }
1100 
1101  allInterfacei++;
1102  }
1103  }
1104 
1105 
1106  // Unmerged processor interfaces
1107  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1108 
1109  //forAllConstIters(unmergedMap, iter)
1110  for (const auto& iter : unmergedMap.csorted())
1111  {
1112  const labelPairList& elems = iter.val();
1113 
1114  // Sort processors in increasing order so both sides walk through in
1115  // same order.
1116  labelPairList procPairs(elems.size());
1117  forAll(elems, i)
1118  {
1119  const labelPair& elem = elems[i];
1120  label procMeshI = elem[0];
1121  label interfacei = elem[1];
1123  (
1124  myMesh,
1125  otherMeshes,
1126  procMeshI
1127  ).interfaces();
1128 
1129  const processorLduInterface& pldui =
1130  refCast<const processorLduInterface>
1131  (
1132  interfaces[interfacei]
1133  );
1134  label myProcNo = pldui.myProcNo();
1135  label nbrProcNo = pldui.neighbProcNo();
1136  procPairs[i] = labelPair
1137  (
1138  min(myProcNo, nbrProcNo),
1139  max(myProcNo, nbrProcNo)
1140  );
1141  }
1142 
1143  labelList order(sortedOrder(procPairs));
1144 
1145  // Count
1146  label n = 0;
1147 
1148  forAll(order, i)
1149  {
1150  const labelPair& elem = elems[order[i]];
1151  label procMeshI = elem[0];
1152  label interfacei = elem[1];
1154  (
1155  myMesh,
1156  otherMeshes,
1157  procMeshI
1158  ).interfaces();
1159 
1160  n += interfaces[interfacei].faceCells().size();
1161  }
1162 
1163  // Size
1164  labelField allFaceCells(n);
1165  labelField allFaceRestrictAddressing(n);
1166  n = 0;
1167 
1168  // Fill
1169  forAll(order, i)
1170  {
1171  const labelPair& elem = elems[order[i]];
1172  label procMeshI = elem[0];
1173  label interfacei = elem[1];
1175  (
1176  myMesh,
1177  otherMeshes,
1178  procMeshI
1179  ).interfaces();
1180 
1181  boundaryMap[procMeshI][interfacei] = allInterfacei;
1182  labelList& bfMap = boundaryFaceMap[procMeshI][interfacei];
1183 
1184  const labelUList& l = interfaces[interfacei].faceCells();
1185  bfMap.setSize(l.size());
1186 
1187  forAll(l, facei)
1188  {
1189  allFaceCells[n] = cellOffsets[procMeshI]+l[facei];
1190  allFaceRestrictAddressing[n] = n;
1191  bfMap[facei] = n;
1192  n++;
1193  }
1194  }
1195 
1196 
1197  // Find out local and remote processor in new communicator
1198 
1199  label neighbProcNo = -1;
1200 
1201  // See what the two processors map onto
1202 
1203  if (iter.key()[0] == myAgglom)
1204  {
1205  if (iter.key()[1] == myAgglom)
1206  {
1208  << "problem procEdge:" << iter.key()
1209  << exit(FatalError);
1210  }
1211 
1212  neighbProcNo = iter.key()[1];
1213  }
1214  else
1215  {
1216  if (iter.key()[1] != myAgglom)
1217  {
1219  << "problem procEdge:" << iter.key()
1220  << exit(FatalError);
1221  }
1222 
1223  neighbProcNo = iter.key()[0];
1224  }
1225 
1226  primitiveInterfaces_.set
1227  (
1228  allInterfacei,
1229  new processorGAMGInterface
1230  (
1231  allInterfacei,
1232  interfaces_,
1233  allFaceCells,
1234  allFaceRestrictAddressing,
1235  comm_,
1236  myAgglom,
1237  neighbProcNo,
1238  tensorField(), // forwardT
1239  Pstream::msgType() // tag
1240  )
1241  );
1242  interfaces_.set(allInterfacei, &primitiveInterfaces_[allInterfacei]);
1243 
1244  if (debug)
1245  {
1246  Pout<< "Created " << interfaces_[allInterfacei].type()
1247  << " interface at " << allInterfacei
1248  << " comm:" << comm_
1249  << " myProcNo:" << myAgglom
1250  << " neighbProcNo:" << neighbProcNo
1251  << " nFaces:" << allFaceCells.size()
1252  << endl;
1253  }
1254 
1255 
1256  allInterfacei++;
1257  }
1258 
1259 
1260  if (allInterfacei != interfaces_.size())
1261  {
1262  FatalErrorInFunction << "allInterfacei:" << allInterfacei
1263  << " interfaces_:" << interfaces_.size() << exit(FatalError);
1264  }
1265 
1266 
1267  if (debug)
1268  {
1269  Pout<< endl
1270  << "Created new lduPrimitiveMesh:" << nl
1271  << " cells:" << this->lduAddr().size() << nl
1272  << " internal face lower:"
1273  << this->lduAddr().lowerAddr().size() << nl
1274  << " internal faces upper:"
1275  << this->lduAddr().upperAddr().size()
1276  << endl;
1277  forAll(interfaces_, i)
1278  {
1279  if (interfaces_.set(i))
1280  {
1281  Pout<< " interface:" << i << " type:" << interfaces_[i].type()
1282  << nl
1283  << " faceCells:" << flatOutput(interfaces_[i].faceCells())
1284  << endl;
1285  }
1286  }
1287 
1288 
1289  Pout<< "Original input meshes:" << endl;
1290  forAll(boundaryMap, procMeshI)
1291  {
1292  const auto& procMesh = mesh(myMesh, otherMeshes, procMeshI);
1293  const lduInterfacePtrsList interfaces = procMesh.interfaces();
1294 
1295  Pout<< " proc:" << procMeshI
1296  << " interfaces:" << interfaces.size() << endl;
1297  forAll(interfaces, inti)
1298  {
1299  if (interfaces.set(inti))
1300  {
1301  Pout<< " int:" << inti
1302  << " type:" << interfaces[inti].type()
1303  << " size:" << interfaces[inti].faceCells().size()
1304  << " maps to:" << boundaryMap[procMeshI][inti]
1305  << endl;
1306  }
1307  }
1308  }
1309  }
1310 
1311 
1312  patchSchedule_ = nonBlockingSchedule<processorGAMGInterface>(interfaces_);
1313 
1314  if (debug)
1315  {
1316  checkUpperTriangular(cellOffsets.last(), lowerAddr_, upperAddr_);
1317  }
1318 }
1319 
1320 
1321 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1322 
1324 (
1325  const lduMesh& myMesh,
1326  const PtrList<lduPrimitiveMesh>& otherMeshes,
1327  const label meshI
1329 {
1330  return (meshI == 0 ? myMesh : otherMeshes[meshI-1]);
1331 }
1332 
1333 
1335 (
1336  const label comm,
1337  const lduMesh& mesh,
1338  PtrList<lduPrimitiveMesh>& otherMeshes
1339 )
1340 {
1341  // Force calculation of schedule (since does parallel comms)
1342  (void)mesh.lduAddr().patchSchedule();
1343 
1344  // Use PstreamBuffers
1345  PstreamBuffers pBufs
1346  (
1349  comm
1350  );
1351 
1352  // Send to master
1353  if (!Pstream::master(comm))
1354  {
1355  const lduAddressing& addressing = mesh.lduAddr();
1356 
1357  lduInterfacePtrsList interfaces(mesh.interfaces());
1358  boolList validInterface(interfaces.size());
1359  forAll(interfaces, intI)
1360  {
1361  validInterface[intI] = interfaces.set(intI);
1362  }
1363 
1364  UOPstream toMaster(UPstream::masterNo(), pBufs);
1365 
1366  toMaster
1367  << addressing.size()
1368  << addressing.lowerAddr()
1369  << addressing.upperAddr()
1370  << validInterface;
1371 
1372  forAll(interfaces, intI)
1373  {
1374  if (interfaces.set(intI))
1375  {
1376  const GAMGInterface& interface = refCast<const GAMGInterface>
1377  (
1378  interfaces[intI]
1379  );
1380 
1381  toMaster << interface.type();
1382  interface.write(toMaster);
1383  }
1384  }
1385  }
1386 
1387  // Wait for finish
1388  pBufs.finishedGathers();
1389 
1390  // Consume
1391  if (Pstream::master(comm))
1392  {
1393  const label nProcs = UPstream::nProcs(comm);
1394 
1395  // Master.
1396  otherMeshes.setSize(nProcs-1);
1397 
1398  for (const int proci : UPstream::subProcs(comm))
1399  {
1400  UIPstream fromProc(proci, pBufs);
1401 
1402  const label nCells = readLabel(fromProc);
1403  labelList lowerAddr(fromProc);
1404  labelList upperAddr(fromProc);
1405  const boolList validInterface(fromProc);
1406 
1407 
1408  // Construct mesh without interfaces
1409  otherMeshes.set
1410  (
1411  proci-1,
1412  new lduPrimitiveMesh
1413  (
1414  nCells,
1415  lowerAddr,
1416  upperAddr,
1417  mesh.comm(),
1418  true // reuse
1419  )
1420  );
1421 
1422  // Construct GAMGInterfaces
1423  lduInterfacePtrsList newInterfaces(validInterface.size());
1424  forAll(validInterface, intI)
1425  {
1426  if (validInterface[intI])
1427  {
1428  word coupleType(fromProc);
1429 
1430  newInterfaces.set
1431  (
1432  intI,
1434  (
1435  coupleType,
1436  intI,
1437  otherMeshes[proci-1].rawInterfaces(),
1438  fromProc
1439  ).ptr()
1440  );
1441  }
1442  }
1443 
1444  otherMeshes[proci-1].addInterfaces
1445  (
1446  newInterfaces,
1447  nonBlockingSchedule<processorGAMGInterface>
1448  (
1449  newInterfaces
1450  )
1451  );
1452  }
1453  }
1454 }
1455 
1456 
1457 // ************************************************************************* //
virtual lduInterfacePtrsList interfaces() const =0
Return a list of pointers for each patch with only those pointing to interfaces being set...
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:48
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void addInterfaces(lduInterfacePtrsList &interfaces, const lduSchedule &ps)
Add interfaces to a mesh constructed without.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
label readLabel(const char *buf)
Parse entire buffer as a label, skipping leading/trailing whitespace.
Definition: label.H:63
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1187
PtrList< const lduInterface > & primitiveInterfaces()
Return a non-const list of primitive interfaces.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
List< lduScheduleEntry > lduSchedule
A List of lduSchedule entries.
Definition: lduSchedule.H:46
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1538
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
virtual label comm() const
Return communicator used for parallel communication.
Definition: fvMesh.H:415
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:741
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1059
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual label comm() const
Return communicator used for parallel communication.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
Simplest concrete lduMesh that stores the addressing needed by lduMatrix.
List< labelListList > labelListListList
List of labelListList.
Definition: labelList.H:41
static const lduMesh & mesh(const lduMesh &mesh0, const PtrList< lduPrimitiveMesh > &otherMeshes, const label meshI)
Select either mesh0 (meshI is 0) or otherMeshes[meshI-1].
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
static void gather(const label agglomComm, const lduMesh &mesh, PtrList< lduPrimitiveMesh > &otherMeshes)
Gather meshes from other processors using agglomComm.
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
defineTypeNameAndDebug(combustionModel, 0)
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:337
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1171
#define WarningInFunction
Report a warning using Foam::Warning.
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
Definition: fvMesh.C:758
static void checkUpperTriangular(const label size, const labelUList &l, const labelUList &u)
Check if in upper-triangular ordering.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static labelListList globalCellCells(const lduMesh &mesh, const globalIndex &globalNumbering)
Calculate global cell-cells.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
"nonBlocking" : (MPI_Isend, MPI_Irecv)
label n
The class contains the addressing required by the lduMatrix: upper, lower and losort.
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:1185
List< label > labelList
A List of labels.
Definition: List.H:62
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:840
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch with only those pointing to interfaces being set...
List< bool > boolList
A List of bools.
Definition: List.H:60
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static labelList upperTriOrder(const label nCells, const labelUList &lower, const labelUList &upper)
Calculate upper-triangular order.
Namespace for OpenFOAM.
label size() const
Return number of equations.
virtual const lduSchedule & patchSchedule() const =0
Return patch field evaluation schedule.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrList.H:366
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127