xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision 78d05565dabe5a0ad9c88b34bc8daea9a6d9da20)
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) {
34462c564dSBarry Smith     PetscCallMPI(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 
152ce78bad3SBarry Smith 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));
212f4f49eeaSPierre 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;
221462c564dSBarry Smith       PetscCallMPI(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 
337ce78bad3SBarry Smith   Level: advanced
338ce78bad3SBarry Smith 
3394b9ad928SBarry Smith   Notes:
340f1580f4eSBarry Smith   After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith   Currently for some matrix implementations only 1 block per processor
3434b9ad928SBarry Smith   is supported.
3444b9ad928SBarry Smith 
345f1580f4eSBarry Smith   You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3464b9ad928SBarry Smith 
347562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3484b9ad928SBarry Smith @*/
349d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
350d71ae5a4SJacob Faibussowitsch {
3514b9ad928SBarry Smith   PetscFunctionBegin;
3520700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
353cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith /*@
3584b9ad928SBarry Smith   PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3594b9ad928SBarry Smith   Jacobi preconditioner.
3604b9ad928SBarry Smith 
361c3339decSBarry Smith   Collective
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith   Input Parameters:
3644b9ad928SBarry Smith + pc     - the preconditioner context
3654b9ad928SBarry Smith . blocks - the number of blocks
3664b9ad928SBarry Smith - lens   - [optional] integer array containing the size of each block
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith   Options Database Key:
3694b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3704b9ad928SBarry Smith 
371ce78bad3SBarry Smith   Level: intermediate
372ce78bad3SBarry Smith 
373f1580f4eSBarry Smith   Note:
3744b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
375f1580f4eSBarry Smith   All processors sharing the `PC` must call this routine with the same data.
3764b9ad928SBarry Smith 
377562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3784b9ad928SBarry Smith @*/
379d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
380d71ae5a4SJacob Faibussowitsch {
3814b9ad928SBarry Smith   PetscFunctionBegin;
3820700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
38308401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
384cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3864b9ad928SBarry Smith }
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith /*@C
3894b9ad928SBarry Smith   PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
390f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
3914b9ad928SBarry Smith 
392ad4df100SBarry Smith   Not Collective
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith   Input Parameter:
3954b9ad928SBarry Smith . pc - the preconditioner context
3964b9ad928SBarry Smith 
397feefa0e1SJacob Faibussowitsch   Output Parameters:
3984b9ad928SBarry Smith + blocks - the number of blocks
3994b9ad928SBarry Smith - lens   - integer array containing the size of each block
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith   Level: intermediate
4024b9ad928SBarry Smith 
403562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4044b9ad928SBarry Smith @*/
405d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
406d71ae5a4SJacob Faibussowitsch {
4074b9ad928SBarry Smith   PetscFunctionBegin;
4080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4094f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
410cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4124b9ad928SBarry Smith }
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith /*@
4154b9ad928SBarry Smith   PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
416f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`,  preconditioner.
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith   Not Collective
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith   Input Parameters:
4214b9ad928SBarry Smith + pc     - the preconditioner context
4224b9ad928SBarry Smith . blocks - the number of blocks
4234b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4244b9ad928SBarry Smith 
425342c94f9SMatthew G. Knepley   Options Database Key:
426342c94f9SMatthew G. Knepley . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
427342c94f9SMatthew G. Knepley 
428ce78bad3SBarry Smith   Level: intermediate
429ce78bad3SBarry Smith 
4304b9ad928SBarry Smith   Note:
4314b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4324b9ad928SBarry Smith 
433562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4344b9ad928SBarry Smith @*/
435d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
436d71ae5a4SJacob Faibussowitsch {
4374b9ad928SBarry Smith   PetscFunctionBegin;
4380700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
43908401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
440cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4424b9ad928SBarry Smith }
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith /*@C
4454b9ad928SBarry Smith   PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
446f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4474b9ad928SBarry Smith 
4484b9ad928SBarry Smith   Not Collective
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith   Input Parameters:
4514b9ad928SBarry Smith + pc     - the preconditioner context
4524b9ad928SBarry Smith . blocks - the number of blocks
4534b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4544b9ad928SBarry Smith 
455ce78bad3SBarry Smith   Level: intermediate
456ce78bad3SBarry Smith 
4574b9ad928SBarry Smith   Note:
4584b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4594b9ad928SBarry Smith 
460562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4614b9ad928SBarry Smith @*/
462d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
463d71ae5a4SJacob Faibussowitsch {
4644b9ad928SBarry Smith   PetscFunctionBegin;
4650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4664f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
467cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4694b9ad928SBarry Smith }
4704b9ad928SBarry Smith 
4714b9ad928SBarry Smith /*MC
4724b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
473f1580f4eSBarry Smith            its own `KSP` object.
4744b9ad928SBarry Smith 
4754b9ad928SBarry Smith    Options Database Keys:
476011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
477011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4784b9ad928SBarry Smith 
479ce78bad3SBarry Smith    Level: beginner
480ce78bad3SBarry Smith 
48195452b02SPatrick Sanan    Notes:
482f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
483468ae2e8SBarry Smith 
48495452b02SPatrick 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.
4854b9ad928SBarry Smith 
486f1580f4eSBarry Smith      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
487d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
4884b9ad928SBarry Smith 
489f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
490f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
4910462cc06SPierre Jolivet          `KSPGetPC()`)
4924b9ad928SBarry Smith 
493f1580f4eSBarry Smith      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
4942210c163SDominic Meiser          performance.  Different block partitioning may lead to additional data transfers
4952210c163SDominic Meiser          between host and GPU that lead to degraded performance.
4962210c163SDominic Meiser 
497011ea8aeSBarry Smith      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
498011ea8aeSBarry Smith 
499562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
500db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
501db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5024b9ad928SBarry Smith M*/
5034b9ad928SBarry Smith 
504d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
505d71ae5a4SJacob Faibussowitsch {
50669a612a9SBarry Smith   PetscMPIInt rank;
5074b9ad928SBarry Smith   PC_BJacobi *jac;
5084b9ad928SBarry Smith 
5094b9ad928SBarry Smith   PetscFunctionBegin;
5104dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
5119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5122fa5cd67SKarl Rupp 
5130a545947SLisandro Dalcin   pc->ops->apply             = NULL;
5147b6e2003SPierre Jolivet   pc->ops->matapply          = NULL;
5150a545947SLisandro Dalcin   pc->ops->applytranspose    = NULL;
5164dbf25a8SPierre Jolivet   pc->ops->matapplytranspose = NULL;
5174b9ad928SBarry Smith   pc->ops->setup             = PCSetUp_BJacobi;
5184b9ad928SBarry Smith   pc->ops->destroy           = PCDestroy_BJacobi;
5194b9ad928SBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_BJacobi;
5204b9ad928SBarry Smith   pc->ops->view              = PCView_BJacobi;
5210a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
5224b9ad928SBarry Smith 
5234b9ad928SBarry Smith   pc->data         = (void *)jac;
5244b9ad928SBarry Smith   jac->n           = -1;
5254b9ad928SBarry Smith   jac->n_local     = -1;
5264b9ad928SBarry Smith   jac->first_local = rank;
5270a545947SLisandro Dalcin   jac->ksp         = NULL;
5280a545947SLisandro Dalcin   jac->g_lens      = NULL;
5290a545947SLisandro Dalcin   jac->l_lens      = NULL;
5300a545947SLisandro Dalcin   jac->psubcomm    = NULL;
5314b9ad928SBarry Smith 
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5384b9ad928SBarry Smith }
5394b9ad928SBarry Smith 
5404b9ad928SBarry Smith /*
5414b9ad928SBarry Smith         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5424b9ad928SBarry Smith */
543d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
544d71ae5a4SJacob Faibussowitsch {
5454b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5464b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5474b9ad928SBarry Smith 
5484b9ad928SBarry Smith   PetscFunctionBegin;
5499566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553e91c6855SBarry Smith }
554e91c6855SBarry Smith 
555d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
556d71ae5a4SJacob Faibussowitsch {
557e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
558e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
559e91c6855SBarry Smith 
560e91c6855SBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5629566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5639566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5649566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5652e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5674b9ad928SBarry Smith }
5684b9ad928SBarry Smith 
569d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
570d71ae5a4SJacob Faibussowitsch {
5714b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5722295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
573539c167fSBarry Smith   KSPConvergedReason reason;
5744b9ad928SBarry Smith 
5754b9ad928SBarry Smith   PetscFunctionBegin;
5769566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5779566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
578ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5804b9ad928SBarry Smith }
5814b9ad928SBarry Smith 
582d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
583d71ae5a4SJacob Faibussowitsch {
5844b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5854b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5866e111a19SKarl Rupp 
5874b9ad928SBarry Smith   PetscFunctionBegin;
5889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
5899566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
590bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
591910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
592910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
5939566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
5949566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
5959566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
5969566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
5979566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
5983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5994b9ad928SBarry Smith }
6004b9ad928SBarry Smith 
6014dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_BJacobi_Singleblock_Private(PC pc, Mat X, Mat Y, PetscBool transpose)
602d71ae5a4SJacob Faibussowitsch {
6037b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6047b6e2003SPierre Jolivet   Mat         sX, sY;
6057b6e2003SPierre Jolivet 
6067b6e2003SPierre Jolivet   PetscFunctionBegin;
6077b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6087b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6097b6e2003SPierre Jolivet      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6109566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6119566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(X, &sX));
6129566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
613*78d05565SPierre Jolivet   if (!transpose) {
614*78d05565SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], sX, sY, 0));
615*78d05565SPierre Jolivet     PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
616*78d05565SPierre Jolivet     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], sX, sY, 0));
617*78d05565SPierre Jolivet   } else {
618*78d05565SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[0], sX, sY, 0));
619*78d05565SPierre Jolivet     PetscCall(KSPMatSolveTranspose(jac->ksp[0], sX, sY));
620*78d05565SPierre Jolivet     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[0], sX, sY, 0));
621*78d05565SPierre Jolivet   }
6224dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
6234dbf25a8SPierre Jolivet }
6244dbf25a8SPierre Jolivet 
6254dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
6264dbf25a8SPierre Jolivet {
6274dbf25a8SPierre Jolivet   PetscFunctionBegin;
6284dbf25a8SPierre Jolivet   PetscCall(PCMatApply_BJacobi_Singleblock_Private(pc, X, Y, PETSC_FALSE));
6294dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
6304dbf25a8SPierre Jolivet }
6314dbf25a8SPierre Jolivet 
6324dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
6334dbf25a8SPierre Jolivet {
6344dbf25a8SPierre Jolivet   PetscFunctionBegin;
6354dbf25a8SPierre Jolivet   PetscCall(PCMatApply_BJacobi_Singleblock_Private(pc, X, Y, PETSC_TRUE));
6363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6377b6e2003SPierre Jolivet }
6387b6e2003SPierre Jolivet 
639d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
640d71ae5a4SJacob Faibussowitsch {
6414b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6424b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
643d9ca1df4SBarry Smith   PetscScalar            *y_array;
644d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6454b9ad928SBarry Smith   PC                      subpc;
6464b9ad928SBarry Smith 
6474b9ad928SBarry Smith   PetscFunctionBegin;
6484b9ad928SBarry Smith   /*
6494b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6504b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6514b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6524b9ad928SBarry Smith     for the sequential solve.
6534b9ad928SBarry Smith   */
6549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6569566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6579566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6584b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
659c3f9dca2SPierre Jolivet   /* note this bypasses the inner KSP and its options completely */
6609566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6619566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6629566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6639566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6674b9ad928SBarry Smith }
6684b9ad928SBarry Smith 
669d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
670d71ae5a4SJacob Faibussowitsch {
6714b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6724b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
673d9ca1df4SBarry Smith   PetscScalar            *y_array;
674d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6754b9ad928SBarry Smith   PC                      subpc;
6764b9ad928SBarry Smith 
6774b9ad928SBarry Smith   PetscFunctionBegin;
6784b9ad928SBarry Smith   /*
6794b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6804b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6814b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6824b9ad928SBarry Smith     for the sequential solve.
6834b9ad928SBarry Smith   */
6849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6869566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6879566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
690c3f9dca2SPierre Jolivet   /* note this bypasses the inner KSP and its options completely */
6914b9ad928SBarry Smith 
6929566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6939566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6944b9ad928SBarry Smith 
69511e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->x));
69611e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->y));
6979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7004b9ad928SBarry Smith }
7014b9ad928SBarry Smith 
702d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
703d71ae5a4SJacob Faibussowitsch {
7044b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
7054b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
706d9ca1df4SBarry Smith   PetscScalar            *y_array;
707d9ca1df4SBarry Smith   const PetscScalar      *x_array;
7084b9ad928SBarry Smith 
7094b9ad928SBarry Smith   PetscFunctionBegin;
7104b9ad928SBarry Smith   /*
7114b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
7124b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
7134b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
7144b9ad928SBarry Smith     for the sequential solve.
7154b9ad928SBarry Smith   */
7169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
7179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
7189566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
7199566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
7209566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7219566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
7229566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
7239566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
7249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
7259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7274b9ad928SBarry Smith }
7284b9ad928SBarry Smith 
729d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
730d71ae5a4SJacob Faibussowitsch {
7314b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
73269a612a9SBarry Smith   PetscInt                m;
7334b9ad928SBarry Smith   KSP                     ksp;
7344b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
735de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7363f6dc190SJunchao Zhang   VecType                 vectype;
737ea41da7aSPierre Jolivet   const char             *prefix;
7384b9ad928SBarry Smith 
7394b9ad928SBarry Smith   PetscFunctionBegin;
7404b9ad928SBarry Smith   if (!pc->setupcalled) {
741c2efdce3SBarry Smith     if (!jac->ksp) {
7423821be0aSBarry Smith       PetscInt nestlevel;
7433821be0aSBarry Smith 
744302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7452fa5cd67SKarl Rupp 
7469566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7473821be0aSBarry Smith       PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
7483821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
7499566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7509566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7519566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7529566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7539566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7549566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7554b9ad928SBarry Smith 
756e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7574b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7584b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7597b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7604dbf25a8SPierre Jolivet       pc->ops->matapplytranspose   = PCMatApplyTranspose_BJacobi_Singleblock;
7614b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7624b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7634b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7644b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7654b9ad928SBarry Smith 
7669566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7674b9ad928SBarry Smith       jac->ksp[0] = ksp;
768c2efdce3SBarry Smith 
7694dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
7704b9ad928SBarry Smith       jac->data = (void *)bjac;
7714b9ad928SBarry Smith     } else {
772c2efdce3SBarry Smith       ksp  = jac->ksp[0];
773c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
774c2efdce3SBarry Smith     }
775c2efdce3SBarry Smith 
776c2efdce3SBarry Smith     /*
777c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
778c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
779c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
780c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
781c2efdce3SBarry Smith       KSPSolve() on the block.
782c2efdce3SBarry Smith     */
7839566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7849566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7859566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7869566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7879566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7889566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
789c2efdce3SBarry Smith   } else {
7904b9ad928SBarry Smith     ksp  = jac->ksp[0];
7914b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7924b9ad928SBarry Smith   }
7939566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
79449517cdeSBarry Smith   if (pc->useAmat) {
7959566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7969566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7974b9ad928SBarry Smith   } else {
7989566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7994b9ad928SBarry Smith   }
8009566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
80190f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
802248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
8039566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
804302a9ddcSMatthew Knepley   }
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8064b9ad928SBarry Smith }
8074b9ad928SBarry Smith 
808d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
809d71ae5a4SJacob Faibussowitsch {
8104b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
8114b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
81269a612a9SBarry Smith   PetscInt               i;
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith   PetscFunctionBegin;
815e91c6855SBarry Smith   if (bjac && bjac->pmat) {
8169566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
81748a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
818e91c6855SBarry Smith   }
8194b9ad928SBarry Smith 
8204b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
8219566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
822e91c6855SBarry Smith     if (bjac && bjac->x) {
8239566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
8249566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
8259566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
8264b9ad928SBarry Smith     }
827e91c6855SBarry Smith   }
8289566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
8299566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
831e91c6855SBarry Smith }
832e91c6855SBarry Smith 
833d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
834d71ae5a4SJacob Faibussowitsch {
835e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
836c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
837e91c6855SBarry Smith   PetscInt               i;
838e91c6855SBarry Smith 
839e91c6855SBarry Smith   PetscFunctionBegin;
8409566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
841c2efdce3SBarry Smith   if (bjac) {
8429566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8439566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8449566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
845c2efdce3SBarry Smith   }
8469566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
84748a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8489566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8492e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8514b9ad928SBarry Smith }
8524b9ad928SBarry Smith 
853d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
854d71ae5a4SJacob Faibussowitsch {
8554b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
85669a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
857539c167fSBarry Smith   KSPConvergedReason reason;
8584b9ad928SBarry Smith 
8594b9ad928SBarry Smith   PetscFunctionBegin;
8604b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8619566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8629566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
863ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8644b9ad928SBarry Smith   }
8653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8664b9ad928SBarry Smith }
8674b9ad928SBarry Smith 
868d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
869d71ae5a4SJacob Faibussowitsch {
8704b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
87169a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8724b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
873d9ca1df4SBarry Smith   PetscScalar           *yin;
874d9ca1df4SBarry Smith   const PetscScalar     *xin;
87558ebbce7SBarry Smith 
8764b9ad928SBarry Smith   PetscFunctionBegin;
8779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8794b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8804b9ad928SBarry Smith     /*
8814b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8824b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8834b9ad928SBarry Smith        the global vector.
8844b9ad928SBarry Smith     */
8859566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8869566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8874b9ad928SBarry Smith 
8889566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8899566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8909566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8919566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
892d11f3a42SBarry Smith 
8939566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8949566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8954b9ad928SBarry Smith   }
8969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8994b9ad928SBarry Smith }
9004b9ad928SBarry Smith 
901f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
902f4d694ddSBarry Smith {
903f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
904f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
905f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
906f4d694ddSBarry Smith   PetscScalar           *yin;
907f4d694ddSBarry Smith   const PetscScalar     *xin;
908f4d694ddSBarry Smith   PC                     subpc;
909f4d694ddSBarry Smith 
910f4d694ddSBarry Smith   PetscFunctionBegin;
911f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
912f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
913f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
9144b9ad928SBarry Smith     /*
915f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
916f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
917f4d694ddSBarry Smith        the global vector.
9184b9ad928SBarry Smith     */
919f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
920f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
921f4d694ddSBarry Smith 
922f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
923f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
924c3f9dca2SPierre Jolivet     /* note this bypasses the inner KSP and its options completely */
925f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
926f4d694ddSBarry Smith     PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
927f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
928f4d694ddSBarry Smith 
929f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
930f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
931f4d694ddSBarry Smith   }
932f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
933f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
935f4d694ddSBarry Smith }
936f4d694ddSBarry Smith 
937f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
938f4d694ddSBarry Smith {
939f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
940f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
941f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
942f4d694ddSBarry Smith   PetscScalar           *yin;
943f4d694ddSBarry Smith   const PetscScalar     *xin;
944f4d694ddSBarry Smith   PC                     subpc;
945f4d694ddSBarry Smith 
946f4d694ddSBarry Smith   PetscFunctionBegin;
947f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
948f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
949f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
950f4d694ddSBarry Smith     /*
951f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
952f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
953f4d694ddSBarry Smith        the global vector.
954f4d694ddSBarry Smith     */
955f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
956f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
957f4d694ddSBarry Smith 
958f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
959f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
960c3f9dca2SPierre Jolivet     /* note this bypasses the inner KSP and its options completely */
961f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
962f4d694ddSBarry Smith     PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
963f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
964f4d694ddSBarry Smith 
965f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
966f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
967f4d694ddSBarry Smith   }
968f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
969f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
971f4d694ddSBarry Smith }
972f4d694ddSBarry Smith 
973d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
974d71ae5a4SJacob Faibussowitsch {
9754b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
97669a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
9774b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
978d9ca1df4SBarry Smith   PetscScalar           *yin;
979d9ca1df4SBarry Smith   const PetscScalar     *xin;
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith   PetscFunctionBegin;
9829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
9839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
9844b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
9854b9ad928SBarry Smith     /*
9864b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
9874b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
9884b9ad928SBarry Smith        the global vector.
9894b9ad928SBarry Smith     */
9909566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9919566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9924b9ad928SBarry Smith 
9939566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
9949566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
9959566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
9969566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
997a7ff49e8SJed Brown 
9989566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
9999566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
10004b9ad928SBarry Smith   }
10019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
10029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
10033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10044b9ad928SBarry Smith }
10054b9ad928SBarry Smith 
1006d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
1007d71ae5a4SJacob Faibussowitsch {
10084b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
100969a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
1010ea41da7aSPierre Jolivet   const char            *prefix;
10114b9ad928SBarry Smith   KSP                    ksp;
10124b9ad928SBarry Smith   Vec                    x, y;
10134b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
10144b9ad928SBarry Smith   PC                     subpc;
10154b9ad928SBarry Smith   IS                     is;
1016434a97beSBrad Aagaard   MatReuse               scall;
10173f6dc190SJunchao Zhang   VecType                vectype;
10184849c82aSBarry Smith   MatNullSpace          *nullsp_mat = NULL, *nullsp_pmat = NULL;
10194b9ad928SBarry Smith 
10204b9ad928SBarry Smith   PetscFunctionBegin;
10219566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10224b9ad928SBarry Smith 
10234b9ad928SBarry Smith   n_local = jac->n_local;
10244b9ad928SBarry Smith 
102549517cdeSBarry Smith   if (pc->useAmat) {
1026ace3abfcSBarry Smith     PetscBool same;
10279566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
102828b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10294b9ad928SBarry Smith   }
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith   if (!pc->setupcalled) {
10323821be0aSBarry Smith     PetscInt nestlevel;
10333821be0aSBarry Smith 
10344b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
1035c2efdce3SBarry Smith 
1036c2efdce3SBarry Smith     if (!jac->ksp) {
1037e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Multiblock;
10384b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Multiblock;
10394b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Multiblock;
10407b6e2003SPierre Jolivet       pc->ops->matapply            = NULL;
10414dbf25a8SPierre Jolivet       pc->ops->matapplytranspose   = NULL;
1042f4d694ddSBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Multiblock;
1043f4d694ddSBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10444b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Multiblock;
10454b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Multiblock;
10464b9ad928SBarry Smith 
10474dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
10489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
10499566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
10514b9ad928SBarry Smith 
10524b9ad928SBarry Smith       jac->data = (void *)bjac;
10539566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
10544b9ad928SBarry Smith 
10554b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
10569566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10573821be0aSBarry Smith         PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
10583821be0aSBarry Smith         PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
10599566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10609566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10619566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
10629566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
10639566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
10649566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10659566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10662fa5cd67SKarl Rupp 
1067c2efdce3SBarry Smith         jac->ksp[i] = ksp;
1068c2efdce3SBarry Smith       }
1069c2efdce3SBarry Smith     } else {
1070c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
1071c2efdce3SBarry Smith     }
10724b9ad928SBarry Smith 
1073c2efdce3SBarry Smith     start = 0;
10749566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
1075c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
10764b9ad928SBarry Smith       m = jac->l_lens[i];
10774b9ad928SBarry Smith       /*
10784b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
10794b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
10804b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
10814b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
10824b9ad928SBarry Smith       KSPSolve() on the block.
10834b9ad928SBarry Smith 
10844b9ad928SBarry Smith       */
10859566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10869566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10879566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
10889566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
10892fa5cd67SKarl Rupp 
10904b9ad928SBarry Smith       bjac->x[i]      = x;
10914b9ad928SBarry Smith       bjac->y[i]      = y;
10924b9ad928SBarry Smith       bjac->starts[i] = start;
10934b9ad928SBarry Smith 
10949566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
10954b9ad928SBarry Smith       bjac->is[i] = is;
10964b9ad928SBarry Smith 
10974b9ad928SBarry Smith       start += m;
10984b9ad928SBarry Smith     }
10994b9ad928SBarry Smith   } else {
11004b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
11014b9ad928SBarry Smith     /*
11024b9ad928SBarry Smith        Destroy the blocks from the previous iteration
11034b9ad928SBarry Smith     */
11044b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
11054849c82aSBarry Smith       PetscCall(MatGetNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
11069566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
11074849c82aSBarry Smith       if (pc->useAmat) {
11084849c82aSBarry Smith         PetscCall(MatGetNullSpaces(n_local, bjac->mat, &nullsp_mat));
11094849c82aSBarry Smith         PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
11104849c82aSBarry Smith       }
11114b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
11122fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
11134b9ad928SBarry Smith   }
11144b9ad928SBarry Smith 
11159566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
11164849c82aSBarry Smith   if (nullsp_pmat) PetscCall(MatRestoreNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
11174849c82aSBarry Smith   if (pc->useAmat) {
11184849c82aSBarry Smith     PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
11194849c82aSBarry Smith     if (nullsp_mat) PetscCall(MatRestoreNullSpaces(n_local, bjac->mat, &nullsp_mat));
11204849c82aSBarry Smith   }
11214b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
11224b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
11239566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
11244b9ad928SBarry Smith 
11254b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
11269566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
112749517cdeSBarry Smith     if (pc->useAmat) {
11289566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11299566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11304b9ad928SBarry Smith     } else {
11319566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11324b9ad928SBarry Smith     }
11339566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
113448a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
113590f1c854SHong Zhang   }
11363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11374b9ad928SBarry Smith }
11385a7e1895SHong Zhang 
11395a7e1895SHong Zhang /*
1140fd0b8932SBarry Smith       These are for a single block with multiple processes
11415a7e1895SHong Zhang */
1142d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1143d71ae5a4SJacob Faibussowitsch {
1144fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1145fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1146fd0b8932SBarry Smith   KSPConvergedReason reason;
1147fd0b8932SBarry Smith 
1148fd0b8932SBarry Smith   PetscFunctionBegin;
11499566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
11509566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1151ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1153fd0b8932SBarry Smith }
1154fd0b8932SBarry Smith 
1155d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1156d71ae5a4SJacob Faibussowitsch {
11575a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11585a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11595a7e1895SHong Zhang 
11605a7e1895SHong Zhang   PetscFunctionBegin;
11619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
11629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
11639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
11649566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1166e91c6855SBarry Smith }
1167e91c6855SBarry Smith 
1168d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1169d71ae5a4SJacob Faibussowitsch {
1170e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1171e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1172e91c6855SBarry Smith 
1173e91c6855SBarry Smith   PetscFunctionBegin;
11749566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
11759566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
11769566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
11779566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11785a7e1895SHong Zhang 
11799566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
11802e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
11813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11825a7e1895SHong Zhang }
11835a7e1895SHong Zhang 
1184d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1185d71ae5a4SJacob Faibussowitsch {
11865a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11875a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1188d9ca1df4SBarry Smith   PetscScalar          *yarray;
1189d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1190539c167fSBarry Smith   KSPConvergedReason    reason;
11915a7e1895SHong Zhang 
11925a7e1895SHong Zhang   PetscFunctionBegin;
11935a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
11949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
11959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
11969566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
11979566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
11985a7e1895SHong Zhang 
11995a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
12009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
12019566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
12029566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
12039566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
12049566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1205ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1206250cf88bSHong Zhang 
12079566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
12089566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
12099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
12109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
12113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12125a7e1895SHong Zhang }
12135a7e1895SHong Zhang 
1214d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1215d71ae5a4SJacob Faibussowitsch {
12167b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
12177b6e2003SPierre Jolivet   KSPConvergedReason reason;
12187b6e2003SPierre Jolivet   Mat                sX, sY;
12197b6e2003SPierre Jolivet   const PetscScalar *x;
12207b6e2003SPierre Jolivet   PetscScalar       *y;
12217b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
12227b6e2003SPierre Jolivet 
12237b6e2003SPierre Jolivet   PetscFunctionBegin;
12247b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
12259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
12269566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
12279566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
12289566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
12299566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
12309566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
12319566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12329566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12339566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
12349566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
12359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12369566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12379566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12399566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
12409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
12419566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12429566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
12439566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1244ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12467b6e2003SPierre Jolivet }
12477b6e2003SPierre Jolivet 
1248d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1249d71ae5a4SJacob Faibussowitsch {
12505a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
12515a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12525a7e1895SHong Zhang   PetscInt              m, n;
1253ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
12545a7e1895SHong Zhang   const char           *prefix;
1255de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
12563f6dc190SJunchao Zhang   VecType               vectype;
12571f62890fSKarl Rupp 
12585a7e1895SHong Zhang   PetscFunctionBegin;
12599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
126008401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12615a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12625a7e1895SHong Zhang   if (!pc->setupcalled) {
12633821be0aSBarry Smith     PetscInt nestlevel;
12643821be0aSBarry Smith 
1265de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
12664dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mpjac));
12675a7e1895SHong Zhang     jac->data = (void *)mpjac;
12685a7e1895SHong Zhang 
12695a7e1895SHong Zhang     /* initialize datastructure mpjac */
12705a7e1895SHong Zhang     if (!jac->psubcomm) {
12715a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
12729566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12739566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12749566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12755a7e1895SHong Zhang     }
12765a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1277306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
12785a7e1895SHong Zhang 
12795a7e1895SHong Zhang     /* Get matrix blocks of pmat */
12809566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12815a7e1895SHong Zhang 
12825a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
12839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
12849566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12853821be0aSBarry Smith     PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
12863821be0aSBarry Smith     PetscCall(KSPSetNestLevel(jac->ksp[0], nestlevel + 1));
12879566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12889566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12899566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12909566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12915a7e1895SHong Zhang 
12929566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
12939566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
12949566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
12959566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
12969566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
12975a7e1895SHong Zhang 
12985a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
12999566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
13009566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
13019566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
13029566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
13039566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
13049566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
13055a7e1895SHong Zhang 
1306fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1307e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
13085a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
13095a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
13107b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1311fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1312306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
13139e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
13145a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
13159566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
13169566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1317fc08c53fSHong Zhang     } else {
13189566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
13195a7e1895SHong Zhang     }
13209566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
13215a7e1895SHong Zhang   }
13225a7e1895SHong Zhang 
132348a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
13243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13255a7e1895SHong Zhang }
1326