domainDecompositionNonConformal.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM 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 OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "domainDecomposition.H"
27 #include "cyclicFvPatch.H"
28 #include "processorCyclicFvPatch.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
41 (
42  const labelPair& procs,
43  const fvPatch& fvp,
44  const fvPatch& nbrFvp,
45  const labelUList& polyFacesPf,
46  const labelUList& nbrPolyFacesPf
47 )
48 {
49  if (polyFacesPf.size() != nbrPolyFacesPf.size())
50  {
52  << "Coupled patches " << fvp.name() << " and "
53  << nbrFvp.name() << " are not the same size"
54  << exit(FatalError);
55  }
56 
57  if (polyFacesPf.size())
58  {
59  for (label i = 1; i < polyFacesPf.size(); ++ i)
60  {
61  if
62  (
63  polyFacesPf[i - 1] > polyFacesPf[i]
64  ? true
65  : polyFacesPf[i - 1] == polyFacesPf[i]
66  ? nbrPolyFacesPf[i - 1] >= nbrPolyFacesPf[i]
67  : false
68  )
69  {
71  << "Coupled patches " << fvp.name()
72  << " and " << nbrFvp.name()
73  << " are not in order";
74 
75  if (procs[0] == procs[1])
76  {
78  << " on process #" << procs[0];
79  }
80 
82  << exit(FatalError);
83  }
84  }
85  }
86 
88  << "Coupled patches " << fvp.name()
89  << " and " << nbrFvp.name()
90  << " are in order" << endl;
91 }
92 
93 
95 (
96  const label& proci,
97  const fvPatch& fvp,
98  const labelList& polyFacesPf
99 )
100 {
101  if (fvp.size())
102  {
103  for (label i = 1; i < fvp.size(); ++ i)
104  {
105  if (polyFacesPf[i - 1] > polyFacesPf[i])
106  {
108  << "Error patch " << fvp.name()
109  << " is not in order";
110 
111  if (proci > 0)
112  {
114  << " on process #" << proci;
115  }
116 
118  << exit(FatalError);
119  }
120  }
121  }
122 
124  << "Error patch " << fvp.name()
125  << " is in order" << endl;
126 }
127 
128 
130 (
131  const fvMesh& completeMesh,
132  const multiDomainDecomposition& regionMeshes
133 )
134 {
135  forAll(completeMesh.boundary(), patchi)
136  {
137  const fvPatch& fvp = completeMesh.boundary()[patchi];
138 
139  // Cyclic patches
140  if
141  (
142  isA<nonConformalCyclicFvPatch>(fvp)
143  && refCast<const nonConformalCoupledFvPatch>(fvp).owner()
144  )
145  {
146  const label nccPatchi = patchi;
147  const label nbrNccPatchi =
148  refCast<const nonConformalCyclicFvPatch>(fvp)
149  .nbrPatchIndex();
150 
152  (
153  {-labelMax, labelMax},
154  completeMesh.boundary()[nccPatchi],
155  completeMesh.boundary()[nbrNccPatchi],
156  completeMesh.polyFacesBf()[nccPatchi],
157  completeMesh.polyFacesBf()[nbrNccPatchi]
158  );
159  }
160 
161  // Mapped patches
162  if
163  (
164  isA<nonConformalMappedWallFvPatch>(fvp)
165  && refCast<const nonConformalMappedWallFvPatch>(fvp).owner()
166  )
167  {
168  const nonConformalMappedWallFvPatch& ncmwFvp =
169  refCast<const nonConformalMappedWallFvPatch>(fvp);
170 
171  const domainDecomposition& nbrDecomposition =
172  regionMeshes[ncmwFvp.nbrRegionName()]();
173 
174  const fvMesh& nbrCompleteMesh = nbrDecomposition.completeMesh();
175 
176  const label ncmwPatchi = patchi;
177  const label nbrNcmwPatchi =
178  nbrCompleteMesh.boundary()[ncmwFvp.nbrPatchName()].index();
179 
181  (
182  {-labelMax, labelMax},
183  completeMesh.boundary()[ncmwPatchi],
184  nbrCompleteMesh.boundary()[nbrNcmwPatchi],
185  completeMesh.polyFacesBf()[ncmwPatchi],
186  nbrCompleteMesh.polyFacesBf()[ncmwPatchi]
187  );
188  }
189 
190  // Error patches
191  if (isA<nonConformalErrorFvPatch>(fvp))
192  {
193  const label ncePatchi = patchi;
194 
196  (
197  -labelMax,
198  completeMesh.boundary()[ncePatchi],
199  completeMesh.polyFacesBf()[ncePatchi]
200  );
201  }
202  }
203 }
204 
205 
207 (
208  const PtrList<fvMesh>& procMeshes,
209  const multiDomainDecomposition& regionMeshes
210 )
211 {
212  forAll(procMeshes, proci)
213  {
214  forAll(procMeshes[proci].boundary(), patchi)
215  {
216  const fvPatch& fvp =
217  procMeshes[proci].boundary()[patchi];
218 
219  // Cyclic patches
220  if
221  (
222  isA<nonConformalCoupledFvPatch>(fvp)
223  && refCast<const nonConformalCoupledFvPatch>(fvp).owner()
224  )
225  {
226  const label nccPatchi = patchi;
227 
228  label nbrProci = -1, nbrNccPatchi = -1;
229  if (isA<cyclicFvPatch>(fvp))
230  {
231  nbrProci = proci;
232  nbrNccPatchi =
233  refCast<const cyclicFvPatch>(fvp).nbrPatchIndex();
234  }
235  else if (isA<processorCyclicFvPatch>(fvp))
236  {
237  typedef processorCyclicFvPatch PcFvp;
238 
239  const PcFvp& pcFvp = refCast<const PcFvp>(fvp);
240 
241  nbrProci = pcFvp.neighbProcNo();
242 
243  const fvBoundaryMesh& nbrFvPatches =
244  procMeshes[nbrProci].boundary();
245 
246  forAll(nbrFvPatches, nbrNccPatchj)
247  {
248  const fvPatch& nbrFvp =
249  nbrFvPatches[nbrNccPatchj];
250 
251  if (isA<PcFvp>(nbrFvp))
252  {
253  const PcFvp& nbrPcFvp =
254  refCast<const PcFvp>(nbrFvp);
255 
256  if
257  (
258  nbrPcFvp.neighbProcNo()
259  == proci
260  && nbrPcFvp.referPatchIndex()
261  == pcFvp.referPatch().nbrPatchIndex()
262  )
263  {
264  nbrNccPatchi = nbrNccPatchj;
265  break;
266  }
267  }
268  }
269 
270  if (nbrNccPatchi == -1)
271  {
273  << "Opposite processor patch not found for "
274  << "patch " << fvp.name() << " on proc #"
275  << proci << exit(FatalError);
276  }
277  }
278  else
279  {
281  << "Non-conformal-coupled type not recognised "
282  << "for patch " << fvp.name() << " on proc #"
283  << proci << exit(FatalError);
284  }
285 
287  (
288  {proci, nbrProci},
289  procMeshes[proci].boundary()[nccPatchi],
290  procMeshes[nbrProci].boundary()[nbrNccPatchi],
291  procMeshes[proci].polyFacesBf()[nccPatchi],
292  procMeshes[nbrProci].polyFacesBf()[nbrNccPatchi]
293  );
294  }
295 
296  // Mapped patches
297  if
298  (
299  isA<nonConformalMappedWallFvPatch>(fvp)
300  && refCast<const nonConformalMappedWallFvPatch>(fvp).owner()
301  )
302  {
303  const label ncmwPatchi = patchi;
304 
305  const nonConformalMappedWallFvPatch& ncmwFvp =
306  refCast<const nonConformalMappedWallFvPatch>(fvp);
307 
308  typedef
310  NcmpfFvsplf;
311 
312  const NcmpfFvsplf& polyFacesPf =
313  refCast<const NcmpfFvsplf>
314  (procMeshes[proci].polyFacesBf()[ncmwPatchi]);
315 
316  forAll(procMeshes, nbrProci)
317  {
318  const fvMesh& nbrProcMesh =
319  regionMeshes[ncmwFvp.nbrRegionName()]()
320  .procMeshes()[nbrProci];
321 
322  if (nbrProcMesh.conformal()) continue;
323 
324  const label nbrNcmwPatchi =
325  nbrProcMesh
326  .boundary()[ncmwFvp.nbrPatchName()]
327  .index();
328 
329  const nonConformalMappedWallFvPatch& nbrNcmwFvp =
330  refCast<const nonConformalMappedWallFvPatch>
331  (nbrProcMesh.boundary()[nbrNcmwPatchi]);
332 
333  const NcmpfFvsplf& nbrPolyFacesPf =
334  refCast<const NcmpfFvsplf>
335  (nbrProcMesh.polyFacesBf()[nbrNcmwPatchi]);
336 
338  (
339  {proci, nbrProci},
340  ncmwFvp,
341  nbrNcmwFvp,
343  (
344  procMeshes[proci].polyFacesBf()[ncmwPatchi],
345  polyFacesPf.procSizes()[nbrProci],
346  polyFacesPf.procOffsets()[nbrProci]
347  ),
349  (
350  nbrProcMesh.polyFacesBf()[nbrNcmwPatchi],
351  nbrPolyFacesPf.procSizes()[proci],
352  nbrPolyFacesPf.procOffsets()[proci]
353  )
354  );
355  }
356  }
357 
358  // Error patches
359  if (isA<nonConformalErrorFvPatch>(fvp))
360  {
361  const label ncePatchi = patchi;
362 
364  (
365  proci,
366  procMeshes[proci].boundary()[ncePatchi],
367  procMeshes[proci].polyFacesBf()[ncePatchi]
368  );
369  }
370  }
371  }
372 }
373 
374 
375 }
376 
377 
378 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
379 
380 bool Foam::domainDecomposition::sortReconstructNonConformalCyclicAddressing_ =
382  (
383  "sortReconstructNonConformalCyclicAddressing",
384  0
385  );
386 
387 
388 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
389 
390 bool Foam::domainDecomposition::completeConformal() const
391 {
392  return completeMesh().conformal();
393 }
394 
395 
396 bool Foam::domainDecomposition::procsConformal() const
397 {
398  forAll(procMeshes_, proci)
399  {
400  if (!procMeshes_[proci].conformal())
401  {
402  return false;
403  }
404  }
405 
406  return true;
407 }
408 
409 
410 Foam::labelList Foam::domainDecomposition::completeFaceAddressing() const
411 {
412  labelList result (completeMesh().nFaces(), -labelMax);
413 
414  forAll(procMeshes_, proci)
415  {
416  forAll(procFaceAddressing()[proci], procFacei)
417  {
418  const label facei =
419  mag(procFaceAddressing()[proci][procFacei]) - 1;
420 
421  result[facei] = result[facei] == -labelMax ? procFacei : -1;
422  }
423  }
424 
425  return result;
426 }
427 
428 
430 Foam::domainDecomposition::nonConformalCyclicProcCyclics() const
431 {
432  List<labelPairLabelTable> result(completeMesh().boundary().size());
433 
434  forAll(procMeshes_, proci)
435  {
436  const fvMesh& procMesh = procMeshes_[proci];
437 
438  forAll(procMesh.boundary(), procPatchi)
439  {
440  const fvPatch& fvp = procMesh.boundary()[procPatchi];
441 
442  if (isA<nonConformalCyclicFvPatch>(fvp))
443  {
444  result[procPatchi].insert
445  (
446  labelPair(proci, proci),
447  procPatchi
448  );
449  }
450 
451  if (isA<nonConformalProcessorCyclicFvPatch>(fvp))
452  {
453  const processorCyclicFvPatch& pcFvp =
454  refCast<const processorCyclicFvPatch>(fvp);
455 
456  result[pcFvp.referPatchIndex()].insert
457  (
458  labelPair(proci, pcFvp.neighbProcNo()),
459  procPatchi
460  );
461  }
462  }
463  }
464 
465  return result;
466 }
467 
468 
470 Foam::domainDecomposition::nonConformalMappedWallProcOffsets
471 (
472  const bool appendSize
473 ) const
474 {
475  PtrList<labelListList> result(completeMesh().boundary().size());
476 
477  const surfaceLabelField::Boundary& polyFacesBf =
478  completeMesh().polyFacesBf();
479 
480  forAll(completeMesh().boundary(), ncmwPatchi)
481  {
482  const fvPatch& fvp = completeMesh().boundary()[ncmwPatchi];
483 
484  if (!isA<nonConformalMappedWallFvPatch>(fvp)) continue;
485 
486  const nonConformalMappedWallFvPatch& ncmwFvp =
487  refCast<const nonConformalMappedWallFvPatch>(fvp);
488 
489  const domainDecomposition& nbrDecomposition =
490  regionMeshes_[ncmwFvp.nbrRegionName()]();
491 
492  const fvMesh& nbrCompleteMesh = nbrDecomposition.completeMesh();
493 
494  const label nbrNcmwPatchi =
495  nbrCompleteMesh.boundary()[ncmwFvp.nbrPatchName()].index();
496 
497  const surfaceLabelField::Boundary& nbrPolyFacesBf =
498  nbrCompleteMesh.polyFacesBf();
499 
500  result.set
501  (
502  ncmwPatchi,
503  new labelListList(nProcs(), labelList(nProcs() + appendSize, 0))
504  );
505 
506  // Determine the number of faces in each processor block of the
507  // mapped patches
508  labelListList& procNbrProcCounts = result[ncmwPatchi];
509  forAll(polyFacesBf[ncmwPatchi], ncmwPatchFacei)
510  {
511  const label facei = polyFacesBf[ncmwPatchi][ncmwPatchFacei];
512  const label celli = completeMesh().faceOwner()[facei];
513  const label proci = cellProc_[celli];
514 
515  const label nbrFacei =
516  nbrPolyFacesBf[nbrNcmwPatchi][ncmwPatchFacei];
517  const label nbrCelli =
518  nbrCompleteMesh.faceOwner()[nbrFacei];
519  const label nbrProci = nbrDecomposition.cellProc()[nbrCelli];
520 
521  procNbrProcCounts[proci][nbrProci] ++;
522  }
523 
524  // Convert the counts to cumulative sums (i.e., offsets)
525  forAll(procNbrProcCounts, proci)
526  {
527  label count = 0;
528 
529  forAll(procNbrProcCounts[proci], nbrProci)
530  {
531  const label count0 = count;
532  count += procNbrProcCounts[proci][nbrProci];
533  procNbrProcCounts[proci][nbrProci] = count0;
534  }
535 
536  if (appendSize)
537  {
538  procNbrProcCounts[proci].last() = count;
539  }
540  }
541  }
542 
543  return result;
544 }
545 
546 
547 void Foam::domainDecomposition::decomposeNonConformalCyclicAddressing
548 (
549  const label nccPatchi,
550  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
551  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
552 ) const
553 {
554  const surfaceLabelField::Boundary& polyFacesBf =
555  completeMesh().polyFacesBf();
556 
557  const nonConformalCyclicFvPatch& nccFvp =
558  refCast<const nonConformalCyclicFvPatch>
559  (
560  completeMesh().boundary()[nccPatchi]
561  );
562 
563  if (!nccFvp.owner()) return;
564 
565  const label nbrNccPatchi = nccFvp.nbrPatchIndex();
566 
567  forAll(polyFacesBf[nccPatchi], nccPatchFacei)
568  {
569  const label facei = polyFacesBf[nccPatchi][nccPatchFacei];
570  const label celli = completeMesh().faceOwner()[facei];
571  const label proci = cellProc_[celli];
572 
573  const label nbrFacei = polyFacesBf[nbrNccPatchi][nccPatchFacei];
574  const label nbrCelli = completeMesh().faceOwner()[nbrFacei];
575  const label nbrProci = cellProc_[nbrCelli];
576 
577  const label procNccPatchi =
578  nonConformalCyclicProcCyclics
579  [nccPatchi][labelPair(proci, nbrProci)];
580  const label nbrProcNccPatchi =
581  nonConformalCyclicProcCyclics
582  [nbrNccPatchi][labelPair(nbrProci, proci)];
583 
584  nonConformalProcFaceAddressingBf[proci][procNccPatchi]
585  .append(nccPatchFacei + 1);
586  nonConformalProcFaceAddressingBf[nbrProci][nbrProcNccPatchi]
587  .append(nccPatchFacei + 1);
588  }
589 }
590 
591 
592 void Foam::domainDecomposition::decomposeNonConformalMappedWallAddressing
593 (
594  const label ncmwPatchi,
595  PtrList<labelListList>& nonConformalMappedWallProcOffsets,
596  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
597 ) const
598 {
599  const surfaceLabelField::Boundary& polyFacesBf =
600  completeMesh().polyFacesBf();
601 
602  const nonConformalMappedWallFvPatch& ncmwFvp =
603  refCast<const nonConformalMappedWallFvPatch>
604  (
605  completeMesh().boundary()[ncmwPatchi]
606  );
607 
608  const domainDecomposition& nbrDecomposition =
609  regionMeshes_[ncmwFvp.nbrRegionName()]();
610 
611  const fvMesh& nbrCompleteMesh = nbrDecomposition.completeMesh();
612 
613  const surfaceLabelField::Boundary& nbrPolyFacesBf =
614  nbrCompleteMesh.polyFacesBf();
615 
616  const label nbrNcmwPatchi =
617  nbrCompleteMesh.boundary()[ncmwFvp.nbrPatchName()].index();
618 
619  // Resize the addressing as necessary
620  labelListList& procOffsets =
621  nonConformalMappedWallProcOffsets[ncmwPatchi];
622  forAll(procOffsets, proci)
623  {
624  nonConformalProcFaceAddressingBf[proci][ncmwPatchi]
625  .resize(procOffsets[proci].last());
626  }
627 
628  // Insert the poly face addressing into the result. Use the
629  // procOffsets array as the index into each processor block.
630  forAll(polyFacesBf[ncmwPatchi], ncmwPatchFacei)
631  {
632  const label facei = polyFacesBf[ncmwPatchi][ncmwPatchFacei];
633  const label celli = completeMesh().faceOwner()[facei];
634  const label proci = cellProc_[celli];
635 
636  const label nbrFacei = nbrPolyFacesBf[nbrNcmwPatchi][ncmwPatchFacei];
637  const label nbrCelli = nbrCompleteMesh.faceOwner()[nbrFacei];
638  const label nbrProci = nbrDecomposition.cellProc()[nbrCelli];
639 
640  nonConformalProcFaceAddressingBf
641  [proci][ncmwPatchi][procOffsets[proci][nbrProci] ++] =
642  ncmwPatchFacei + 1;
643  }
644 }
645 
646 
647 void Foam::domainDecomposition::decomposeNonConformalErrorAddressing
648 (
649  const label ncePatchi,
650  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
651 ) const
652 {
653  const surfaceLabelField::Boundary& polyFacesBf =
654  completeMesh().polyFacesBf();
655 
656  forAll(polyFacesBf[ncePatchi], ncePatchFacei)
657  {
658  const label facei = polyFacesBf[ncePatchi][ncePatchFacei];
659  const label celli = completeMesh().faceOwner()[facei];
660  const label proci = cellProc_[celli];
661 
662  nonConformalProcFaceAddressingBf[proci][ncePatchi]
663  .append(ncePatchFacei + 1);
664  }
665 }
666 
667 
668 void Foam::domainDecomposition::reconstructNonConformalCyclicAddressing
669 (
670  const label nccPatchi,
671  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
672  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
673 ) const
674 {
675  const nonConformalCyclicFvPatch& nccFvp =
676  refCast<const nonConformalCyclicFvPatch>
677  (
678  completeMesh().boundary()[nccPatchi]
679  );
680 
681  if (!nccFvp.owner()) return;
682 
683  const label nbrNccPatchi = nccFvp.nbrPatchIndex();
684 
685  label nccPatchFacei = 0;
686  labelPairLabelTable procNccPatchFaceis;
688  (
689  labelPairLabelTable,
690  nonConformalCyclicProcCyclics[nccPatchi],
691  iter
692  )
693  {
694  procNccPatchFaceis.insert(iter.key(), 0);
695  }
696 
697  while (true)
698  {
699  labelPair procNbrProc(labelMax, labelMax);
700  labelPair faceNbrFace(labelMax, labelMax);
701 
702  forAllConstIter(labelPairLabelTable, procNccPatchFaceis, iter)
703  {
704  const label proci = iter.key().first();
705  const label nbrProci = iter.key().second();
706 
707  const labelPair procNbrProcStar(proci, nbrProci);
708  const labelPair nbrProcProcStar(nbrProci, proci);
709 
710  const label procNccPatchi =
711  nonConformalCyclicProcCyclics[nccPatchi][procNbrProcStar];
712  const label nbrProcNccPatchi =
713  nonConformalCyclicProcCyclics[nbrNccPatchi][nbrProcProcStar];
714 
715  const label procNccPatchFacei = iter();
716  const label size =
717  procMeshes_[proci].polyFacesBf()[procNccPatchi].size();
718 
719  if (procNccPatchFacei >= size) continue;
720 
721  const label procFacei =
722  procMeshes_[proci].polyFacesBf()
723  [procNccPatchi][procNccPatchFacei];
724  const label nbrProcFacei =
725  procMeshes_[nbrProci].polyFacesBf()
726  [nbrProcNccPatchi][procNccPatchFacei];
727 
728  const labelPair faceNbrFaceStar
729  (
730  procFaceAddressing_[proci][procFacei] - 1,
731  procFaceAddressing_[nbrProci][nbrProcFacei] - 1
732  );
733 
734  if (faceNbrFace > faceNbrFaceStar)
735  {
736  procNbrProc = procNbrProcStar;
737  faceNbrFace = faceNbrFaceStar;
738  }
739  }
740 
741  if (faceNbrFace == labelPair(labelMax, labelMax))
742  {
743  break;
744  }
745  else
746  {
747  const label proci = procNbrProc.first();
748  const label nbrProci = procNbrProc.second();
749 
750  const labelPair nbrProcProc(nbrProci, proci);
751 
752  const label procNccPatchi =
753  nonConformalCyclicProcCyclics[nccPatchi][procNbrProc];
754  const label nbrProcNccPatchi =
755  nonConformalCyclicProcCyclics[nbrNccPatchi][nbrProcProc];
756 
757  nonConformalProcFaceAddressingBf[proci][procNccPatchi]
758  .append(nccPatchFacei + 1);
759  nonConformalProcFaceAddressingBf[nbrProci][nbrProcNccPatchi]
760  .append(nccPatchFacei + 1);
761 
762  nccPatchFacei ++;
763  procNccPatchFaceis[procNbrProc] ++;
764  }
765  }
766 }
767 
768 
769 void Foam::domainDecomposition::sortReconstructNonConformalCyclicAddressing
770 (
771  const label nccPatchi,
772  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
773  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
774 ) const
775 {
776  const nonConformalCyclicFvPatch& nccFvp =
777  refCast<const nonConformalCyclicFvPatch>
778  (
779  completeMesh().boundary()[nccPatchi]
780  );
781 
782  if (!nccFvp.owner()) return;
783 
784  const label nbrNccPatchi = nccFvp.nbrPatchIndex();
785 
786  // Resize the relevant patches
788  (
789  labelPairLabelTable,
790  nonConformalCyclicProcCyclics[nccPatchi],
791  iter
792  )
793  {
794  const label proci = iter.key().first();
795  const label nbrProci = iter.key().second();
796 
797  const label procNccPatchi = iter();
798  const label nbrProcNccPatchi =
799  nonConformalCyclicProcCyclics
800  [nbrNccPatchi][labelPair(nbrProci, proci)];
801 
802  const labelUList& procPolyFacesPf =
803  procMeshes_[proci].polyFacesBf()[procNccPatchi];
804  const labelUList& nbrProcPolyFacesPf =
805  procMeshes_[nbrProci].polyFacesBf()[nbrProcNccPatchi];
806 
807  nonConformalProcFaceAddressingBf[proci][procNccPatchi]
808  .resize(procPolyFacesPf.size());
809  nonConformalProcFaceAddressingBf[nbrProci][nbrProcNccPatchi]
810  .resize(nbrProcPolyFacesPf.size());
811  }
812 
813  // Obtain references to the original patches
814  const fvPatch& origFvp = nccFvp.origPatch();
815  const fvPatch& nbrOrigFvp = nccFvp.nbrPatch().origPatch();
816 
817  // Create a list of "indices", containing every label relevant to each
818  // non-conformal face
819  DynamicList<FixedList<label, 7>> indices;
821  (
822  labelPairLabelTable,
823  nonConformalCyclicProcCyclics[nccPatchi],
824  iter
825  )
826  {
827  const label proci = iter.key().first();
828  const label nbrProci = iter.key().second();
829 
830  const label procNccPatchi = iter();
831  const label nbrProcNccPatchi =
832  nonConformalCyclicProcCyclics
833  [nbrNccPatchi][labelPair(nbrProci, proci)];
834 
835  const labelUList& procPolyFacesPf =
836  procMeshes_[proci].polyFacesBf()[procNccPatchi];
837  const labelUList& nbrProcPolyFacesPf =
838  procMeshes_[nbrProci].polyFacesBf()[nbrProcNccPatchi];
839 
840  forAll(procPolyFacesPf, procPatchFacei)
841  {
842  const label procPolyFacei =
843  procPolyFacesPf[procPatchFacei];
844  const label nbrProcPolyFacei =
845  nbrProcPolyFacesPf[procPatchFacei];
846 
847  const label completePolyFacei =
848  procFaceAddressing_[proci][procPolyFacei] - 1;
849  const label nbcCompletePolyFacei =
850  procFaceAddressing_[nbrProci][nbrProcPolyFacei] - 1;
851 
852  indices.append
853  ({
854  proci,
855  nbrProci,
856  procNccPatchi,
857  nbrProcNccPatchi,
858  completePolyFacei - origFvp.start(),
859  nbcCompletePolyFacei - nbrOrigFvp.start(),
860  procPatchFacei
861  });
862  }
863  }
864 
865  // Sort the indices by the owner poly face, then by the neighbour poly face
867  (
868  indices,
869  [](const FixedList<label, 7>& a, const FixedList<label, 7>& b)
870  {
871  return labelPair(a[4], a[5]) < labelPair(b[4], b[5]);
872  }
873  );
874 
875  // Unpack into the addressing
876  forAll(indices, i)
877  {
878  const label proci = indices[i][0];
879  const label nbrProci = indices[i][1];
880  const label procNccPatchi = indices[i][2];
881  const label nbrProcNccPatchi = indices[i][3];
882  //const label completePolyPatchFacei = indices[i][4];
883  //const label nbrCompletePolyPatchFacei = indices[i][5];
884  const label procPatchFacei = indices[i][6];
885 
886  nonConformalProcFaceAddressingBf
887  [proci][procNccPatchi][procPatchFacei] = i + 1;
888  nonConformalProcFaceAddressingBf
889  [nbrProci][nbrProcNccPatchi][procPatchFacei] = i + 1;
890  }
891 }
892 
893 
894 void Foam::domainDecomposition::reconstructNonConformalMappedWallAddressing
895 (
896  const label ncmwPatchi,
897  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
898 ) const
899 {
900  const nonConformalMappedWallFvPatch& ncmwFvp =
901  refCast<const nonConformalMappedWallFvPatch>
902  (
903  completeMesh().boundary()[ncmwPatchi]
904  );
905 
906  const bool owner = ncmwFvp.owner();
907 
908  const domainDecomposition& nbrDecomposition =
909  regionMeshes_[ncmwFvp.nbrRegionName()]();
910 
911  const label nbrNcmwPatchi =
912  nbrDecomposition.completeMesh()
913  .boundary()[ncmwFvp.nbrPatchName()]
914  .index();
915 
916  auto calcProcOffsets = []
917  (
918  const domainDecomposition& meshes,
919  const label ncmwPatchi
920  )
921  {
922  labelListList result(meshes.nProcs());
923 
924  forAll(meshes.procMeshes(), proci)
925  {
926  typedef
927  nonConformalMappedPolyFacesFvsPatchLabelField
928  NcmpfFvsplf;
929 
930  const NcmpfFvsplf& polyFacesPf =
931  refCast<const NcmpfFvsplf>
932  (meshes.procMeshes()[proci].polyFacesBf()[ncmwPatchi]);
933 
934  result[proci] = polyFacesPf.procOffsets();
935  result[proci].append(polyFacesPf.size());
936  }
937 
938  return result;
939  };
940  const labelListList procOffsets =
941  calcProcOffsets(*this, ncmwPatchi);
942  const labelListList nbrProcOffsets =
943  calcProcOffsets(nbrDecomposition, nbrNcmwPatchi);
944 
945  label ncmwPatchFacei = 0;
946  labelListList procNcmwPatchFaceis(nProcs(), labelList(nProcs(), 0));
947 
948  while (true)
949  {
950  labelPair procNbrProc(labelMax, labelMax);
951  labelPair faceNbrFace(labelMax, labelMax);
952 
953  forAll(procNcmwPatchFaceis, proci)
954  {
955  forAll(procNcmwPatchFaceis[proci], nbrProci)
956  {
957  const labelPair procNbrProcStar(proci, nbrProci);
958 
959  const label procNcmwPatchFacei =
960  procNcmwPatchFaceis[proci][nbrProci];
961  const label size =
962  procOffsets[proci][nbrProci + 1]
963  - procOffsets[proci][nbrProci];
964 
965  if (procNcmwPatchFacei >= size) continue;
966 
967  const label procFacei =
968  procMeshes_[proci].polyFacesBf()
969  [ncmwPatchi]
970  [procNcmwPatchFacei + procOffsets[proci][nbrProci]];
971  const label nbrProcFacei =
972  nbrDecomposition.procMeshes()[nbrProci].polyFacesBf()
973  [nbrNcmwPatchi]
974  [procNcmwPatchFacei + nbrProcOffsets[nbrProci][proci]];
975 
976  const labelPair faceNbrFaceStar
977  (
978  procFaceAddressing_[proci][procFacei] - 1,
979  nbrDecomposition
980  .procFaceAddressing_[nbrProci][nbrProcFacei] - 1
981  );
982 
983  if
984  (
985  owner
986  ? faceNbrFace > faceNbrFaceStar
987  : reverse(faceNbrFace) > reverse(faceNbrFaceStar)
988  )
989  {
990  procNbrProc = procNbrProcStar;
991  faceNbrFace = faceNbrFaceStar;
992  }
993  }
994  }
995 
996  if (faceNbrFace == labelPair(labelMax, labelMax))
997  {
998  break;
999  }
1000  else
1001  {
1002  const label proci = procNbrProc.first();
1003  const label nbrProci = procNbrProc.second();
1004 
1005  nonConformalProcFaceAddressingBf[proci][ncmwPatchi]
1006  .append(ncmwPatchFacei + 1);
1007 
1008  ncmwPatchFacei ++;
1009  procNcmwPatchFaceis[proci][nbrProci] ++;
1010  }
1011  }
1012 }
1013 
1014 
1015 void Foam::domainDecomposition::reconstructNonConformalErrorAddressing
1016 (
1017  const label ncePatchi,
1018  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
1019 ) const
1020 {
1021  label ncePatchFacei = 0;
1022  labelList procNcePatchFaceis(nProcs(), 0);
1023 
1024  while (true)
1025  {
1026  label facei = labelMax, proci = labelMax;
1027 
1028  forAll(procNcePatchFaceis, procStari)
1029  {
1030  const label size =
1031  procMeshes_[procStari].polyFacesBf()[ncePatchi].size();
1032 
1033  if (procNcePatchFaceis[procStari] >= size) continue;
1034 
1035  const label procFacei =
1036  procMeshes_[procStari].polyFacesBf()
1037  [ncePatchi][procNcePatchFaceis[procStari]];
1038 
1039  const label faceStari =
1040  procFaceAddressing_[procStari][procFacei] - 1;
1041 
1042  if (facei > faceStari)
1043  {
1044  facei = faceStari;
1045  proci = procStari;
1046  }
1047  }
1048 
1049  if (facei == labelMax)
1050  {
1051  break;
1052  }
1053  else
1054  {
1055  nonConformalProcFaceAddressingBf[proci][ncePatchi]
1056  .append(ncePatchFacei + 1);
1057 
1058  ncePatchFacei ++;
1059  procNcePatchFaceis[proci] ++;
1060  }
1061  }
1062 }
1063 
1064 
1066 Foam::domainDecomposition::nonConformalProcFaceAddressingBf() const
1067 {
1068  validateComplete();
1069  validateProcs();
1070 
1071  // Build non-conformal finite volume face addressing for each processor
1072  List<List<DynamicList<label>>> result(nProcs());
1073  forAll(result, proci)
1074  {
1075  result[proci].resize
1076  (
1077  procMeshes_[proci].boundary().size()
1078  );
1079  }
1080 
1081  if (completeConformal() && procsConformal())
1082  {
1083  // Nothing to do
1084  }
1085  else if (!completeConformal())
1086  {
1087  const List<labelPairLabelTable> nonConformalCyclicProcCyclics =
1088  this->nonConformalCyclicProcCyclics();
1089 
1090  PtrList<labelListList> nonConformalMappedWallProcOffsets =
1091  this->nonConformalMappedWallProcOffsets(true);
1092 
1093  // Decompose non-conformal addressing
1094  forAll(completeMesh().boundary(), ncPatchi)
1095  {
1096  const fvPatch& fvp = completeMesh().boundary()[ncPatchi];
1097 
1098  if (isA<nonConformalCyclicFvPatch>(fvp))
1099  {
1100  decomposeNonConformalCyclicAddressing
1101  (
1102  ncPatchi,
1103  nonConformalCyclicProcCyclics,
1104  result
1105  );
1106  }
1107 
1108  if (isA<nonConformalMappedWallFvPatch>(fvp))
1109  {
1110  decomposeNonConformalMappedWallAddressing
1111  (
1112  ncPatchi,
1113  nonConformalMappedWallProcOffsets,
1114  result
1115  );
1116  }
1117 
1118  if (isA<nonConformalErrorFvPatch>(fvp))
1119  {
1120  decomposeNonConformalErrorAddressing(ncPatchi, result);
1121  }
1122  }
1123  }
1124  else // if (!procsConformal())
1125  {
1126  const List<labelPairLabelTable> nonConformalCyclicProcCyclics =
1127  this->nonConformalCyclicProcCyclics();
1128 
1129  // Reconstruct non-conformal addressing
1130  forAll(completeMesh().boundary(), ncPatchi)
1131  {
1132  const fvPatch& fvp = completeMesh().boundary()[ncPatchi];
1133 
1134  if (isA<nonConformalCyclicFvPatch>(fvp))
1135  {
1136  if (sortReconstructNonConformalCyclicAddressing_)
1137  {
1138  sortReconstructNonConformalCyclicAddressing
1139  (
1140  ncPatchi,
1141  nonConformalCyclicProcCyclics,
1142  result
1143  );
1144  }
1145  else
1146  {
1147  reconstructNonConformalCyclicAddressing
1148  (
1149  ncPatchi,
1150  nonConformalCyclicProcCyclics,
1151  result
1152  );
1153  }
1154  }
1155 
1156  if (isA<nonConformalMappedWallFvPatch>(fvp))
1157  {
1158  reconstructNonConformalMappedWallAddressing(ncPatchi, result);
1159  }
1160 
1161  if (isA<nonConformalErrorFvPatch>(fvp))
1162  {
1163  reconstructNonConformalErrorAddressing(ncPatchi, result);
1164  }
1165  }
1166  }
1167 
1168  return result;
1169 }
1170 
1171 
1172 void Foam::domainDecomposition::unconformComplete()
1173 {
1174  surfaceLabelField::Boundary polyFacesBf
1175  (
1177  completeMesh().polyFacesBf()
1178  );
1179  surfaceVectorField Sf(completeMesh().Sf().cloneUnSliced());
1180  surfaceVectorField Cf(completeMesh().Cf().cloneUnSliced());
1181 
1182  forAll(procMeshes_, proci)
1183  {
1184  const fvMesh& procMesh = procMeshes_[proci];
1185 
1186  const surfaceLabelField::Boundary& faceAddressingBf =
1187  procFaceAddressingBf()[proci];
1188 
1189  forAll(procMesh.boundary(), procNcPatchi)
1190  {
1191  const fvPatch& fvp = procMesh.boundary()[procNcPatchi];
1192 
1193  if (!isA<nonConformalFvPatch>(fvp)) continue;
1194 
1195  const label completeNcPatchi =
1196  isA<processorCyclicFvPatch>(fvp)
1197  ? refCast<const processorCyclicFvPatch>(fvp)
1198  .referPatchIndex()
1199  : procNcPatchi;
1200 
1201  const label size =
1202  max
1203  (
1204  max(mag(faceAddressingBf[procNcPatchi])),
1205  polyFacesBf[completeNcPatchi].size()
1206  );
1207 
1208  polyFacesBf[completeNcPatchi].resize(size, -1);
1209  polyFacesBf[completeNcPatchi].labelField::rmap
1210  (
1211  mag
1212  (
1213  labelField
1214  (
1215  procFaceAddressing_[proci],
1216  procMesh.polyFacesBf()[procNcPatchi]
1217  )
1218  ) - 1,
1219  mag(faceAddressingBf[procNcPatchi]) - 1
1220  );
1221 
1222  // Set dummy data for the face geometry. This should not be
1223  // used during reconstruction.
1224  Sf.boundaryFieldRef()[completeNcPatchi].resize(size, Zero);
1225  Cf.boundaryFieldRef()[completeNcPatchi].resize(size, Zero);
1226  }
1227  }
1228 
1229  completeMesh_->unconform
1230  (
1231  polyFacesBf,
1232  Sf,
1233  Cf,
1234  NullObjectRef<surfaceScalarField>(),
1235  false
1236  );
1237 
1238  completeMesh_->setPolyFacesBfInstance(procMeshes_[0].polyFacesBfInstance());
1239 
1240  if (debug)
1241  {
1242  checkCompleteMeshOrdering(completeMesh(), regionMeshes_);
1243  }
1244 }
1245 
1246 
1247 void Foam::domainDecomposition::unconformProcs()
1248 {
1249  const labelList completeFaceAddressing =
1250  this->completeFaceAddressing();
1251 
1252  const PtrList<labelListList> nonConformalMappedWallProcOffsets =
1253  this->nonConformalMappedWallProcOffsets(false);
1254 
1255  forAll(procMeshes_, proci)
1256  {
1257  fvMesh& procMesh = procMeshes_[proci];
1258 
1259  const surfaceLabelField::Boundary& faceAddressingBf =
1260  procFaceAddressingBf()[proci];
1261 
1262  surfaceLabelField::Boundary polyFacesBf
1263  (
1265  procMesh.polyFacesBf()
1266  );
1267  surfaceVectorField Sf(procMesh.Sf().cloneUnSliced());
1268  surfaceVectorField Cf(procMesh.Cf().cloneUnSliced());
1269 
1270  forAll(procMesh.boundary(), procNcPatchi)
1271  {
1272  const fvPatch& fvp = procMesh.boundary()[procNcPatchi];
1273 
1274  if (!isA<nonConformalFvPatch>(fvp)) continue;
1275 
1276  const label completeNcPatchi =
1277  isA<processorCyclicFvPatch>(fvp)
1278  ? refCast<const processorCyclicFvPatch>(fvp)
1279  .referPatchIndex()
1280  : procNcPatchi;
1281 
1282  polyFacesBf[procNcPatchi] =
1283  labelField
1284  (
1285  completeFaceAddressing,
1286  labelField
1287  (
1288  completeMesh().polyFacesBf()[completeNcPatchi],
1289  mag(faceAddressingBf[procNcPatchi]) - 1
1290  )
1291  );
1292 
1293  if (isA<nonConformalMappedWallFvPatch>(fvp))
1294  {
1295  typedef
1296  nonConformalMappedPolyFacesFvsPatchLabelField
1297  NcmpfFvsplf;
1298 
1299  refCast<NcmpfFvsplf>(polyFacesBf[procNcPatchi]).procOffsets() =
1300  nonConformalMappedWallProcOffsets[procNcPatchi][proci];
1301  }
1302 
1303  const label size = polyFacesBf[procNcPatchi].size();
1304 
1305  // Set dummy data for the face geometry. This should not be
1306  // used during decomposition.
1307  Sf.boundaryFieldRef()[procNcPatchi].resize(size, Zero);
1308  Cf.boundaryFieldRef()[procNcPatchi].resize(size, Zero);
1309  }
1310 
1311  procMesh.unconform
1312  (
1313  polyFacesBf,
1314  Sf,
1315  Cf,
1316  NullObjectRef<surfaceScalarField>(),
1317  false
1318  );
1319 
1320  procMesh.setPolyFacesBfInstance(completeMesh().polyFacesBfInstance());
1321  }
1322 
1323  // Check ordering
1324  if (debug)
1325  {
1326  checkProcMeshesOrdering(procMeshes(), regionMeshes_);
1327  }
1328 }
1329 
1330 
1331 void Foam::domainDecomposition::unconform()
1332 {
1333  if (completeConformal() && procsConformal())
1334  {
1335  // Nothing to do
1336  }
1337  else if (!completeConformal() && procsConformal())
1338  {
1339  unconformProcs();
1340  }
1341  else if (completeConformal() && !procsConformal())
1342  {
1343  unconformComplete();
1344  }
1345  else // if (!completeConformal() && !procsConformal())
1346  {
1347  // Assume everything is consistent and do nothing
1348  }
1349 }
1350 
1351 
1352 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A List obtained as a section of another List.
Definition: SubList.H:56
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
Automatic domain decomposition class for finite-volume meshes.
const fvMesh & completeMesh() const
Access the global mesh.
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:880
bool conformal() const
Return whether the fvMesh is conformal with the polyMesh.
Definition: fvMesh.C:874
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:138
virtual const word & name() const
Return name.
Definition: fvPatch.H:126
const word & nbrPatchName() const
Name of the patch to map from.
const word & nbrRegionName() const
Name of the region to map from.
Wall fv patch which can do non-conformal mapping of values from another potentially non-globally conf...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
volScalarField & b
Definition: createFields.H:25
#define InfoInFunction
Report an information message using Foam::Info.
int optimisationSwitch(const char *name, const int defaultValue=0)
Lookup optimisation switch or add default value.
Definition: debug.C:255
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
void checkNonConformalErrorPatchOrdering(const label &proci, const fvPatch &fvp, const labelList &polyFacesPf)
void checkProcMeshesOrdering(const PtrList< fvMesh > &procMeshes, const multiDomainDecomposition &regionMeshes)
void reverse(UList< T > &, const label n)
Definition: UListI.H:334
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
error FatalError
void stableSort(UList< T > &)
Definition: UList.C:129
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
static const label labelMax
Definition: label.H:62
SurfaceField< vector > surfaceVectorField
UList< label > labelUList
Definition: UList.H:65
void checkCompleteMeshOrdering(const fvMesh &completeMesh, const multiDomainDecomposition &regionMeshes)
void checkNonConformalCoupledPatchOrdering(const labelPair &procs, const fvPatch &fvp, const fvPatch &nbrFvp, const labelUList &polyFacesPf, const labelUList &nbrPolyFacesPf)
faceListList boundary(nPatches)