xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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) {
201b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPushTab(sviewer));
2029566063dSJacob Faibussowitsch           PetscCall(KSPView(jac->ksp[0], sviewer));
203b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPopTab(sviewer));
2044b9ad928SBarry Smith         }
2059566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
206e4de9e1dSBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2079566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
208e4fa1014SBarry Smith       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
2099566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
210e4fa1014SBarry Smith         if (!mpjac->psubcomm->color) {
211b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPushTab(sviewer));
212*f4f49eeaSPierre Jolivet           PetscCall(KSPView(*jac->ksp, sviewer));
213b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPopTab(sviewer));
214cbe18068SHong Zhang         }
2159566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2169530cbd7SBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
218e4de9e1dSBarry Smith       }
219e4de9e1dSBarry Smith     } else {
2204814766eSBarry Smith       PetscInt n_global;
2211c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2229566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2239566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
2249566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2259566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
226b4025f61SBarry Smith       PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local));
227dd2fa690SBarry Smith       for (i = 0; i < jac->n_local; i++) {
228b4025f61SBarry Smith         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2299566063dSJacob Faibussowitsch         PetscCall(KSPView(jac->ksp[i], sviewer));
230b4025f61SBarry Smith         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
2314b9ad928SBarry Smith       }
2329566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
2349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2354b9ad928SBarry Smith     }
2364b9ad928SBarry Smith   } else if (isstring) {
23763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2389566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2399566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2409566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
241d9884438SBarry Smith   } else if (isdraw) {
242d9884438SBarry Smith     PetscDraw draw;
243d9884438SBarry Smith     char      str[25];
244d9884438SBarry Smith     PetscReal x, y, bottom, h;
245d9884438SBarry Smith 
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2479566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
24863a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2499566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
250d9884438SBarry Smith     bottom = y - h;
2519566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
252d9884438SBarry Smith     /* warning the communicator on viewer is different then on ksp in parallel */
2539566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2549566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2554b9ad928SBarry Smith   }
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2574b9ad928SBarry Smith }
2584b9ad928SBarry Smith 
259d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
260d71ae5a4SJacob Faibussowitsch {
261feb237baSPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2624b9ad928SBarry Smith 
2634b9ad928SBarry Smith   PetscFunctionBegin;
26428b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2654b9ad928SBarry Smith 
2664b9ad928SBarry Smith   if (n_local) *n_local = jac->n_local;
2674b9ad928SBarry Smith   if (first_local) *first_local = jac->first_local;
268020d6619SPierre Jolivet   if (ksp) *ksp = jac->ksp;
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2704b9ad928SBarry Smith }
2714b9ad928SBarry Smith 
272f9663b93SPierre Jolivet static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt *lens)
273d71ae5a4SJacob Faibussowitsch {
2744b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2754b9ad928SBarry Smith 
2764b9ad928SBarry Smith   PetscFunctionBegin;
2772472a847SBarry 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");
2784b9ad928SBarry Smith   jac->n = blocks;
2790a545947SLisandro Dalcin   if (!lens) jac->g_lens = NULL;
2802fa5cd67SKarl Rupp   else {
2819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2829566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2834b9ad928SBarry Smith   }
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2854b9ad928SBarry Smith }
2864b9ad928SBarry Smith 
287d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
288d71ae5a4SJacob Faibussowitsch {
2894b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2904b9ad928SBarry Smith 
2914b9ad928SBarry Smith   PetscFunctionBegin;
2924b9ad928SBarry Smith   *blocks = jac->n;
2934b9ad928SBarry Smith   if (lens) *lens = jac->g_lens;
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2954b9ad928SBarry Smith }
2964b9ad928SBarry Smith 
297d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[])
298d71ae5a4SJacob Faibussowitsch {
2994b9ad928SBarry Smith   PC_BJacobi *jac;
3004b9ad928SBarry Smith 
3014b9ad928SBarry Smith   PetscFunctionBegin;
3024b9ad928SBarry Smith   jac = (PC_BJacobi *)pc->data;
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith   jac->n_local = blocks;
3050a545947SLisandro Dalcin   if (!lens) jac->l_lens = NULL;
3062fa5cd67SKarl Rupp   else {
3079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3089566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3094b9ad928SBarry Smith   }
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3114b9ad928SBarry Smith }
3124b9ad928SBarry Smith 
313d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
314d71ae5a4SJacob Faibussowitsch {
3154b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith   PetscFunctionBegin;
3184b9ad928SBarry Smith   *blocks = jac->n_local;
3194b9ad928SBarry Smith   if (lens) *lens = jac->l_lens;
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3214b9ad928SBarry Smith }
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith /*@C
324f1580f4eSBarry Smith   PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3254b9ad928SBarry Smith   this processor.
3264b9ad928SBarry Smith 
3276da92b7fSHong Zhang   Not Collective
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith   Input Parameter:
3304b9ad928SBarry Smith . pc - the preconditioner context
3314b9ad928SBarry Smith 
3324b9ad928SBarry Smith   Output Parameters:
3330298fd71SBarry Smith + n_local     - the number of blocks on this processor, or NULL
3340298fd71SBarry Smith . first_local - the global number of the first block on this processor, or NULL
3354b9ad928SBarry Smith - ksp         - the array of KSP contexts
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith   Notes:
338f1580f4eSBarry Smith   After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith   Currently for some matrix implementations only 1 block per processor
3414b9ad928SBarry Smith   is supported.
3424b9ad928SBarry Smith 
343f1580f4eSBarry Smith   You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3444b9ad928SBarry Smith 
345feefa0e1SJacob Faibussowitsch   Fortran Notes:
346f1580f4eSBarry Smith   You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
347f1580f4eSBarry Smith 
348f1580f4eSBarry Smith   You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
349f1580f4eSBarry Smith   `KSP` array must be.
350196cc216SBarry Smith 
3514b9ad928SBarry Smith   Level: advanced
3524b9ad928SBarry Smith 
353562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3544b9ad928SBarry Smith @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
356d71ae5a4SJacob Faibussowitsch {
3574b9ad928SBarry Smith   PetscFunctionBegin;
3580700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
359cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3614b9ad928SBarry Smith }
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith /*@
3644b9ad928SBarry Smith   PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3654b9ad928SBarry Smith   Jacobi preconditioner.
3664b9ad928SBarry Smith 
367c3339decSBarry Smith   Collective
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith   Input Parameters:
3704b9ad928SBarry Smith + pc     - the preconditioner context
3714b9ad928SBarry Smith . blocks - the number of blocks
3724b9ad928SBarry Smith - lens   - [optional] integer array containing the size of each block
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith   Options Database Key:
3754b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3764b9ad928SBarry Smith 
377f1580f4eSBarry Smith   Note:
3784b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
379f1580f4eSBarry Smith   All processors sharing the `PC` must call this routine with the same data.
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith   Level: intermediate
3824b9ad928SBarry Smith 
383562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3844b9ad928SBarry Smith @*/
385d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
386d71ae5a4SJacob Faibussowitsch {
3874b9ad928SBarry Smith   PetscFunctionBegin;
3880700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
38908401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
390cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3924b9ad928SBarry Smith }
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith /*@C
3954b9ad928SBarry Smith   PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
396f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
3974b9ad928SBarry Smith 
398ad4df100SBarry Smith   Not Collective
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith   Input Parameter:
4014b9ad928SBarry Smith . pc - the preconditioner context
4024b9ad928SBarry Smith 
403feefa0e1SJacob Faibussowitsch   Output Parameters:
4044b9ad928SBarry Smith + blocks - the number of blocks
4054b9ad928SBarry Smith - lens   - integer array containing the size of each block
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith   Level: intermediate
4084b9ad928SBarry Smith 
409562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4104b9ad928SBarry Smith @*/
411d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
412d71ae5a4SJacob Faibussowitsch {
4134b9ad928SBarry Smith   PetscFunctionBegin;
4140700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4154f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
416cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4184b9ad928SBarry Smith }
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith /*@
4214b9ad928SBarry Smith   PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
422f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`,  preconditioner.
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith   Not Collective
4254b9ad928SBarry Smith 
4264b9ad928SBarry Smith   Input Parameters:
4274b9ad928SBarry Smith + pc     - the preconditioner context
4284b9ad928SBarry Smith . blocks - the number of blocks
4294b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4304b9ad928SBarry Smith 
431342c94f9SMatthew G. Knepley   Options Database Key:
432342c94f9SMatthew G. Knepley . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
433342c94f9SMatthew G. Knepley 
4344b9ad928SBarry Smith   Note:
4354b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith   Level: intermediate
4384b9ad928SBarry Smith 
439562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4404b9ad928SBarry Smith @*/
441d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
442d71ae5a4SJacob Faibussowitsch {
4434b9ad928SBarry Smith   PetscFunctionBegin;
4440700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
44508401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
446cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4484b9ad928SBarry Smith }
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith /*@C
4514b9ad928SBarry Smith   PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
452f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4534b9ad928SBarry Smith 
4544b9ad928SBarry Smith   Not Collective
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith   Input Parameters:
4574b9ad928SBarry Smith + pc     - the preconditioner context
4584b9ad928SBarry Smith . blocks - the number of blocks
4594b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith   Note:
4624b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith   Level: intermediate
4654b9ad928SBarry Smith 
466562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4674b9ad928SBarry Smith @*/
468d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
469d71ae5a4SJacob Faibussowitsch {
4704b9ad928SBarry Smith   PetscFunctionBegin;
4710700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4724f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
473cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4754b9ad928SBarry Smith }
4764b9ad928SBarry Smith 
4774b9ad928SBarry Smith /*MC
4784b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
479f1580f4eSBarry Smith            its own `KSP` object.
4804b9ad928SBarry Smith 
4814b9ad928SBarry Smith    Options Database Keys:
482011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
483011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4844b9ad928SBarry Smith 
48595452b02SPatrick Sanan    Notes:
486f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
487468ae2e8SBarry Smith 
48895452b02SPatrick 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.
4894b9ad928SBarry Smith 
490f1580f4eSBarry Smith      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
491d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
4924b9ad928SBarry Smith 
493f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
494f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
4950462cc06SPierre Jolivet          `KSPGetPC()`)
4964b9ad928SBarry Smith 
497f1580f4eSBarry Smith      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
4982210c163SDominic Meiser          performance.  Different block partitioning may lead to additional data transfers
4992210c163SDominic Meiser          between host and GPU that lead to degraded performance.
5002210c163SDominic Meiser 
501011ea8aeSBarry Smith      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
502011ea8aeSBarry Smith 
5034b9ad928SBarry Smith    Level: beginner
5044b9ad928SBarry Smith 
505562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
506db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
507db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5084b9ad928SBarry Smith M*/
5094b9ad928SBarry Smith 
510d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
511d71ae5a4SJacob Faibussowitsch {
51269a612a9SBarry Smith   PetscMPIInt rank;
5134b9ad928SBarry Smith   PC_BJacobi *jac;
5144b9ad928SBarry Smith 
5154b9ad928SBarry Smith   PetscFunctionBegin;
5164dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
5179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5182fa5cd67SKarl Rupp 
5190a545947SLisandro Dalcin   pc->ops->apply           = NULL;
5207b6e2003SPierre Jolivet   pc->ops->matapply        = NULL;
5210a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
5224b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_BJacobi;
5234b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_BJacobi;
5244b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
5254b9ad928SBarry Smith   pc->ops->view            = PCView_BJacobi;
5260a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith   pc->data         = (void *)jac;
5294b9ad928SBarry Smith   jac->n           = -1;
5304b9ad928SBarry Smith   jac->n_local     = -1;
5314b9ad928SBarry Smith   jac->first_local = rank;
5320a545947SLisandro Dalcin   jac->ksp         = NULL;
5330a545947SLisandro Dalcin   jac->g_lens      = NULL;
5340a545947SLisandro Dalcin   jac->l_lens      = NULL;
5350a545947SLisandro Dalcin   jac->psubcomm    = NULL;
5364b9ad928SBarry Smith 
5379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5434b9ad928SBarry Smith }
5444b9ad928SBarry Smith 
5454b9ad928SBarry Smith /*
5464b9ad928SBarry Smith         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5474b9ad928SBarry Smith */
548d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
549d71ae5a4SJacob Faibussowitsch {
5504b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5514b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5524b9ad928SBarry Smith 
5534b9ad928SBarry Smith   PetscFunctionBegin;
5549566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
5573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
558e91c6855SBarry Smith }
559e91c6855SBarry Smith 
560d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
561d71ae5a4SJacob Faibussowitsch {
562e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
563e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
564e91c6855SBarry Smith 
565e91c6855SBarry Smith   PetscFunctionBegin;
5669566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5679566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5702e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5724b9ad928SBarry Smith }
5734b9ad928SBarry Smith 
574d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
575d71ae5a4SJacob Faibussowitsch {
5764b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5772295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
578539c167fSBarry Smith   KSPConvergedReason reason;
5794b9ad928SBarry Smith 
5804b9ad928SBarry Smith   PetscFunctionBegin;
5819566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5829566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
583ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5854b9ad928SBarry Smith }
5864b9ad928SBarry Smith 
587d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
588d71ae5a4SJacob Faibussowitsch {
5894b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5904b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5916e111a19SKarl Rupp 
5924b9ad928SBarry Smith   PetscFunctionBegin;
5939566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
5949566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
595bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
596910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
597910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
5989566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
5999566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
6009566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
6019566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
6029566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
6033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6044b9ad928SBarry Smith }
6054b9ad928SBarry Smith 
606d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
607d71ae5a4SJacob Faibussowitsch {
6087b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6097b6e2003SPierre Jolivet   Mat         sX, sY;
6107b6e2003SPierre Jolivet 
6117b6e2003SPierre Jolivet   PetscFunctionBegin;
6127b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6137b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6147b6e2003SPierre Jolivet      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6159566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6169566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(X, &sX));
6179566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
6189566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
6193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6207b6e2003SPierre Jolivet }
6217b6e2003SPierre Jolivet 
622d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
623d71ae5a4SJacob Faibussowitsch {
6244b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6254b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
626d9ca1df4SBarry Smith   PetscScalar            *y_array;
627d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6284b9ad928SBarry Smith   PC                      subpc;
6294b9ad928SBarry Smith 
6304b9ad928SBarry Smith   PetscFunctionBegin;
6314b9ad928SBarry Smith   /*
6324b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6334b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6344b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6354b9ad928SBarry Smith     for the sequential solve.
6364b9ad928SBarry Smith   */
6379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6399566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6409566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6414b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
6424b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6439566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6449566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6459566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6469566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6504b9ad928SBarry Smith }
6514b9ad928SBarry Smith 
652d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
653d71ae5a4SJacob Faibussowitsch {
6544b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6554b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
656d9ca1df4SBarry Smith   PetscScalar            *y_array;
657d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6584b9ad928SBarry Smith   PC                      subpc;
6594b9ad928SBarry Smith 
6604b9ad928SBarry Smith   PetscFunctionBegin;
6614b9ad928SBarry Smith   /*
6624b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6634b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6644b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6654b9ad928SBarry Smith     for the sequential solve.
6664b9ad928SBarry Smith   */
6679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6699566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6709566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6714b9ad928SBarry Smith 
6724b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
6734b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6744b9ad928SBarry Smith 
6759566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6769566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6774b9ad928SBarry Smith 
67811e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->x));
67911e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->y));
6809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6834b9ad928SBarry Smith }
6844b9ad928SBarry Smith 
685d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
686d71ae5a4SJacob Faibussowitsch {
6874b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6884b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
689d9ca1df4SBarry Smith   PetscScalar            *y_array;
690d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6914b9ad928SBarry Smith 
6924b9ad928SBarry Smith   PetscFunctionBegin;
6934b9ad928SBarry Smith   /*
6944b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6954b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6964b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6974b9ad928SBarry Smith     for the sequential solve.
6984b9ad928SBarry Smith   */
6999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
7009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
7019566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
7029566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
7039566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7049566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
7059566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
7069566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
7079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
7089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7104b9ad928SBarry Smith }
7114b9ad928SBarry Smith 
712d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
713d71ae5a4SJacob Faibussowitsch {
7144b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
71569a612a9SBarry Smith   PetscInt                m;
7164b9ad928SBarry Smith   KSP                     ksp;
7174b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
718de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7193f6dc190SJunchao Zhang   VecType                 vectype;
720ea41da7aSPierre Jolivet   const char             *prefix;
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith   PetscFunctionBegin;
7234b9ad928SBarry Smith   if (!pc->setupcalled) {
724c2efdce3SBarry Smith     if (!jac->ksp) {
7253821be0aSBarry Smith       PetscInt nestlevel;
7263821be0aSBarry Smith 
727302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7282fa5cd67SKarl Rupp 
7299566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7303821be0aSBarry Smith       PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
7313821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
7329566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7339566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7349566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7359566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7369566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7379566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7384b9ad928SBarry Smith 
739e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7404b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7414b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7427b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7434b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7444b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7454b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7464b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7474b9ad928SBarry Smith 
7489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7494b9ad928SBarry Smith       jac->ksp[0] = ksp;
750c2efdce3SBarry Smith 
7514dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
7524b9ad928SBarry Smith       jac->data = (void *)bjac;
7534b9ad928SBarry Smith     } else {
754c2efdce3SBarry Smith       ksp  = jac->ksp[0];
755c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
756c2efdce3SBarry Smith     }
757c2efdce3SBarry Smith 
758c2efdce3SBarry Smith     /*
759c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
760c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
761c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
762c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
763c2efdce3SBarry Smith       KSPSolve() on the block.
764c2efdce3SBarry Smith     */
7659566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7669566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7679566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7689566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7699566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7709566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
771c2efdce3SBarry Smith   } else {
7724b9ad928SBarry Smith     ksp  = jac->ksp[0];
7734b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7744b9ad928SBarry Smith   }
7759566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
77649517cdeSBarry Smith   if (pc->useAmat) {
7779566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7789566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7794b9ad928SBarry Smith   } else {
7809566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7814b9ad928SBarry Smith   }
7829566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
78390f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
784248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
7859566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
786302a9ddcSMatthew Knepley   }
7873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7884b9ad928SBarry Smith }
7894b9ad928SBarry Smith 
790d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
791d71ae5a4SJacob Faibussowitsch {
7924b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
7934b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
79469a612a9SBarry Smith   PetscInt               i;
7954b9ad928SBarry Smith 
7964b9ad928SBarry Smith   PetscFunctionBegin;
797e91c6855SBarry Smith   if (bjac && bjac->pmat) {
7989566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
79948a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
800e91c6855SBarry Smith   }
8014b9ad928SBarry Smith 
8024b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
8039566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
804e91c6855SBarry Smith     if (bjac && bjac->x) {
8059566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
8069566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
8079566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
8084b9ad928SBarry Smith     }
809e91c6855SBarry Smith   }
8109566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
8119566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
8123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
813e91c6855SBarry Smith }
814e91c6855SBarry Smith 
815d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
816d71ae5a4SJacob Faibussowitsch {
817e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
818c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
819e91c6855SBarry Smith   PetscInt               i;
820e91c6855SBarry Smith 
821e91c6855SBarry Smith   PetscFunctionBegin;
8229566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
823c2efdce3SBarry Smith   if (bjac) {
8249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8259566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8269566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
827c2efdce3SBarry Smith   }
8289566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
82948a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8309566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8312e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8334b9ad928SBarry Smith }
8344b9ad928SBarry Smith 
835d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
836d71ae5a4SJacob Faibussowitsch {
8374b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
83869a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
839539c167fSBarry Smith   KSPConvergedReason reason;
8404b9ad928SBarry Smith 
8414b9ad928SBarry Smith   PetscFunctionBegin;
8424b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8439566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8449566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
845ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8464b9ad928SBarry Smith   }
8473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8484b9ad928SBarry Smith }
8494b9ad928SBarry Smith 
850d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
851d71ae5a4SJacob Faibussowitsch {
8524b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
85369a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8544b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
855d9ca1df4SBarry Smith   PetscScalar           *yin;
856d9ca1df4SBarry Smith   const PetscScalar     *xin;
85758ebbce7SBarry Smith 
8584b9ad928SBarry Smith   PetscFunctionBegin;
8599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8614b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8624b9ad928SBarry Smith     /*
8634b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8644b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8654b9ad928SBarry Smith        the global vector.
8664b9ad928SBarry Smith     */
8679566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8689566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8694b9ad928SBarry Smith 
8709566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8719566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8729566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8739566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
874d11f3a42SBarry Smith 
8759566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8769566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8774b9ad928SBarry Smith   }
8789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8814b9ad928SBarry Smith }
8824b9ad928SBarry Smith 
883f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
884f4d694ddSBarry Smith {
885f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
886f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
887f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
888f4d694ddSBarry Smith   PetscScalar           *yin;
889f4d694ddSBarry Smith   const PetscScalar     *xin;
890f4d694ddSBarry Smith   PC                     subpc;
891f4d694ddSBarry Smith 
892f4d694ddSBarry Smith   PetscFunctionBegin;
893f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
894f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
895f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
8964b9ad928SBarry Smith     /*
897f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
898f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
899f4d694ddSBarry Smith        the global vector.
9004b9ad928SBarry Smith     */
901f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
902f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
903f4d694ddSBarry Smith 
904f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
905f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
906f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
907f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
908f4d694ddSBarry Smith     PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
909f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
910f4d694ddSBarry Smith 
911f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
912f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
913f4d694ddSBarry Smith   }
914f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
915f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
917f4d694ddSBarry Smith }
918f4d694ddSBarry Smith 
919f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
920f4d694ddSBarry Smith {
921f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
922f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
923f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
924f4d694ddSBarry Smith   PetscScalar           *yin;
925f4d694ddSBarry Smith   const PetscScalar     *xin;
926f4d694ddSBarry Smith   PC                     subpc;
927f4d694ddSBarry Smith 
928f4d694ddSBarry Smith   PetscFunctionBegin;
929f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
930f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
931f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
932f4d694ddSBarry Smith     /*
933f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
934f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
935f4d694ddSBarry Smith        the global vector.
936f4d694ddSBarry Smith     */
937f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
938f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
939f4d694ddSBarry Smith 
940f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
941f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
942f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
943f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
944f4d694ddSBarry Smith     PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
945f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
946f4d694ddSBarry Smith 
947f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
948f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
949f4d694ddSBarry Smith   }
950f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
951f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
953f4d694ddSBarry Smith }
954f4d694ddSBarry Smith 
955d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
956d71ae5a4SJacob Faibussowitsch {
9574b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
95869a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
9594b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
960d9ca1df4SBarry Smith   PetscScalar           *yin;
961d9ca1df4SBarry Smith   const PetscScalar     *xin;
9624b9ad928SBarry Smith 
9634b9ad928SBarry Smith   PetscFunctionBegin;
9649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
9659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
9664b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
9674b9ad928SBarry Smith     /*
9684b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
9694b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
9704b9ad928SBarry Smith        the global vector.
9714b9ad928SBarry Smith     */
9729566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9739566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9744b9ad928SBarry Smith 
9759566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
9769566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
9779566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
9789566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
979a7ff49e8SJed Brown 
9809566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
9819566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
9824b9ad928SBarry Smith   }
9839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
9849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
9853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9864b9ad928SBarry Smith }
9874b9ad928SBarry Smith 
988d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
989d71ae5a4SJacob Faibussowitsch {
9904b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
99169a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
992ea41da7aSPierre Jolivet   const char            *prefix;
9934b9ad928SBarry Smith   KSP                    ksp;
9944b9ad928SBarry Smith   Vec                    x, y;
9954b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
9964b9ad928SBarry Smith   PC                     subpc;
9974b9ad928SBarry Smith   IS                     is;
998434a97beSBrad Aagaard   MatReuse               scall;
9993f6dc190SJunchao Zhang   VecType                vectype;
10004849c82aSBarry Smith   MatNullSpace          *nullsp_mat = NULL, *nullsp_pmat = NULL;
10014b9ad928SBarry Smith 
10024b9ad928SBarry Smith   PetscFunctionBegin;
10039566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10044b9ad928SBarry Smith 
10054b9ad928SBarry Smith   n_local = jac->n_local;
10064b9ad928SBarry Smith 
100749517cdeSBarry Smith   if (pc->useAmat) {
1008ace3abfcSBarry Smith     PetscBool same;
10099566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
101028b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10114b9ad928SBarry Smith   }
10124b9ad928SBarry Smith 
10134b9ad928SBarry Smith   if (!pc->setupcalled) {
10143821be0aSBarry Smith     PetscInt nestlevel;
10153821be0aSBarry Smith 
10164b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
1017c2efdce3SBarry Smith 
1018c2efdce3SBarry Smith     if (!jac->ksp) {
1019e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Multiblock;
10204b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Multiblock;
10214b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Multiblock;
10227b6e2003SPierre Jolivet       pc->ops->matapply            = NULL;
1023f4d694ddSBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Multiblock;
1024f4d694ddSBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10254b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Multiblock;
10264b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Multiblock;
10274b9ad928SBarry Smith 
10284dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
10299566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
10309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10319566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith       jac->data = (void *)bjac;
10349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
10354b9ad928SBarry Smith 
10364b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
10379566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10383821be0aSBarry Smith         PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
10393821be0aSBarry Smith         PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
10409566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10419566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10429566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
10439566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
10449566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
10459566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10469566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10472fa5cd67SKarl Rupp 
1048c2efdce3SBarry Smith         jac->ksp[i] = ksp;
1049c2efdce3SBarry Smith       }
1050c2efdce3SBarry Smith     } else {
1051c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
1052c2efdce3SBarry Smith     }
10534b9ad928SBarry Smith 
1054c2efdce3SBarry Smith     start = 0;
10559566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
1056c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
10574b9ad928SBarry Smith       m = jac->l_lens[i];
10584b9ad928SBarry Smith       /*
10594b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
10604b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
10614b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
10624b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
10634b9ad928SBarry Smith       KSPSolve() on the block.
10644b9ad928SBarry Smith 
10654b9ad928SBarry Smith       */
10669566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10679566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10689566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
10699566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
10702fa5cd67SKarl Rupp 
10714b9ad928SBarry Smith       bjac->x[i]      = x;
10724b9ad928SBarry Smith       bjac->y[i]      = y;
10734b9ad928SBarry Smith       bjac->starts[i] = start;
10744b9ad928SBarry Smith 
10759566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
10764b9ad928SBarry Smith       bjac->is[i] = is;
10774b9ad928SBarry Smith 
10784b9ad928SBarry Smith       start += m;
10794b9ad928SBarry Smith     }
10804b9ad928SBarry Smith   } else {
10814b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
10824b9ad928SBarry Smith     /*
10834b9ad928SBarry Smith        Destroy the blocks from the previous iteration
10844b9ad928SBarry Smith     */
10854b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
10864849c82aSBarry Smith       PetscCall(MatGetNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
10879566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
10884849c82aSBarry Smith       if (pc->useAmat) {
10894849c82aSBarry Smith         PetscCall(MatGetNullSpaces(n_local, bjac->mat, &nullsp_mat));
10904849c82aSBarry Smith         PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
10914849c82aSBarry Smith       }
10924b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
10932fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
10944b9ad928SBarry Smith   }
10954b9ad928SBarry Smith 
10969566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
10974849c82aSBarry Smith   if (nullsp_pmat) PetscCall(MatRestoreNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
10984849c82aSBarry Smith   if (pc->useAmat) {
10994849c82aSBarry Smith     PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
11004849c82aSBarry Smith     if (nullsp_mat) PetscCall(MatRestoreNullSpaces(n_local, bjac->mat, &nullsp_mat));
11014849c82aSBarry Smith   }
11024b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
11034b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
11049566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
11054b9ad928SBarry Smith 
11064b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
11079566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
110849517cdeSBarry Smith     if (pc->useAmat) {
11099566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11109566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11114b9ad928SBarry Smith     } else {
11129566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11134b9ad928SBarry Smith     }
11149566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
111548a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
111690f1c854SHong Zhang   }
11173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11184b9ad928SBarry Smith }
11195a7e1895SHong Zhang 
11205a7e1895SHong Zhang /*
1121fd0b8932SBarry Smith       These are for a single block with multiple processes
11225a7e1895SHong Zhang */
1123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1124d71ae5a4SJacob Faibussowitsch {
1125fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1126fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1127fd0b8932SBarry Smith   KSPConvergedReason reason;
1128fd0b8932SBarry Smith 
1129fd0b8932SBarry Smith   PetscFunctionBegin;
11309566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
11319566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1132ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1134fd0b8932SBarry Smith }
1135fd0b8932SBarry Smith 
1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1137d71ae5a4SJacob Faibussowitsch {
11385a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11395a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11405a7e1895SHong Zhang 
11415a7e1895SHong Zhang   PetscFunctionBegin;
11429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
11439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
11449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
11459566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1147e91c6855SBarry Smith }
1148e91c6855SBarry Smith 
1149d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1150d71ae5a4SJacob Faibussowitsch {
1151e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1152e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1153e91c6855SBarry Smith 
1154e91c6855SBarry Smith   PetscFunctionBegin;
11559566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
11569566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
11579566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
11589566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11595a7e1895SHong Zhang 
11609566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
11612e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
11623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11635a7e1895SHong Zhang }
11645a7e1895SHong Zhang 
1165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1166d71ae5a4SJacob Faibussowitsch {
11675a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11685a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1169d9ca1df4SBarry Smith   PetscScalar          *yarray;
1170d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1171539c167fSBarry Smith   KSPConvergedReason    reason;
11725a7e1895SHong Zhang 
11735a7e1895SHong Zhang   PetscFunctionBegin;
11745a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
11759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
11769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
11779566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
11789566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
11795a7e1895SHong Zhang 
11805a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
11819566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11829566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
11839566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
11849566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11859566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1186ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1187250cf88bSHong Zhang 
11889566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
11899566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
11909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
11919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
11923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11935a7e1895SHong Zhang }
11945a7e1895SHong Zhang 
1195d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1196d71ae5a4SJacob Faibussowitsch {
11977b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
11987b6e2003SPierre Jolivet   KSPConvergedReason reason;
11997b6e2003SPierre Jolivet   Mat                sX, sY;
12007b6e2003SPierre Jolivet   const PetscScalar *x;
12017b6e2003SPierre Jolivet   PetscScalar       *y;
12027b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
12037b6e2003SPierre Jolivet 
12047b6e2003SPierre Jolivet   PetscFunctionBegin;
12057b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
12069566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
12079566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
12089566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
12099566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
12109566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
12119566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
12129566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12139566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12149566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
12159566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
12169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12179566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12189566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
12219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
12229566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12239566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
12249566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1225ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12277b6e2003SPierre Jolivet }
12287b6e2003SPierre Jolivet 
1229d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1230d71ae5a4SJacob Faibussowitsch {
12315a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
12325a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12335a7e1895SHong Zhang   PetscInt              m, n;
1234ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
12355a7e1895SHong Zhang   const char           *prefix;
1236de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
12373f6dc190SJunchao Zhang   VecType               vectype;
12381f62890fSKarl Rupp 
12395a7e1895SHong Zhang   PetscFunctionBegin;
12409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
124108401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12425a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12435a7e1895SHong Zhang   if (!pc->setupcalled) {
12443821be0aSBarry Smith     PetscInt nestlevel;
12453821be0aSBarry Smith 
1246de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
12474dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mpjac));
12485a7e1895SHong Zhang     jac->data = (void *)mpjac;
12495a7e1895SHong Zhang 
12505a7e1895SHong Zhang     /* initialize datastructure mpjac */
12515a7e1895SHong Zhang     if (!jac->psubcomm) {
12525a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
12539566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12549566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12559566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12565a7e1895SHong Zhang     }
12575a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1258306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
12595a7e1895SHong Zhang 
12605a7e1895SHong Zhang     /* Get matrix blocks of pmat */
12619566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12625a7e1895SHong Zhang 
12635a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
12649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
12659566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12663821be0aSBarry Smith     PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
12673821be0aSBarry Smith     PetscCall(KSPSetNestLevel(jac->ksp[0], nestlevel + 1));
12689566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12699566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12709566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12719566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12725a7e1895SHong Zhang 
12739566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
12749566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
12759566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
12769566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
12779566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
12785a7e1895SHong Zhang 
12795a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
12809566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
12819566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
12829566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
12839566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
12849566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
12859566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
12865a7e1895SHong Zhang 
1287fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1288e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
12895a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
12905a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
12917b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1292fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1293306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
12949e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
12955a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
12969566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
12979566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1298fc08c53fSHong Zhang     } else {
12999566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
13005a7e1895SHong Zhang     }
13019566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
13025a7e1895SHong Zhang   }
13035a7e1895SHong Zhang 
130448a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
13053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13065a7e1895SHong Zhang }
1307