xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    Defines a block Jacobi preconditioner.
34b9ad928SBarry Smith */
400e125f8SBarry Smith 
500e125f8SBarry Smith #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/
64b9ad928SBarry Smith 
76849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat);
86849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat);
95a7e1895SHong Zhang static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
104b9ad928SBarry Smith 
11d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi(PC pc)
12d71ae5a4SJacob Faibussowitsch {
134b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
144b9ad928SBarry Smith   Mat         mat = pc->mat, pmat = pc->pmat;
15976e8c5aSLisandro Dalcin   PetscBool   hasop;
1669a612a9SBarry Smith   PetscInt    N, M, start, i, sum, end;
1769a612a9SBarry Smith   PetscInt    bs, i_start = -1, i_end = -1;
1869a612a9SBarry Smith   PetscMPIInt rank, size;
194b9ad928SBarry Smith 
204b9ad928SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
249566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(pc->pmat, &bs));
254b9ad928SBarry Smith 
265a7e1895SHong Zhang   if (jac->n > 0 && jac->n < size) {
279566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiproc(pc));
283ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
295a7e1895SHong Zhang   }
305a7e1895SHong Zhang 
31f1580f4eSBarry Smith   /*    Determines the number of blocks assigned to each processor */
324b9ad928SBarry Smith   /*   local block count  given */
334b9ad928SBarry Smith   if (jac->n_local > 0 && jac->n < 0) {
341c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
354b9ad928SBarry Smith     if (jac->l_lens) { /* check that user set these correctly */
364b9ad928SBarry Smith       sum = 0;
374b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
3808401ef6SPierre Jolivet         PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
394b9ad928SBarry Smith         sum += jac->l_lens[i];
404b9ad928SBarry Smith       }
4108401ef6SPierre Jolivet       PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly");
424b9ad928SBarry Smith     } else {
439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
442fa5cd67SKarl Rupp       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
454b9ad928SBarry Smith     }
464b9ad928SBarry Smith   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
474b9ad928SBarry Smith     /* global blocks given: determine which ones are local */
484b9ad928SBarry Smith     if (jac->g_lens) {
494b9ad928SBarry Smith       /* check if the g_lens is has valid entries */
504b9ad928SBarry Smith       for (i = 0; i < jac->n; i++) {
517827d75bSBarry Smith         PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed");
5208401ef6SPierre Jolivet         PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
534b9ad928SBarry Smith       }
544b9ad928SBarry Smith       if (size == 1) {
554b9ad928SBarry Smith         jac->n_local = jac->n;
569566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
579566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local));
584b9ad928SBarry Smith         /* check that user set these correctly */
594b9ad928SBarry Smith         sum = 0;
604b9ad928SBarry Smith         for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i];
6108401ef6SPierre Jolivet         PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly");
624b9ad928SBarry Smith       } else {
639566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end));
64aaa8cc7dSPierre Jolivet         /* loop over blocks determining first one owned by me */
654b9ad928SBarry Smith         sum = 0;
664b9ad928SBarry Smith         for (i = 0; i < jac->n + 1; i++) {
679371c9d4SSatish Balay           if (sum == start) {
689371c9d4SSatish Balay             i_start = i;
699371c9d4SSatish Balay             goto start_1;
709371c9d4SSatish Balay           }
714b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
724b9ad928SBarry Smith         }
73e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
744b9ad928SBarry Smith       start_1:
754b9ad928SBarry Smith         for (i = i_start; i < jac->n + 1; i++) {
769371c9d4SSatish Balay           if (sum == end) {
779371c9d4SSatish Balay             i_end = i;
789371c9d4SSatish Balay             goto end_1;
799371c9d4SSatish Balay           }
804b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
814b9ad928SBarry Smith         }
82e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
834b9ad928SBarry Smith       end_1:
844b9ad928SBarry Smith         jac->n_local = i_end - i_start;
859566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
869566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local));
874b9ad928SBarry Smith       }
884b9ad928SBarry Smith     } else { /* no global blocks given, determine then using default layout */
894b9ad928SBarry Smith       jac->n_local = jac->n / size + ((jac->n % size) > rank);
909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
914b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
924b9ad928SBarry Smith         jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs;
937827d75bSBarry Smith         PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given");
944b9ad928SBarry Smith       }
954b9ad928SBarry Smith     }
964b9ad928SBarry Smith   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
974b9ad928SBarry Smith     jac->n       = size;
984b9ad928SBarry Smith     jac->n_local = 1;
999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->l_lens));
1004b9ad928SBarry Smith     jac->l_lens[0] = M;
1017271b979SBarry Smith   } else { /* jac->n > 0 && jac->n_local > 0 */
1027271b979SBarry Smith     if (!jac->l_lens) {
1039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
1042fa5cd67SKarl Rupp       for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
1057271b979SBarry Smith     }
1064b9ad928SBarry Smith   }
10708401ef6SPierre Jolivet   PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors");
1084b9ad928SBarry Smith 
109f1580f4eSBarry Smith   /*    Determines mat and pmat */
1109566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
111976e8c5aSLisandro Dalcin   if (!hasop && size == 1) {
1124b9ad928SBarry Smith     mat  = pc->mat;
1134b9ad928SBarry Smith     pmat = pc->pmat;
1144b9ad928SBarry Smith   } else {
11549517cdeSBarry Smith     if (pc->useAmat) {
11649517cdeSBarry Smith       /* use block from Amat matrix, not Pmat for local MatMult() */
1179566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->mat, &mat));
1184b9ad928SBarry Smith     }
11949517cdeSBarry Smith     if (pc->pmat != pc->mat || !pc->useAmat) {
1209566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat));
1212fa5cd67SKarl Rupp     } else pmat = mat;
1224b9ad928SBarry Smith   }
1234b9ad928SBarry Smith 
124f1580f4eSBarry Smith   /*
1254b9ad928SBarry Smith      Setup code depends on the number of blocks
1264b9ad928SBarry Smith   */
127cc1d6079SHong Zhang   if (jac->n_local == 1) {
1289566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat));
1294b9ad928SBarry Smith   } else {
1309566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat));
1314b9ad928SBarry Smith   }
1323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith /* Default destroy, if it has never been setup */
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi(PC pc)
137d71ae5a4SJacob Faibussowitsch {
1384b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL));
1452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL));
1462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL));
1472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1504b9ad928SBarry Smith }
1514b9ad928SBarry Smith 
152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject)
153d71ae5a4SJacob Faibussowitsch {
1544b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
155248bfaf0SJed Brown   PetscInt    blocks, i;
156ace3abfcSBarry Smith   PetscBool   flg;
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   PetscFunctionBegin;
159d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
1619566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
1629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
1639566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
164248bfaf0SJed Brown   if (jac->ksp) {
165248bfaf0SJed Brown     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
166248bfaf0SJed Brown      * unless we had already been called. */
16748a46eb9SPierre Jolivet     for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
168248bfaf0SJed Brown   }
169d0609cedSBarry Smith   PetscOptionsHeadEnd();
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1714b9ad928SBarry Smith }
1724b9ad928SBarry Smith 
1739804daf3SBarry Smith #include <petscdraw.h>
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer)
175d71ae5a4SJacob Faibussowitsch {
1764b9ad928SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
177cbe18068SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
17869a612a9SBarry Smith   PetscMPIInt           rank;
17969a612a9SBarry Smith   PetscInt              i;
180d9884438SBarry Smith   PetscBool             iascii, isstring, isdraw;
1814b9ad928SBarry Smith   PetscViewer           sviewer;
182020d6619SPierre Jolivet   PetscViewerFormat     format;
183020d6619SPierre Jolivet   const char           *prefix;
1844b9ad928SBarry Smith 
1854b9ad928SBarry Smith   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
18932077d6dSBarry Smith   if (iascii) {
19048a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
19163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of blocks = %" PetscInt_FMT "\n", jac->n));
1929566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
1939566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
194020d6619SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
1969566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
1979566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
198933995ecSHong Zhang       if (jac->ksp && !jac->psubcomm) {
1999566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
200dd400576SPatrick Sanan         if (rank == 0) {
2019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2029566063dSJacob Faibussowitsch           PetscCall(KSPView(jac->ksp[0], sviewer));
2039566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2044b9ad928SBarry Smith         }
2059566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2069566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2079566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
208e4de9e1dSBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
210e4fa1014SBarry Smith       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
2119566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
212e4fa1014SBarry Smith         if (!mpjac->psubcomm->color) {
2139566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2149566063dSJacob Faibussowitsch           PetscCall(KSPView(*(jac->ksp), sviewer));
2159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
216cbe18068SHong Zhang         }
2179566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2189566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2199566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2209530cbd7SBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2224b9ad928SBarry Smith       } else {
2239566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
224e4de9e1dSBarry Smith       }
225e4de9e1dSBarry Smith     } else {
2264814766eSBarry Smith       PetscInt n_global;
2271c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2289566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
2309371c9d4SSatish Balay       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local));
2319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2329566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
233dd2fa690SBarry Smith       for (i = 0; i < jac->n_local; i++) {
23463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2359566063dSJacob Faibussowitsch         PetscCall(KSPView(jac->ksp[i], sviewer));
2369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2374b9ad928SBarry Smith       }
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2424b9ad928SBarry Smith     }
2434b9ad928SBarry Smith   } else if (isstring) {
24463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2459566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2469566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2479566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
248d9884438SBarry Smith   } else if (isdraw) {
249d9884438SBarry Smith     PetscDraw draw;
250d9884438SBarry Smith     char      str[25];
251d9884438SBarry Smith     PetscReal x, y, bottom, h;
252d9884438SBarry Smith 
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2549566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
25563a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2569566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
257d9884438SBarry Smith     bottom = y - h;
2589566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
259d9884438SBarry Smith     /* warning the communicator on viewer is different then on ksp in parallel */
2609566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2624b9ad928SBarry Smith   }
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2644b9ad928SBarry Smith }
2654b9ad928SBarry Smith 
266d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
267d71ae5a4SJacob Faibussowitsch {
268feb237baSPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith   PetscFunctionBegin;
27128b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2724b9ad928SBarry Smith 
2734b9ad928SBarry Smith   if (n_local) *n_local = jac->n_local;
2744b9ad928SBarry Smith   if (first_local) *first_local = jac->first_local;
275020d6619SPierre Jolivet   if (ksp) *ksp = jac->ksp;
2763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2774b9ad928SBarry Smith }
2784b9ad928SBarry Smith 
279f9663b93SPierre Jolivet static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt *lens)
280d71ae5a4SJacob Faibussowitsch {
2814b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith   PetscFunctionBegin;
2842472a847SBarry 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");
2854b9ad928SBarry Smith   jac->n = blocks;
2860a545947SLisandro Dalcin   if (!lens) jac->g_lens = NULL;
2872fa5cd67SKarl Rupp   else {
2889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2904b9ad928SBarry Smith   }
2913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2924b9ad928SBarry Smith }
2934b9ad928SBarry Smith 
294d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
295d71ae5a4SJacob Faibussowitsch {
2964b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith   PetscFunctionBegin;
2994b9ad928SBarry Smith   *blocks = jac->n;
3004b9ad928SBarry Smith   if (lens) *lens = jac->g_lens;
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3024b9ad928SBarry Smith }
3034b9ad928SBarry Smith 
304d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[])
305d71ae5a4SJacob Faibussowitsch {
3064b9ad928SBarry Smith   PC_BJacobi *jac;
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith   PetscFunctionBegin;
3094b9ad928SBarry Smith   jac = (PC_BJacobi *)pc->data;
3104b9ad928SBarry Smith 
3114b9ad928SBarry Smith   jac->n_local = blocks;
3120a545947SLisandro Dalcin   if (!lens) jac->l_lens = NULL;
3132fa5cd67SKarl Rupp   else {
3149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3159566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3164b9ad928SBarry Smith   }
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3184b9ad928SBarry Smith }
3194b9ad928SBarry Smith 
320d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
321d71ae5a4SJacob Faibussowitsch {
3224b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
3254b9ad928SBarry Smith   *blocks = jac->n_local;
3264b9ad928SBarry Smith   if (lens) *lens = jac->l_lens;
3273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3284b9ad928SBarry Smith }
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith /*@C
331f1580f4eSBarry Smith   PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3324b9ad928SBarry Smith   this processor.
3334b9ad928SBarry Smith 
3346da92b7fSHong Zhang   Not Collective
3354b9ad928SBarry Smith 
3364b9ad928SBarry Smith   Input Parameter:
3374b9ad928SBarry Smith . pc - the preconditioner context
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith   Output Parameters:
3400298fd71SBarry Smith + n_local     - the number of blocks on this processor, or NULL
3410298fd71SBarry Smith . first_local - the global number of the first block on this processor, or NULL
3424b9ad928SBarry Smith - ksp         - the array of KSP contexts
3434b9ad928SBarry Smith 
3444b9ad928SBarry Smith   Notes:
345f1580f4eSBarry Smith   After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3464b9ad928SBarry Smith 
3474b9ad928SBarry Smith   Currently for some matrix implementations only 1 block per processor
3484b9ad928SBarry Smith   is supported.
3494b9ad928SBarry Smith 
350f1580f4eSBarry Smith   You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3514b9ad928SBarry Smith 
352feefa0e1SJacob Faibussowitsch   Fortran Notes:
353f1580f4eSBarry Smith   You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
354f1580f4eSBarry Smith 
355f1580f4eSBarry Smith   You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
356f1580f4eSBarry Smith   `KSP` array must be.
357196cc216SBarry Smith 
3584b9ad928SBarry Smith   Level: advanced
3594b9ad928SBarry Smith 
360*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3614b9ad928SBarry Smith @*/
362d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
363d71ae5a4SJacob Faibussowitsch {
3644b9ad928SBarry Smith   PetscFunctionBegin;
3650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
366cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3684b9ad928SBarry Smith }
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith /*@
3714b9ad928SBarry Smith   PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3724b9ad928SBarry Smith   Jacobi preconditioner.
3734b9ad928SBarry Smith 
374c3339decSBarry Smith   Collective
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith   Input Parameters:
3774b9ad928SBarry Smith + pc     - the preconditioner context
3784b9ad928SBarry Smith . blocks - the number of blocks
3794b9ad928SBarry Smith - lens   - [optional] integer array containing the size of each block
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith   Options Database Key:
3824b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3834b9ad928SBarry Smith 
384f1580f4eSBarry Smith   Note:
3854b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
386f1580f4eSBarry Smith   All processors sharing the `PC` must call this routine with the same data.
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith   Level: intermediate
3894b9ad928SBarry Smith 
390*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3914b9ad928SBarry Smith @*/
392d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
393d71ae5a4SJacob Faibussowitsch {
3944b9ad928SBarry Smith   PetscFunctionBegin;
3950700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39608401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
397cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3994b9ad928SBarry Smith }
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith /*@C
4024b9ad928SBarry Smith   PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
403f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4044b9ad928SBarry Smith 
405ad4df100SBarry Smith   Not Collective
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith   Input Parameter:
4084b9ad928SBarry Smith . pc - the preconditioner context
4094b9ad928SBarry Smith 
410feefa0e1SJacob Faibussowitsch   Output Parameters:
4114b9ad928SBarry Smith + blocks - the number of blocks
4124b9ad928SBarry Smith - lens   - integer array containing the size of each block
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith   Level: intermediate
4154b9ad928SBarry Smith 
416*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4174b9ad928SBarry Smith @*/
418d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
419d71ae5a4SJacob Faibussowitsch {
4204b9ad928SBarry Smith   PetscFunctionBegin;
4210700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4224f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
423cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4254b9ad928SBarry Smith }
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith /*@
4284b9ad928SBarry Smith   PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
429f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`,  preconditioner.
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith   Not Collective
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith   Input Parameters:
4344b9ad928SBarry Smith + pc     - the preconditioner context
4354b9ad928SBarry Smith . blocks - the number of blocks
4364b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4374b9ad928SBarry Smith 
438342c94f9SMatthew G. Knepley   Options Database Key:
439342c94f9SMatthew G. Knepley . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
440342c94f9SMatthew G. Knepley 
4414b9ad928SBarry Smith   Note:
4424b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith   Level: intermediate
4454b9ad928SBarry Smith 
446*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4474b9ad928SBarry Smith @*/
448d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
449d71ae5a4SJacob Faibussowitsch {
4504b9ad928SBarry Smith   PetscFunctionBegin;
4510700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
45208401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
453cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4554b9ad928SBarry Smith }
4564b9ad928SBarry Smith 
4574b9ad928SBarry Smith /*@C
4584b9ad928SBarry Smith   PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
459f1580f4eSBarry Smith   Jacobi, `PCBJACOBI`, preconditioner.
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith   Not Collective
4624b9ad928SBarry Smith 
4634b9ad928SBarry Smith   Input Parameters:
4644b9ad928SBarry Smith + pc     - the preconditioner context
4654b9ad928SBarry Smith . blocks - the number of blocks
4664b9ad928SBarry Smith - lens   - [optional] integer array containing size of each block
4674b9ad928SBarry Smith 
4684b9ad928SBarry Smith   Note:
4694b9ad928SBarry Smith   Currently only a limited number of blocking configurations are supported.
4704b9ad928SBarry Smith 
4714b9ad928SBarry Smith   Level: intermediate
4724b9ad928SBarry Smith 
473*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4744b9ad928SBarry Smith @*/
475d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
476d71ae5a4SJacob Faibussowitsch {
4774b9ad928SBarry Smith   PetscFunctionBegin;
4780700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4794f572ea9SToby Isaac   PetscAssertPointer(blocks, 2);
480cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4824b9ad928SBarry Smith }
4834b9ad928SBarry Smith 
4844b9ad928SBarry Smith /*MC
4854b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
486f1580f4eSBarry Smith            its own `KSP` object.
4874b9ad928SBarry Smith 
4884b9ad928SBarry Smith    Options Database Keys:
489011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
490011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4914b9ad928SBarry Smith 
49295452b02SPatrick Sanan    Notes:
493f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
494468ae2e8SBarry Smith 
49595452b02SPatrick 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.
4964b9ad928SBarry Smith 
497f1580f4eSBarry Smith      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
498d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
4994b9ad928SBarry Smith 
500f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
501f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
5020462cc06SPierre Jolivet          `KSPGetPC()`)
5034b9ad928SBarry Smith 
504f1580f4eSBarry Smith      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
5052210c163SDominic Meiser          performance.  Different block partitioning may lead to additional data transfers
5062210c163SDominic Meiser          between host and GPU that lead to degraded performance.
5072210c163SDominic Meiser 
508011ea8aeSBarry Smith      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
509011ea8aeSBarry Smith 
5104b9ad928SBarry Smith    Level: beginner
5114b9ad928SBarry Smith 
512*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
513db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
514db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5154b9ad928SBarry Smith M*/
5164b9ad928SBarry Smith 
517d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
518d71ae5a4SJacob Faibussowitsch {
51969a612a9SBarry Smith   PetscMPIInt rank;
5204b9ad928SBarry Smith   PC_BJacobi *jac;
5214b9ad928SBarry Smith 
5224b9ad928SBarry Smith   PetscFunctionBegin;
5234dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
5249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5252fa5cd67SKarl Rupp 
5260a545947SLisandro Dalcin   pc->ops->apply           = NULL;
5277b6e2003SPierre Jolivet   pc->ops->matapply        = NULL;
5280a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
5294b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_BJacobi;
5304b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_BJacobi;
5314b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
5324b9ad928SBarry Smith   pc->ops->view            = PCView_BJacobi;
5330a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
5344b9ad928SBarry Smith 
5354b9ad928SBarry Smith   pc->data         = (void *)jac;
5364b9ad928SBarry Smith   jac->n           = -1;
5374b9ad928SBarry Smith   jac->n_local     = -1;
5384b9ad928SBarry Smith   jac->first_local = rank;
5390a545947SLisandro Dalcin   jac->ksp         = NULL;
5400a545947SLisandro Dalcin   jac->g_lens      = NULL;
5410a545947SLisandro Dalcin   jac->l_lens      = NULL;
5420a545947SLisandro Dalcin   jac->psubcomm    = NULL;
5434b9ad928SBarry Smith 
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5504b9ad928SBarry Smith }
5514b9ad928SBarry Smith 
5524b9ad928SBarry Smith /*
5534b9ad928SBarry Smith         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5544b9ad928SBarry Smith */
555d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
556d71ae5a4SJacob Faibussowitsch {
5574b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5584b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5594b9ad928SBarry Smith 
5604b9ad928SBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
5643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
565e91c6855SBarry Smith }
566e91c6855SBarry Smith 
567d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
568d71ae5a4SJacob Faibussowitsch {
569e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
570e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
571e91c6855SBarry Smith 
572e91c6855SBarry Smith   PetscFunctionBegin;
5739566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5749566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5759566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5769566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5772e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5794b9ad928SBarry Smith }
5804b9ad928SBarry Smith 
581d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
582d71ae5a4SJacob Faibussowitsch {
5834b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5842295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
585539c167fSBarry Smith   KSPConvergedReason reason;
5864b9ad928SBarry Smith 
5874b9ad928SBarry Smith   PetscFunctionBegin;
5889566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5899566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
590ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5924b9ad928SBarry Smith }
5934b9ad928SBarry Smith 
594d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
595d71ae5a4SJacob Faibussowitsch {
5964b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5974b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5986e111a19SKarl Rupp 
5994b9ad928SBarry Smith   PetscFunctionBegin;
6009566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
6019566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
602bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
603910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
604910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6059566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6069566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
6079566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
6089566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
6099566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
6103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6114b9ad928SBarry Smith }
6124b9ad928SBarry Smith 
613d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
614d71ae5a4SJacob Faibussowitsch {
6157b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6167b6e2003SPierre Jolivet   Mat         sX, sY;
6177b6e2003SPierre Jolivet 
6187b6e2003SPierre Jolivet   PetscFunctionBegin;
6197b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6207b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6217b6e2003SPierre Jolivet      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6229566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6239566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(X, &sX));
6249566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
6259566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6277b6e2003SPierre Jolivet }
6287b6e2003SPierre Jolivet 
629d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
630d71ae5a4SJacob Faibussowitsch {
6314b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6324b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
633d9ca1df4SBarry Smith   PetscScalar            *y_array;
634d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6354b9ad928SBarry Smith   PC                      subpc;
6364b9ad928SBarry Smith 
6374b9ad928SBarry Smith   PetscFunctionBegin;
6384b9ad928SBarry Smith   /*
6394b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6404b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6414b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6424b9ad928SBarry Smith     for the sequential solve.
6434b9ad928SBarry Smith   */
6449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6469566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6479566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6484b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
6494b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6509566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6519566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6529566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6539566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6574b9ad928SBarry Smith }
6584b9ad928SBarry Smith 
659d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
660d71ae5a4SJacob Faibussowitsch {
6614b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6624b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
663d9ca1df4SBarry Smith   PetscScalar            *y_array;
664d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6654b9ad928SBarry Smith   PC                      subpc;
6664b9ad928SBarry Smith 
6674b9ad928SBarry Smith   PetscFunctionBegin;
6684b9ad928SBarry Smith   /*
6694b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6704b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6714b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6724b9ad928SBarry Smith     for the sequential solve.
6734b9ad928SBarry Smith   */
6749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6769566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6779566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6784b9ad928SBarry Smith 
6794b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
6804b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6814b9ad928SBarry Smith 
6829566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6839566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6844b9ad928SBarry Smith 
68511e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->x));
68611e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->y));
6879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6904b9ad928SBarry Smith }
6914b9ad928SBarry Smith 
692d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
693d71ae5a4SJacob Faibussowitsch {
6944b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6954b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
696d9ca1df4SBarry Smith   PetscScalar            *y_array;
697d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6984b9ad928SBarry Smith 
6994b9ad928SBarry Smith   PetscFunctionBegin;
7004b9ad928SBarry Smith   /*
7014b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
7024b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
7034b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
7044b9ad928SBarry Smith     for the sequential solve.
7054b9ad928SBarry Smith   */
7069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
7079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
7089566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
7099566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
7109566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7119566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
7129566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
7139566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
7149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
7159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
7163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7174b9ad928SBarry Smith }
7184b9ad928SBarry Smith 
719d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
720d71ae5a4SJacob Faibussowitsch {
7214b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
72269a612a9SBarry Smith   PetscInt                m;
7234b9ad928SBarry Smith   KSP                     ksp;
7244b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
725de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7263f6dc190SJunchao Zhang   VecType                 vectype;
727ea41da7aSPierre Jolivet   const char             *prefix;
7284b9ad928SBarry Smith 
7294b9ad928SBarry Smith   PetscFunctionBegin;
7304b9ad928SBarry Smith   if (!pc->setupcalled) {
731c2efdce3SBarry Smith     if (!jac->ksp) {
7323821be0aSBarry Smith       PetscInt nestlevel;
7333821be0aSBarry Smith 
734302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7352fa5cd67SKarl Rupp 
7369566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7373821be0aSBarry Smith       PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
7383821be0aSBarry Smith       PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
7399566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7409566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7419566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7429566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7439566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7449566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7454b9ad928SBarry Smith 
746e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7474b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7484b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7497b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7504b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7514b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7524b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7534b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7544b9ad928SBarry Smith 
7559566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7564b9ad928SBarry Smith       jac->ksp[0] = ksp;
757c2efdce3SBarry Smith 
7584dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
7594b9ad928SBarry Smith       jac->data = (void *)bjac;
7604b9ad928SBarry Smith     } else {
761c2efdce3SBarry Smith       ksp  = jac->ksp[0];
762c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
763c2efdce3SBarry Smith     }
764c2efdce3SBarry Smith 
765c2efdce3SBarry Smith     /*
766c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
767c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
768c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
769c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
770c2efdce3SBarry Smith       KSPSolve() on the block.
771c2efdce3SBarry Smith     */
7729566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7739566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7749566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7759566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7769566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7779566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
778c2efdce3SBarry Smith   } else {
7794b9ad928SBarry Smith     ksp  = jac->ksp[0];
7804b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7814b9ad928SBarry Smith   }
7829566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
78349517cdeSBarry Smith   if (pc->useAmat) {
7849566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7859566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7864b9ad928SBarry Smith   } else {
7879566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7884b9ad928SBarry Smith   }
7899566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
79090f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
791248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
7929566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
793302a9ddcSMatthew Knepley   }
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7954b9ad928SBarry Smith }
7964b9ad928SBarry Smith 
797d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
798d71ae5a4SJacob Faibussowitsch {
7994b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
8004b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
80169a612a9SBarry Smith   PetscInt               i;
8024b9ad928SBarry Smith 
8034b9ad928SBarry Smith   PetscFunctionBegin;
804e91c6855SBarry Smith   if (bjac && bjac->pmat) {
8059566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
80648a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
807e91c6855SBarry Smith   }
8084b9ad928SBarry Smith 
8094b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
8109566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
811e91c6855SBarry Smith     if (bjac && bjac->x) {
8129566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
8139566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
8149566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
8154b9ad928SBarry Smith     }
816e91c6855SBarry Smith   }
8179566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
8189566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
8193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
820e91c6855SBarry Smith }
821e91c6855SBarry Smith 
822d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
823d71ae5a4SJacob Faibussowitsch {
824e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
825c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
826e91c6855SBarry Smith   PetscInt               i;
827e91c6855SBarry Smith 
828e91c6855SBarry Smith   PetscFunctionBegin;
8299566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
830c2efdce3SBarry Smith   if (bjac) {
8319566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8329566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8339566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
834c2efdce3SBarry Smith   }
8359566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
83648a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8379566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8382e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8404b9ad928SBarry Smith }
8414b9ad928SBarry Smith 
842d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
843d71ae5a4SJacob Faibussowitsch {
8444b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
84569a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
846539c167fSBarry Smith   KSPConvergedReason reason;
8474b9ad928SBarry Smith 
8484b9ad928SBarry Smith   PetscFunctionBegin;
8494b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8509566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8519566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
852ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8534b9ad928SBarry Smith   }
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8554b9ad928SBarry Smith }
8564b9ad928SBarry Smith 
857d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
858d71ae5a4SJacob Faibussowitsch {
8594b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
86069a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8614b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
862d9ca1df4SBarry Smith   PetscScalar           *yin;
863d9ca1df4SBarry Smith   const PetscScalar     *xin;
86458ebbce7SBarry Smith 
8654b9ad928SBarry Smith   PetscFunctionBegin;
8669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8684b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8694b9ad928SBarry Smith     /*
8704b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8714b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8724b9ad928SBarry Smith        the global vector.
8734b9ad928SBarry Smith     */
8749566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8759566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8764b9ad928SBarry Smith 
8779566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8789566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8799566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8809566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
881d11f3a42SBarry Smith 
8829566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8839566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8844b9ad928SBarry Smith   }
8859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8884b9ad928SBarry Smith }
8894b9ad928SBarry Smith 
890f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
891f4d694ddSBarry Smith {
892f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
893f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
894f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
895f4d694ddSBarry Smith   PetscScalar           *yin;
896f4d694ddSBarry Smith   const PetscScalar     *xin;
897f4d694ddSBarry Smith   PC                     subpc;
898f4d694ddSBarry Smith 
899f4d694ddSBarry Smith   PetscFunctionBegin;
900f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
901f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
902f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
9034b9ad928SBarry Smith     /*
904f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
905f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
906f4d694ddSBarry Smith        the global vector.
9074b9ad928SBarry Smith     */
908f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
909f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
910f4d694ddSBarry Smith 
911f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
912f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
913f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
914f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
915f4d694ddSBarry Smith     PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
916f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
917f4d694ddSBarry Smith 
918f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
919f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
920f4d694ddSBarry Smith   }
921f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
922f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
924f4d694ddSBarry Smith }
925f4d694ddSBarry Smith 
926f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
927f4d694ddSBarry Smith {
928f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
929f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
930f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
931f4d694ddSBarry Smith   PetscScalar           *yin;
932f4d694ddSBarry Smith   const PetscScalar     *xin;
933f4d694ddSBarry Smith   PC                     subpc;
934f4d694ddSBarry Smith 
935f4d694ddSBarry Smith   PetscFunctionBegin;
936f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
937f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
938f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
939f4d694ddSBarry Smith     /*
940f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
941f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
942f4d694ddSBarry Smith        the global vector.
943f4d694ddSBarry Smith     */
944f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
945f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
946f4d694ddSBarry Smith 
947f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
948f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
949f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
950f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
951f4d694ddSBarry Smith     PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
952f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
953f4d694ddSBarry Smith 
954f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
955f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
956f4d694ddSBarry Smith   }
957f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
958f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
960f4d694ddSBarry Smith }
961f4d694ddSBarry Smith 
962d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
963d71ae5a4SJacob Faibussowitsch {
9644b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
96569a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
9664b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
967d9ca1df4SBarry Smith   PetscScalar           *yin;
968d9ca1df4SBarry Smith   const PetscScalar     *xin;
9694b9ad928SBarry Smith 
9704b9ad928SBarry Smith   PetscFunctionBegin;
9719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
9729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
9734b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
9744b9ad928SBarry Smith     /*
9754b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
9764b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
9774b9ad928SBarry Smith        the global vector.
9784b9ad928SBarry Smith     */
9799566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9809566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9814b9ad928SBarry Smith 
9829566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
9839566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
9849566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
9859566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
986a7ff49e8SJed Brown 
9879566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
9889566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
9894b9ad928SBarry Smith   }
9909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
9919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
9923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9934b9ad928SBarry Smith }
9944b9ad928SBarry Smith 
995d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
996d71ae5a4SJacob Faibussowitsch {
9974b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
99869a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
999ea41da7aSPierre Jolivet   const char            *prefix;
10004b9ad928SBarry Smith   KSP                    ksp;
10014b9ad928SBarry Smith   Vec                    x, y;
10024b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
10034b9ad928SBarry Smith   PC                     subpc;
10044b9ad928SBarry Smith   IS                     is;
1005434a97beSBrad Aagaard   MatReuse               scall;
10063f6dc190SJunchao Zhang   VecType                vectype;
10074b9ad928SBarry Smith 
10084b9ad928SBarry Smith   PetscFunctionBegin;
10099566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10104b9ad928SBarry Smith 
10114b9ad928SBarry Smith   n_local = jac->n_local;
10124b9ad928SBarry Smith 
101349517cdeSBarry Smith   if (pc->useAmat) {
1014ace3abfcSBarry Smith     PetscBool same;
10159566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
101628b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10174b9ad928SBarry Smith   }
10184b9ad928SBarry Smith 
10194b9ad928SBarry Smith   if (!pc->setupcalled) {
10203821be0aSBarry Smith     PetscInt nestlevel;
10213821be0aSBarry Smith 
10224b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
1023c2efdce3SBarry Smith 
1024c2efdce3SBarry Smith     if (!jac->ksp) {
1025e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Multiblock;
10264b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Multiblock;
10274b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Multiblock;
10287b6e2003SPierre Jolivet       pc->ops->matapply            = NULL;
1029f4d694ddSBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Multiblock;
1030f4d694ddSBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10314b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Multiblock;
10324b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Multiblock;
10334b9ad928SBarry Smith 
10344dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
10359566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
10369566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10379566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
10384b9ad928SBarry Smith 
10394b9ad928SBarry Smith       jac->data = (void *)bjac;
10409566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
10414b9ad928SBarry Smith 
10424b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
10439566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10443821be0aSBarry Smith         PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
10453821be0aSBarry Smith         PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
10469566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10479566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10489566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
10499566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
10509566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
10519566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10529566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10532fa5cd67SKarl Rupp 
1054c2efdce3SBarry Smith         jac->ksp[i] = ksp;
1055c2efdce3SBarry Smith       }
1056c2efdce3SBarry Smith     } else {
1057c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
1058c2efdce3SBarry Smith     }
10594b9ad928SBarry Smith 
1060c2efdce3SBarry Smith     start = 0;
10619566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
1062c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
10634b9ad928SBarry Smith       m = jac->l_lens[i];
10644b9ad928SBarry Smith       /*
10654b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
10664b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
10674b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
10684b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
10694b9ad928SBarry Smith       KSPSolve() on the block.
10704b9ad928SBarry Smith 
10714b9ad928SBarry Smith       */
10729566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10739566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10749566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
10759566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
10762fa5cd67SKarl Rupp 
10774b9ad928SBarry Smith       bjac->x[i]      = x;
10784b9ad928SBarry Smith       bjac->y[i]      = y;
10794b9ad928SBarry Smith       bjac->starts[i] = start;
10804b9ad928SBarry Smith 
10819566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
10824b9ad928SBarry Smith       bjac->is[i] = is;
10834b9ad928SBarry Smith 
10844b9ad928SBarry Smith       start += m;
10854b9ad928SBarry Smith     }
10864b9ad928SBarry Smith   } else {
10874b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
10884b9ad928SBarry Smith     /*
10894b9ad928SBarry Smith        Destroy the blocks from the previous iteration
10904b9ad928SBarry Smith     */
10914b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
10929566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
109348a46eb9SPierre Jolivet       if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
10944b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
10952fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
10964b9ad928SBarry Smith   }
10974b9ad928SBarry Smith 
10989566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
109948a46eb9SPierre Jolivet   if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
11004b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
11014b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
11029566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
11034b9ad928SBarry Smith 
11044b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
11059566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
110649517cdeSBarry Smith     if (pc->useAmat) {
11079566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11089566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11094b9ad928SBarry Smith     } else {
11109566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11114b9ad928SBarry Smith     }
11129566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
111348a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
111490f1c854SHong Zhang   }
11153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11164b9ad928SBarry Smith }
11175a7e1895SHong Zhang 
11185a7e1895SHong Zhang /*
1119fd0b8932SBarry Smith       These are for a single block with multiple processes
11205a7e1895SHong Zhang */
1121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1122d71ae5a4SJacob Faibussowitsch {
1123fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1124fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1125fd0b8932SBarry Smith   KSPConvergedReason reason;
1126fd0b8932SBarry Smith 
1127fd0b8932SBarry Smith   PetscFunctionBegin;
11289566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
11299566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1130ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1132fd0b8932SBarry Smith }
1133fd0b8932SBarry Smith 
1134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1135d71ae5a4SJacob Faibussowitsch {
11365a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11375a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11385a7e1895SHong Zhang 
11395a7e1895SHong Zhang   PetscFunctionBegin;
11409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
11419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
11429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
11439566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1145e91c6855SBarry Smith }
1146e91c6855SBarry Smith 
1147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1148d71ae5a4SJacob Faibussowitsch {
1149e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1150e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1151e91c6855SBarry Smith 
1152e91c6855SBarry Smith   PetscFunctionBegin;
11539566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
11549566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
11559566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
11569566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11575a7e1895SHong Zhang 
11589566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
11592e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
11603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11615a7e1895SHong Zhang }
11625a7e1895SHong Zhang 
1163d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1164d71ae5a4SJacob Faibussowitsch {
11655a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11665a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1167d9ca1df4SBarry Smith   PetscScalar          *yarray;
1168d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1169539c167fSBarry Smith   KSPConvergedReason    reason;
11705a7e1895SHong Zhang 
11715a7e1895SHong Zhang   PetscFunctionBegin;
11725a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
11739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
11749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
11759566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
11769566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
11775a7e1895SHong Zhang 
11785a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
11799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11809566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
11819566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
11829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11839566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1184ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1185250cf88bSHong Zhang 
11869566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
11879566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
11889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
11899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
11903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11915a7e1895SHong Zhang }
11925a7e1895SHong Zhang 
1193d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1194d71ae5a4SJacob Faibussowitsch {
11957b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
11967b6e2003SPierre Jolivet   KSPConvergedReason reason;
11977b6e2003SPierre Jolivet   Mat                sX, sY;
11987b6e2003SPierre Jolivet   const PetscScalar *x;
11997b6e2003SPierre Jolivet   PetscScalar       *y;
12007b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
12017b6e2003SPierre Jolivet 
12027b6e2003SPierre Jolivet   PetscFunctionBegin;
12037b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
12049566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
12059566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
12069566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
12079566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
12089566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
12099566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
12109566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12119566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12129566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
12139566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
12149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12159566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12169566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
12199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
12209566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12219566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
12229566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1223ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12257b6e2003SPierre Jolivet }
12267b6e2003SPierre Jolivet 
1227d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1228d71ae5a4SJacob Faibussowitsch {
12295a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
12305a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12315a7e1895SHong Zhang   PetscInt              m, n;
1232ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
12335a7e1895SHong Zhang   const char           *prefix;
1234de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
12353f6dc190SJunchao Zhang   VecType               vectype;
12361f62890fSKarl Rupp 
12375a7e1895SHong Zhang   PetscFunctionBegin;
12389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
123908401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12405a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12415a7e1895SHong Zhang   if (!pc->setupcalled) {
12423821be0aSBarry Smith     PetscInt nestlevel;
12433821be0aSBarry Smith 
1244de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
12454dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mpjac));
12465a7e1895SHong Zhang     jac->data = (void *)mpjac;
12475a7e1895SHong Zhang 
12485a7e1895SHong Zhang     /* initialize datastructure mpjac */
12495a7e1895SHong Zhang     if (!jac->psubcomm) {
12505a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
12519566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12529566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12539566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12545a7e1895SHong Zhang     }
12555a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1256306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
12575a7e1895SHong Zhang 
12585a7e1895SHong Zhang     /* Get matrix blocks of pmat */
12599566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12605a7e1895SHong Zhang 
12615a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
12629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
12639566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12643821be0aSBarry Smith     PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
12653821be0aSBarry Smith     PetscCall(KSPSetNestLevel(jac->ksp[0], nestlevel + 1));
12669566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12679566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12689566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12699566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12705a7e1895SHong Zhang 
12719566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
12729566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
12739566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
12749566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
12759566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
12765a7e1895SHong Zhang 
12775a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
12789566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
12799566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
12809566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
12819566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
12829566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
12839566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
12845a7e1895SHong Zhang 
1285fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1286e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
12875a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
12885a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
12897b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1290fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1291306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
12929e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
12935a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
12949566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
12959566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1296fc08c53fSHong Zhang     } else {
12979566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
12985a7e1895SHong Zhang     }
12999566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
13005a7e1895SHong Zhang   }
13015a7e1895SHong Zhang 
130248a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
13033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13045a7e1895SHong Zhang }
1305