xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision 4849c82aea8bbb0f1c72cc370030898e85d96b88)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    Defines a block Jacobi preconditioner.
34b9ad928SBarry Smith */
400e125f8SBarry Smith 
500e125f8SBarry Smith #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/
64b9ad928SBarry Smith 
76849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat);
86849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat);
95a7e1895SHong Zhang static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
104b9ad928SBarry Smith 
11d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi(PC pc)
12d71ae5a4SJacob Faibussowitsch {
134b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
144b9ad928SBarry Smith   Mat         mat = pc->mat, pmat = pc->pmat;
15976e8c5aSLisandro Dalcin   PetscBool   hasop;
1669a612a9SBarry Smith   PetscInt    N, M, start, i, sum, end;
1769a612a9SBarry Smith   PetscInt    bs, i_start = -1, i_end = -1;
1869a612a9SBarry Smith   PetscMPIInt rank, size;
194b9ad928SBarry Smith 
204b9ad928SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
249566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(pc->pmat, &bs));
254b9ad928SBarry Smith 
265a7e1895SHong Zhang   if (jac->n > 0 && jac->n < size) {
279566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiproc(pc));
283ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
295a7e1895SHong Zhang   }
305a7e1895SHong Zhang 
31f1580f4eSBarry Smith   /*    Determines the number of blocks assigned to each processor */
324b9ad928SBarry Smith   /*   local block count  given */
334b9ad928SBarry Smith   if (jac->n_local > 0 && jac->n < 0) {
341c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
354b9ad928SBarry Smith     if (jac->l_lens) { /* check that user set these correctly */
364b9ad928SBarry Smith       sum = 0;
374b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
3808401ef6SPierre Jolivet         PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
394b9ad928SBarry Smith         sum += jac->l_lens[i];
404b9ad928SBarry Smith       }
4108401ef6SPierre Jolivet       PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly");
424b9ad928SBarry Smith     } else {
439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
442fa5cd67SKarl Rupp       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
454b9ad928SBarry Smith     }
464b9ad928SBarry Smith   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
474b9ad928SBarry Smith     /* global blocks given: determine which ones are local */
484b9ad928SBarry Smith     if (jac->g_lens) {
494b9ad928SBarry Smith       /* check if the g_lens is has valid entries */
504b9ad928SBarry Smith       for (i = 0; i < jac->n; i++) {
517827d75bSBarry Smith         PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed");
5208401ef6SPierre Jolivet         PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
534b9ad928SBarry Smith       }
544b9ad928SBarry Smith       if (size == 1) {
554b9ad928SBarry Smith         jac->n_local = jac->n;
569566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
579566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local));
584b9ad928SBarry Smith         /* check that user set these correctly */
594b9ad928SBarry Smith         sum = 0;
604b9ad928SBarry Smith         for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i];
6108401ef6SPierre Jolivet         PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly");
624b9ad928SBarry Smith       } else {
639566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end));
64aaa8cc7dSPierre Jolivet         /* loop over blocks determining first one owned by me */
654b9ad928SBarry Smith         sum = 0;
664b9ad928SBarry Smith         for (i = 0; i < jac->n + 1; i++) {
679371c9d4SSatish Balay           if (sum == start) {
689371c9d4SSatish Balay             i_start = i;
699371c9d4SSatish Balay             goto start_1;
709371c9d4SSatish Balay           }
714b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
724b9ad928SBarry Smith         }
73e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
744b9ad928SBarry Smith       start_1:
754b9ad928SBarry Smith         for (i = i_start; i < jac->n + 1; i++) {
769371c9d4SSatish Balay           if (sum == end) {
779371c9d4SSatish Balay             i_end = i;
789371c9d4SSatish Balay             goto end_1;
799371c9d4SSatish Balay           }
804b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
814b9ad928SBarry Smith         }
82e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
834b9ad928SBarry Smith       end_1:
844b9ad928SBarry Smith         jac->n_local = i_end - i_start;
859566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
869566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local));
874b9ad928SBarry Smith       }
884b9ad928SBarry Smith     } else { /* no global blocks given, determine then using default layout */
894b9ad928SBarry Smith       jac->n_local = jac->n / size + ((jac->n % size) > rank);
909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
914b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
924b9ad928SBarry Smith         jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs;
937827d75bSBarry Smith         PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given");
944b9ad928SBarry Smith       }
954b9ad928SBarry Smith     }
964b9ad928SBarry Smith   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
974b9ad928SBarry Smith     jac->n       = size;
984b9ad928SBarry Smith     jac->n_local = 1;
999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->l_lens));
1004b9ad928SBarry Smith     jac->l_lens[0] = M;
1017271b979SBarry Smith   } else { /* jac->n > 0 && jac->n_local > 0 */
1027271b979SBarry Smith     if (!jac->l_lens) {
1039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
1042fa5cd67SKarl Rupp       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
1057271b979SBarry Smith     }
1064b9ad928SBarry Smith   }
10708401ef6SPierre Jolivet   PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors");
1084b9ad928SBarry Smith 
109f1580f4eSBarry Smith   /*    Determines mat and pmat */
1109566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
111976e8c5aSLisandro Dalcin   if (!hasop && size == 1) {
1124b9ad928SBarry Smith     mat  = pc->mat;
1134b9ad928SBarry Smith     pmat = pc->pmat;
1144b9ad928SBarry Smith   } else {
11549517cdeSBarry Smith     if (pc->useAmat) {
11649517cdeSBarry Smith       /* use block from Amat matrix, not Pmat for local MatMult() */
1179566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->mat, &mat));
1184b9ad928SBarry Smith     }
11949517cdeSBarry Smith     if (pc->pmat != pc->mat || !pc->useAmat) {
1209566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat));
1212fa5cd67SKarl Rupp     } else pmat = mat;
1224b9ad928SBarry Smith   }
1234b9ad928SBarry Smith 
124f1580f4eSBarry Smith   /*
1254b9ad928SBarry Smith      Setup code depends on the number of blocks
1264b9ad928SBarry Smith   */
127cc1d6079SHong Zhang   if (jac->n_local == 1) {
1289566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat));
1294b9ad928SBarry Smith   } else {
1309566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat));
1314b9ad928SBarry Smith   }
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith /* Default destroy, if it has never been setup */
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi(PC pc)
137d71ae5a4SJacob Faibussowitsch {
1384b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL));
1452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL));
1462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL));
1472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1504b9ad928SBarry Smith }
1514b9ad928SBarry Smith 
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject)
153d71ae5a4SJacob Faibussowitsch {
1544b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
155248bfaf0SJed Brown   PetscInt    blocks, i;
156ace3abfcSBarry Smith   PetscBool   flg;
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   PetscFunctionBegin;
159d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
1619566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
1629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
1639566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
164248bfaf0SJed Brown   if (jac->ksp) {
165248bfaf0SJed Brown     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
166248bfaf0SJed Brown      * unless we had already been called. */
16748a46eb9SPierre Jolivet     for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
168248bfaf0SJed Brown   }
169d0609cedSBarry Smith   PetscOptionsHeadEnd();
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1714b9ad928SBarry Smith }
1724b9ad928SBarry Smith 
1739804daf3SBarry Smith #include <petscdraw.h>
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer)
175d71ae5a4SJacob Faibussowitsch {
1764b9ad928SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
177cbe18068SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
17869a612a9SBarry Smith   PetscMPIInt           rank;
17969a612a9SBarry Smith   PetscInt              i;
180d9884438SBarry Smith   PetscBool             iascii, isstring, isdraw;
1814b9ad928SBarry Smith   PetscViewer           sviewer;
182020d6619SPierre Jolivet   PetscViewerFormat     format;
183020d6619SPierre Jolivet   const char           *prefix;
1844b9ad928SBarry Smith 
1854b9ad928SBarry Smith   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
18932077d6dSBarry Smith   if (iascii) {
19048a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
19163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of blocks = %" PetscInt_FMT "\n", jac->n));
1929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
1939566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
194020d6619SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
1969566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
1979566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
198933995ecSHong Zhang       if (jac->ksp && !jac->psubcomm) {
1999566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
200dd400576SPatrick Sanan         if (rank == 0) {
2019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2029566063dSJacob Faibussowitsch           PetscCall(KSPView(jac->ksp[0], sviewer));
2039566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2044b9ad928SBarry Smith         }
2059566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2069566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
207e4de9e1dSBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
209e4fa1014SBarry Smith       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
2109566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
211e4fa1014SBarry Smith         if (!mpjac->psubcomm->color) {
2129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2139566063dSJacob Faibussowitsch           PetscCall(KSPView(*(jac->ksp), sviewer));
2149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
215cbe18068SHong Zhang         }
2169566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2179566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2189530cbd7SBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2199566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2204b9ad928SBarry Smith       } else {
2219566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
222e4de9e1dSBarry Smith       }
223e4de9e1dSBarry Smith     } else {
2244814766eSBarry Smith       PetscInt n_global;
2251c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2269566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2279566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
2289371c9d4SSatish Balay       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2309566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
231dd2fa690SBarry Smith       for (i = 0; i < jac->n_local; i++) {
23263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2339566063dSJacob Faibussowitsch         PetscCall(KSPView(jac->ksp[i], sviewer));
2349566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2354b9ad928SBarry Smith       }
2369566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2404b9ad928SBarry Smith     }
2414b9ad928SBarry Smith   } else if (isstring) {
24263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2449566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2459566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
246d9884438SBarry Smith   } else if (isdraw) {
247d9884438SBarry Smith     PetscDraw draw;
248d9884438SBarry Smith     char      str[25];
249d9884438SBarry Smith     PetscReal x, y, bottom, h;
250d9884438SBarry Smith 
2519566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2529566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
25363a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2549566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
255d9884438SBarry Smith     bottom = y - h;
2569566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
257d9884438SBarry Smith     /* warning the communicator on viewer is different then on ksp in parallel */
2589566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2599566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2604b9ad928SBarry Smith   }
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2624b9ad928SBarry Smith }
2634b9ad928SBarry Smith 
264d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
265d71ae5a4SJacob Faibussowitsch {
266feb237baSPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2674b9ad928SBarry Smith 
2684b9ad928SBarry Smith   PetscFunctionBegin;
26928b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith   if (n_local) *n_local = jac->n_local;
2724b9ad928SBarry Smith   if (first_local) *first_local = jac->first_local;
273020d6619SPierre Jolivet   if (ksp) *ksp = jac->ksp;
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2754b9ad928SBarry Smith }
2764b9ad928SBarry Smith 
277f9663b93SPierre Jolivet static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt *lens)
278d71ae5a4SJacob Faibussowitsch {
2794b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith   PetscFunctionBegin;
2822472a847SBarry Smith   PetscCheck(pc->setupcalled <= 0 || jac->n == blocks, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
2834b9ad928SBarry Smith   jac->n = blocks;
2840a545947SLisandro Dalcin   if (!lens) jac->g_lens = NULL;
2852fa5cd67SKarl Rupp   else {
2869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2879566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2884b9ad928SBarry Smith   }
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2904b9ad928SBarry Smith }
2914b9ad928SBarry Smith 
292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
293d71ae5a4SJacob Faibussowitsch {
2944b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2954b9ad928SBarry Smith 
2964b9ad928SBarry Smith   PetscFunctionBegin;
2974b9ad928SBarry Smith   *blocks = jac->n;
2984b9ad928SBarry Smith   if (lens) *lens = jac->g_lens;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3004b9ad928SBarry Smith }
3014b9ad928SBarry Smith 
302d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[])
303d71ae5a4SJacob Faibussowitsch {
3044b9ad928SBarry Smith   PC_BJacobi *jac;
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith   PetscFunctionBegin;
3074b9ad928SBarry Smith   jac = (PC_BJacobi *)pc->data;
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith   jac->n_local = blocks;
3100a545947SLisandro Dalcin   if (!lens) jac->l_lens = NULL;
3112fa5cd67SKarl Rupp   else {
3129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3139566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3144b9ad928SBarry Smith   }
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3164b9ad928SBarry Smith }
3174b9ad928SBarry Smith 
318d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
319d71ae5a4SJacob Faibussowitsch {
3204b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith   PetscFunctionBegin;
3234b9ad928SBarry Smith   *blocks = jac->n_local;
3244b9ad928SBarry Smith   if (lens) *lens = jac->l_lens;
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3264b9ad928SBarry Smith }
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith /*@C
329f1580f4eSBarry Smith   PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3304b9ad928SBarry Smith   this processor.
3314b9ad928SBarry Smith 
3326da92b7fSHong Zhang   Not Collective
3334b9ad928SBarry Smith 
3344b9ad928SBarry Smith   Input Parameter:
3354b9ad928SBarry Smith . pc - the preconditioner context
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith   Output Parameters:
3380298fd71SBarry Smith + n_local     - the number of blocks on this processor, or NULL
3390298fd71SBarry Smith . first_local - the global number of the first block on this processor, or NULL
3404b9ad928SBarry Smith - ksp         - the array of KSP contexts
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith   Notes:
343f1580f4eSBarry Smith   After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith   Currently for some matrix implementations only 1 block per processor
3464b9ad928SBarry Smith   is supported.
3474b9ad928SBarry Smith 
348f1580f4eSBarry Smith   You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3494b9ad928SBarry Smith 
350feefa0e1SJacob Faibussowitsch   Fortran Notes:
351f1580f4eSBarry Smith   You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
352f1580f4eSBarry Smith 
353f1580f4eSBarry Smith   You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
354f1580f4eSBarry Smith   `KSP` array must be.
355196cc216SBarry Smith 
3564b9ad928SBarry Smith   Level: advanced
3574b9ad928SBarry Smith 
358562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3594b9ad928SBarry Smith @*/
360d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
361d71ae5a4SJacob Faibussowitsch {
3624b9ad928SBarry Smith   PetscFunctionBegin;
3630700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
364cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3664b9ad928SBarry Smith }
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith /*@
3694b9ad928SBarry Smith   PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3704b9ad928SBarry Smith   Jacobi preconditioner.
3714b9ad928SBarry Smith 
372c3339decSBarry Smith   Collective
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith   Input Parameters:
3754b9ad928SBarry Smith + pc     - the preconditioner context
3764b9ad928SBarry Smith . blocks - the number of blocks
3774b9ad928SBarry Smith - lens   - [optional] integer array containing the size of each block
3784b9ad928SBarry Smith 
3794b9ad928SBarry Smith   Options Database Key:
3804b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3814b9ad928SBarry Smith 
382f1580f4eSBarry Smith   Note:
3834b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
384f1580f4eSBarry Smith   All processors sharing the `PC` must call this routine with the same data.
3854b9ad928SBarry Smith 
3864b9ad928SBarry Smith   Level: intermediate
3874b9ad928SBarry Smith 
388562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3894b9ad928SBarry Smith @*/
390d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
391d71ae5a4SJacob Faibussowitsch {
3924b9ad928SBarry Smith   PetscFunctionBegin;
3930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39408401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
395cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3974b9ad928SBarry Smith }
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith /*@C
4004b9ad928SBarry Smith   PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
401f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4024b9ad928SBarry Smith 
403ad4df100SBarry Smith   Not Collective
4044b9ad928SBarry Smith 
4054b9ad928SBarry Smith   Input Parameter:
4064b9ad928SBarry Smith . pc - the preconditioner context
4074b9ad928SBarry Smith 
408feefa0e1SJacob Faibussowitsch   Output Parameters:
4094b9ad928SBarry Smith + blocks - the number of blocks
4104b9ad928SBarry Smith - lens   - integer array containing the size of each block
4114b9ad928SBarry Smith 
4124b9ad928SBarry Smith   Level: intermediate
4134b9ad928SBarry Smith 
414562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4154b9ad928SBarry Smith @*/
416d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
417d71ae5a4SJacob Faibussowitsch {
4184b9ad928SBarry Smith   PetscFunctionBegin;
4190700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4204f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
421cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4234b9ad928SBarry Smith }
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith /*@
4264b9ad928SBarry Smith   PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
427f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`,  preconditioner.
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith   Not Collective
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith   Input Parameters:
4324b9ad928SBarry Smith + pc     - the preconditioner context
4334b9ad928SBarry Smith . blocks - the number of blocks
4344b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4354b9ad928SBarry Smith 
436342c94f9SMatthew G. Knepley   Options Database Key:
437342c94f9SMatthew G. Knepley . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
438342c94f9SMatthew G. Knepley 
4394b9ad928SBarry Smith   Note:
4404b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4414b9ad928SBarry Smith 
4424b9ad928SBarry Smith   Level: intermediate
4434b9ad928SBarry Smith 
444562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4454b9ad928SBarry Smith @*/
446d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
447d71ae5a4SJacob Faibussowitsch {
4484b9ad928SBarry Smith   PetscFunctionBegin;
4490700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
45008401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
451cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4534b9ad928SBarry Smith }
4544b9ad928SBarry Smith 
4554b9ad928SBarry Smith /*@C
4564b9ad928SBarry Smith   PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
457f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4584b9ad928SBarry Smith 
4594b9ad928SBarry Smith   Not Collective
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith   Input Parameters:
4624b9ad928SBarry Smith + pc     - the preconditioner context
4634b9ad928SBarry Smith . blocks - the number of blocks
4644b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4654b9ad928SBarry Smith 
4664b9ad928SBarry Smith   Note:
4674b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith   Level: intermediate
4704b9ad928SBarry Smith 
471562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4724b9ad928SBarry Smith @*/
473d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
474d71ae5a4SJacob Faibussowitsch {
4754b9ad928SBarry Smith   PetscFunctionBegin;
4760700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4774f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
478cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4804b9ad928SBarry Smith }
4814b9ad928SBarry Smith 
4824b9ad928SBarry Smith /*MC
4834b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
484f1580f4eSBarry Smith            its own `KSP` object.
4854b9ad928SBarry Smith 
4864b9ad928SBarry Smith    Options Database Keys:
487011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
488011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4894b9ad928SBarry Smith 
49095452b02SPatrick Sanan    Notes:
491f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
492468ae2e8SBarry Smith 
49395452b02SPatrick Sanan     Each processor can have one or more blocks, or a single block can be shared by several processes. Defaults to one block per processor.
4944b9ad928SBarry Smith 
495f1580f4eSBarry Smith      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
496d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
4974b9ad928SBarry Smith 
498f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
499f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
5000462cc06SPierre Jolivet          `KSPGetPC()`)
5014b9ad928SBarry Smith 
502f1580f4eSBarry Smith      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
5032210c163SDominic Meiser          performance.  Different block partitioning may lead to additional data transfers
5042210c163SDominic Meiser          between host and GPU that lead to degraded performance.
5052210c163SDominic Meiser 
506011ea8aeSBarry Smith      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
507011ea8aeSBarry Smith 
5084b9ad928SBarry Smith    Level: beginner
5094b9ad928SBarry Smith 
510562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
511db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
512db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5134b9ad928SBarry Smith M*/
5144b9ad928SBarry Smith 
515d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
516d71ae5a4SJacob Faibussowitsch {
51769a612a9SBarry Smith   PetscMPIInt rank;
5184b9ad928SBarry Smith   PC_BJacobi *jac;
5194b9ad928SBarry Smith 
5204b9ad928SBarry Smith   PetscFunctionBegin;
5214dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
5229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5232fa5cd67SKarl Rupp 
5240a545947SLisandro Dalcin   pc->ops->apply           = NULL;
5257b6e2003SPierre Jolivet   pc->ops->matapply        = NULL;
5260a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
5274b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_BJacobi;
5284b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_BJacobi;
5294b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
5304b9ad928SBarry Smith   pc->ops->view            = PCView_BJacobi;
5310a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
5324b9ad928SBarry Smith 
5334b9ad928SBarry Smith   pc->data         = (void *)jac;
5344b9ad928SBarry Smith   jac->n           = -1;
5354b9ad928SBarry Smith   jac->n_local     = -1;
5364b9ad928SBarry Smith   jac->first_local = rank;
5370a545947SLisandro Dalcin   jac->ksp         = NULL;
5380a545947SLisandro Dalcin   jac->g_lens      = NULL;
5390a545947SLisandro Dalcin   jac->l_lens      = NULL;
5400a545947SLisandro Dalcin   jac->psubcomm    = NULL;
5414b9ad928SBarry Smith 
5429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5484b9ad928SBarry Smith }
5494b9ad928SBarry Smith 
5504b9ad928SBarry Smith /*
5514b9ad928SBarry Smith         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5524b9ad928SBarry Smith */
553d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
554d71ae5a4SJacob Faibussowitsch {
5554b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5564b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5574b9ad928SBarry Smith 
5584b9ad928SBarry Smith   PetscFunctionBegin;
5599566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563e91c6855SBarry Smith }
564e91c6855SBarry Smith 
565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
566d71ae5a4SJacob Faibussowitsch {
567e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
568e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
569e91c6855SBarry Smith 
570e91c6855SBarry Smith   PetscFunctionBegin;
5719566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5729566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5739566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5749566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5752e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5774b9ad928SBarry Smith }
5784b9ad928SBarry Smith 
579d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
580d71ae5a4SJacob Faibussowitsch {
5814b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5822295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
583539c167fSBarry Smith   KSPConvergedReason reason;
5844b9ad928SBarry Smith 
5854b9ad928SBarry Smith   PetscFunctionBegin;
5869566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5879566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
588ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5904b9ad928SBarry Smith }
5914b9ad928SBarry Smith 
592d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
593d71ae5a4SJacob Faibussowitsch {
5944b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5954b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5966e111a19SKarl Rupp 
5974b9ad928SBarry Smith   PetscFunctionBegin;
5989566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
5999566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
600bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
601910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
602910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6039566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6049566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
6059566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
6069566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
6079566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6094b9ad928SBarry Smith }
6104b9ad928SBarry Smith 
611d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
612d71ae5a4SJacob Faibussowitsch {
6137b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6147b6e2003SPierre Jolivet   Mat         sX, sY;
6157b6e2003SPierre Jolivet 
6167b6e2003SPierre Jolivet   PetscFunctionBegin;
6177b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6187b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6197b6e2003SPierre Jolivet      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6209566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6219566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(X, &sX));
6229566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
6239566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
6243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6257b6e2003SPierre Jolivet }
6267b6e2003SPierre Jolivet 
627d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
628d71ae5a4SJacob Faibussowitsch {
6294b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6304b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
631d9ca1df4SBarry Smith   PetscScalar            *y_array;
632d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6334b9ad928SBarry Smith   PC                      subpc;
6344b9ad928SBarry Smith 
6354b9ad928SBarry Smith   PetscFunctionBegin;
6364b9ad928SBarry Smith   /*
6374b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6384b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6394b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6404b9ad928SBarry Smith     for the sequential solve.
6414b9ad928SBarry Smith   */
6429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6449566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6459566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6464b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
6474b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6489566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6499566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6509566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6519566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6554b9ad928SBarry Smith }
6564b9ad928SBarry Smith 
657d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
658d71ae5a4SJacob Faibussowitsch {
6594b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6604b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
661d9ca1df4SBarry Smith   PetscScalar            *y_array;
662d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6634b9ad928SBarry Smith   PC                      subpc;
6644b9ad928SBarry Smith 
6654b9ad928SBarry Smith   PetscFunctionBegin;
6664b9ad928SBarry Smith   /*
6674b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6684b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6694b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6704b9ad928SBarry Smith     for the sequential solve.
6714b9ad928SBarry Smith   */
6729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6749566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6759566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6764b9ad928SBarry Smith 
6774b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
6784b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6794b9ad928SBarry Smith 
6809566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6819566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6824b9ad928SBarry Smith 
68311e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->x));
68411e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->y));
6859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6884b9ad928SBarry Smith }
6894b9ad928SBarry Smith 
690d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
691d71ae5a4SJacob Faibussowitsch {
6924b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6934b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
694d9ca1df4SBarry Smith   PetscScalar            *y_array;
695d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6964b9ad928SBarry Smith 
6974b9ad928SBarry Smith   PetscFunctionBegin;
6984b9ad928SBarry Smith   /*
6994b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
7004b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
7014b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
7024b9ad928SBarry Smith     for the sequential solve.
7034b9ad928SBarry Smith   */
7049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
7059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
7069566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
7079566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
7089566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7099566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
7109566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
7119566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
7129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
7139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
7143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7154b9ad928SBarry Smith }
7164b9ad928SBarry Smith 
717d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
718d71ae5a4SJacob Faibussowitsch {
7194b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
72069a612a9SBarry Smith   PetscInt                m;
7214b9ad928SBarry Smith   KSP                     ksp;
7224b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
723de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7243f6dc190SJunchao Zhang   VecType                 vectype;
725ea41da7aSPierre Jolivet   const char             *prefix;
7264b9ad928SBarry Smith 
7274b9ad928SBarry Smith   PetscFunctionBegin;
7284b9ad928SBarry Smith   if (!pc->setupcalled) {
729c2efdce3SBarry Smith     if (!jac->ksp) {
7303821be0aSBarry Smith       PetscInt nestlevel;
7313821be0aSBarry Smith 
732302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7332fa5cd67SKarl Rupp 
7349566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7353821be0aSBarry Smith       PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
7363821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
7379566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7389566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7399566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7409566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7419566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7429566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7434b9ad928SBarry Smith 
744e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7454b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7464b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7477b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7484b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7494b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7504b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7514b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7524b9ad928SBarry Smith 
7539566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7544b9ad928SBarry Smith       jac->ksp[0] = ksp;
755c2efdce3SBarry Smith 
7564dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
7574b9ad928SBarry Smith       jac->data = (void *)bjac;
7584b9ad928SBarry Smith     } else {
759c2efdce3SBarry Smith       ksp  = jac->ksp[0];
760c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
761c2efdce3SBarry Smith     }
762c2efdce3SBarry Smith 
763c2efdce3SBarry Smith     /*
764c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
765c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
766c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
767c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
768c2efdce3SBarry Smith       KSPSolve() on the block.
769c2efdce3SBarry Smith     */
7709566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7719566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7729566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7739566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7749566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7759566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
776c2efdce3SBarry Smith   } else {
7774b9ad928SBarry Smith     ksp  = jac->ksp[0];
7784b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7794b9ad928SBarry Smith   }
7809566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
78149517cdeSBarry Smith   if (pc->useAmat) {
7829566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7839566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7844b9ad928SBarry Smith   } else {
7859566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7864b9ad928SBarry Smith   }
7879566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
78890f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
789248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
7909566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
791302a9ddcSMatthew Knepley   }
7923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7934b9ad928SBarry Smith }
7944b9ad928SBarry Smith 
795d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
796d71ae5a4SJacob Faibussowitsch {
7974b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
7984b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
79969a612a9SBarry Smith   PetscInt               i;
8004b9ad928SBarry Smith 
8014b9ad928SBarry Smith   PetscFunctionBegin;
802e91c6855SBarry Smith   if (bjac && bjac->pmat) {
8039566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
80448a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
805e91c6855SBarry Smith   }
8064b9ad928SBarry Smith 
8074b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
8089566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
809e91c6855SBarry Smith     if (bjac && bjac->x) {
8109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
8119566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
8129566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
8134b9ad928SBarry Smith     }
814e91c6855SBarry Smith   }
8159566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
8169566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
8173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
818e91c6855SBarry Smith }
819e91c6855SBarry Smith 
820d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
821d71ae5a4SJacob Faibussowitsch {
822e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
823c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
824e91c6855SBarry Smith   PetscInt               i;
825e91c6855SBarry Smith 
826e91c6855SBarry Smith   PetscFunctionBegin;
8279566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
828c2efdce3SBarry Smith   if (bjac) {
8299566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8309566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8319566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
832c2efdce3SBarry Smith   }
8339566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
83448a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8359566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8362e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8384b9ad928SBarry Smith }
8394b9ad928SBarry Smith 
840d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
841d71ae5a4SJacob Faibussowitsch {
8424b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
84369a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
844539c167fSBarry Smith   KSPConvergedReason reason;
8454b9ad928SBarry Smith 
8464b9ad928SBarry Smith   PetscFunctionBegin;
8474b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8489566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8499566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
850ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8514b9ad928SBarry Smith   }
8523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8534b9ad928SBarry Smith }
8544b9ad928SBarry Smith 
855d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
856d71ae5a4SJacob Faibussowitsch {
8574b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
85869a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8594b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
860d9ca1df4SBarry Smith   PetscScalar           *yin;
861d9ca1df4SBarry Smith   const PetscScalar     *xin;
86258ebbce7SBarry Smith 
8634b9ad928SBarry Smith   PetscFunctionBegin;
8649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8664b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8674b9ad928SBarry Smith     /*
8684b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8694b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8704b9ad928SBarry Smith        the global vector.
8714b9ad928SBarry Smith     */
8729566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8739566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8744b9ad928SBarry Smith 
8759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8769566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8779566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
879d11f3a42SBarry Smith 
8809566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8819566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8824b9ad928SBarry Smith   }
8839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8864b9ad928SBarry Smith }
8874b9ad928SBarry Smith 
888f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
889f4d694ddSBarry Smith {
890f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
891f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
892f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
893f4d694ddSBarry Smith   PetscScalar           *yin;
894f4d694ddSBarry Smith   const PetscScalar     *xin;
895f4d694ddSBarry Smith   PC                     subpc;
896f4d694ddSBarry Smith 
897f4d694ddSBarry Smith   PetscFunctionBegin;
898f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
899f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
900f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
9014b9ad928SBarry Smith     /*
902f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
903f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
904f4d694ddSBarry Smith        the global vector.
9054b9ad928SBarry Smith     */
906f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
907f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
908f4d694ddSBarry Smith 
909f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
910f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
911f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
912f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
913f4d694ddSBarry Smith     PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
914f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
915f4d694ddSBarry Smith 
916f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
917f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
918f4d694ddSBarry Smith   }
919f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
920f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
922f4d694ddSBarry Smith }
923f4d694ddSBarry Smith 
924f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
925f4d694ddSBarry Smith {
926f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
927f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
928f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
929f4d694ddSBarry Smith   PetscScalar           *yin;
930f4d694ddSBarry Smith   const PetscScalar     *xin;
931f4d694ddSBarry Smith   PC                     subpc;
932f4d694ddSBarry Smith 
933f4d694ddSBarry Smith   PetscFunctionBegin;
934f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
935f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
936f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
937f4d694ddSBarry Smith     /*
938f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
939f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
940f4d694ddSBarry Smith        the global vector.
941f4d694ddSBarry Smith     */
942f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
943f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
944f4d694ddSBarry Smith 
945f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
946f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
947f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
948f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
949f4d694ddSBarry Smith     PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
950f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
951f4d694ddSBarry Smith 
952f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
953f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
954f4d694ddSBarry Smith   }
955f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
956f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
958f4d694ddSBarry Smith }
959f4d694ddSBarry Smith 
960d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
961d71ae5a4SJacob Faibussowitsch {
9624b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
96369a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
9644b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
965d9ca1df4SBarry Smith   PetscScalar           *yin;
966d9ca1df4SBarry Smith   const PetscScalar     *xin;
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith   PetscFunctionBegin;
9699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
9709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
9714b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
9724b9ad928SBarry Smith     /*
9734b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
9744b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
9754b9ad928SBarry Smith        the global vector.
9764b9ad928SBarry Smith     */
9779566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9789566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9794b9ad928SBarry Smith 
9809566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
9819566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
9829566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
9839566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
984a7ff49e8SJed Brown 
9859566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
9869566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
9874b9ad928SBarry Smith   }
9889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
9899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9914b9ad928SBarry Smith }
9924b9ad928SBarry Smith 
993d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
994d71ae5a4SJacob Faibussowitsch {
9954b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
99669a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
997ea41da7aSPierre Jolivet   const char            *prefix;
9984b9ad928SBarry Smith   KSP                    ksp;
9994b9ad928SBarry Smith   Vec                    x, y;
10004b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
10014b9ad928SBarry Smith   PC                     subpc;
10024b9ad928SBarry Smith   IS                     is;
1003434a97beSBrad Aagaard   MatReuse               scall;
10043f6dc190SJunchao Zhang   VecType                vectype;
1005*4849c82aSBarry Smith   MatNullSpace          *nullsp_mat = NULL, *nullsp_pmat = NULL;
10064b9ad928SBarry Smith 
10074b9ad928SBarry Smith   PetscFunctionBegin;
10089566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10094b9ad928SBarry Smith 
10104b9ad928SBarry Smith   n_local = jac->n_local;
10114b9ad928SBarry Smith 
101249517cdeSBarry Smith   if (pc->useAmat) {
1013ace3abfcSBarry Smith     PetscBool same;
10149566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
101528b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10164b9ad928SBarry Smith   }
10174b9ad928SBarry Smith 
10184b9ad928SBarry Smith   if (!pc->setupcalled) {
10193821be0aSBarry Smith     PetscInt nestlevel;
10203821be0aSBarry Smith 
10214b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
1022c2efdce3SBarry Smith 
1023c2efdce3SBarry Smith     if (!jac->ksp) {
1024e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Multiblock;
10254b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Multiblock;
10264b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Multiblock;
10277b6e2003SPierre Jolivet       pc->ops->matapply            = NULL;
1028f4d694ddSBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Multiblock;
1029f4d694ddSBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10304b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Multiblock;
10314b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Multiblock;
10324b9ad928SBarry Smith 
10334dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
10349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
10359566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10369566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
10374b9ad928SBarry Smith 
10384b9ad928SBarry Smith       jac->data = (void *)bjac;
10399566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
10404b9ad928SBarry Smith 
10414b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
10429566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10433821be0aSBarry Smith         PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
10443821be0aSBarry Smith         PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
10459566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10469566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10479566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
10489566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
10499566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
10509566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10519566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10522fa5cd67SKarl Rupp 
1053c2efdce3SBarry Smith         jac->ksp[i] = ksp;
1054c2efdce3SBarry Smith       }
1055c2efdce3SBarry Smith     } else {
1056c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
1057c2efdce3SBarry Smith     }
10584b9ad928SBarry Smith 
1059c2efdce3SBarry Smith     start = 0;
10609566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
1061c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
10624b9ad928SBarry Smith       m = jac->l_lens[i];
10634b9ad928SBarry Smith       /*
10644b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
10654b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
10664b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
10674b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
10684b9ad928SBarry Smith       KSPSolve() on the block.
10694b9ad928SBarry Smith 
10704b9ad928SBarry Smith       */
10719566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10729566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10739566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
10749566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
10752fa5cd67SKarl Rupp 
10764b9ad928SBarry Smith       bjac->x[i]      = x;
10774b9ad928SBarry Smith       bjac->y[i]      = y;
10784b9ad928SBarry Smith       bjac->starts[i] = start;
10794b9ad928SBarry Smith 
10809566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
10814b9ad928SBarry Smith       bjac->is[i] = is;
10824b9ad928SBarry Smith 
10834b9ad928SBarry Smith       start += m;
10844b9ad928SBarry Smith     }
10854b9ad928SBarry Smith   } else {
10864b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
10874b9ad928SBarry Smith     /*
10884b9ad928SBarry Smith        Destroy the blocks from the previous iteration
10894b9ad928SBarry Smith     */
10904b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1091*4849c82aSBarry Smith       PetscCall(MatGetNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
10929566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
1093*4849c82aSBarry Smith       if (pc->useAmat) {
1094*4849c82aSBarry Smith         PetscCall(MatGetNullSpaces(n_local, bjac->mat, &nullsp_mat));
1095*4849c82aSBarry Smith         PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
1096*4849c82aSBarry Smith       }
10974b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
10982fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
10994b9ad928SBarry Smith   }
11004b9ad928SBarry Smith 
11019566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
1102*4849c82aSBarry Smith   if (nullsp_pmat) PetscCall(MatRestoreNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
1103*4849c82aSBarry Smith   if (pc->useAmat) {
1104*4849c82aSBarry Smith     PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
1105*4849c82aSBarry Smith     if (nullsp_mat) PetscCall(MatRestoreNullSpaces(n_local, bjac->mat, &nullsp_mat));
1106*4849c82aSBarry Smith   }
11074b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
11084b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
11099566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
11104b9ad928SBarry Smith 
11114b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
11129566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
111349517cdeSBarry Smith     if (pc->useAmat) {
11149566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11159566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11164b9ad928SBarry Smith     } else {
11179566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11184b9ad928SBarry Smith     }
11199566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
112048a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
112190f1c854SHong Zhang   }
11223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11234b9ad928SBarry Smith }
11245a7e1895SHong Zhang 
11255a7e1895SHong Zhang /*
1126fd0b8932SBarry Smith       These are for a single block with multiple processes
11275a7e1895SHong Zhang */
1128d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1129d71ae5a4SJacob Faibussowitsch {
1130fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1131fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1132fd0b8932SBarry Smith   KSPConvergedReason reason;
1133fd0b8932SBarry Smith 
1134fd0b8932SBarry Smith   PetscFunctionBegin;
11359566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
11369566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1137ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1139fd0b8932SBarry Smith }
1140fd0b8932SBarry Smith 
1141d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1142d71ae5a4SJacob Faibussowitsch {
11435a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11445a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11455a7e1895SHong Zhang 
11465a7e1895SHong Zhang   PetscFunctionBegin;
11479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
11489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
11499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
11509566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1152e91c6855SBarry Smith }
1153e91c6855SBarry Smith 
1154d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1155d71ae5a4SJacob Faibussowitsch {
1156e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1157e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1158e91c6855SBarry Smith 
1159e91c6855SBarry Smith   PetscFunctionBegin;
11609566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
11619566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
11629566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
11639566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11645a7e1895SHong Zhang 
11659566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
11662e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
11673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11685a7e1895SHong Zhang }
11695a7e1895SHong Zhang 
1170d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1171d71ae5a4SJacob Faibussowitsch {
11725a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11735a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1174d9ca1df4SBarry Smith   PetscScalar          *yarray;
1175d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1176539c167fSBarry Smith   KSPConvergedReason    reason;
11775a7e1895SHong Zhang 
11785a7e1895SHong Zhang   PetscFunctionBegin;
11795a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
11809566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
11819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
11829566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
11839566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
11845a7e1895SHong Zhang 
11855a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
11869566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11879566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
11889566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
11899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11909566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1191ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1192250cf88bSHong Zhang 
11939566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
11949566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
11959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
11969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
11973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11985a7e1895SHong Zhang }
11995a7e1895SHong Zhang 
1200d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1201d71ae5a4SJacob Faibussowitsch {
12027b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
12037b6e2003SPierre Jolivet   KSPConvergedReason reason;
12047b6e2003SPierre Jolivet   Mat                sX, sY;
12057b6e2003SPierre Jolivet   const PetscScalar *x;
12067b6e2003SPierre Jolivet   PetscScalar       *y;
12077b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
12087b6e2003SPierre Jolivet 
12097b6e2003SPierre Jolivet   PetscFunctionBegin;
12107b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
12119566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
12129566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
12139566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
12149566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
12159566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
12169566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
12179566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12189566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12199566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
12209566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
12219566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12229566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12239566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
12269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
12279566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12289566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
12299566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1230ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12327b6e2003SPierre Jolivet }
12337b6e2003SPierre Jolivet 
1234d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1235d71ae5a4SJacob Faibussowitsch {
12365a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
12375a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12385a7e1895SHong Zhang   PetscInt              m, n;
1239ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
12405a7e1895SHong Zhang   const char           *prefix;
1241de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
12423f6dc190SJunchao Zhang   VecType               vectype;
12431f62890fSKarl Rupp 
12445a7e1895SHong Zhang   PetscFunctionBegin;
12459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
124608401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12475a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12485a7e1895SHong Zhang   if (!pc->setupcalled) {
12493821be0aSBarry Smith     PetscInt nestlevel;
12503821be0aSBarry Smith 
1251de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
12524dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mpjac));
12535a7e1895SHong Zhang     jac->data = (void *)mpjac;
12545a7e1895SHong Zhang 
12555a7e1895SHong Zhang     /* initialize datastructure mpjac */
12565a7e1895SHong Zhang     if (!jac->psubcomm) {
12575a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
12589566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12599566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12609566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12615a7e1895SHong Zhang     }
12625a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1263306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
12645a7e1895SHong Zhang 
12655a7e1895SHong Zhang     /* Get matrix blocks of pmat */
12669566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12675a7e1895SHong Zhang 
12685a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
12699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
12709566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12713821be0aSBarry Smith     PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
12723821be0aSBarry Smith     PetscCall(KSPSetNestLevel(jac->ksp[0], nestlevel + 1));
12739566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12749566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12759566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12769566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12775a7e1895SHong Zhang 
12789566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
12799566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
12809566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
12819566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
12829566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
12835a7e1895SHong Zhang 
12845a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
12859566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
12869566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
12879566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
12889566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
12899566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
12909566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
12915a7e1895SHong Zhang 
1292fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1293e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
12945a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
12955a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
12967b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1297fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1298306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
12999e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
13005a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
13019566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
13029566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1303fc08c53fSHong Zhang     } else {
13049566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
13055a7e1895SHong Zhang     }
13069566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
13075a7e1895SHong Zhang   }
13085a7e1895SHong Zhang 
130948a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
13103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13115a7e1895SHong Zhang }
1312