xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith    Defines a block Jacobi preconditioner.
44b9ad928SBarry Smith */
500e125f8SBarry Smith 
600e125f8SBarry Smith #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/
74b9ad928SBarry Smith 
86849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat);
96849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat);
105a7e1895SHong Zhang static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
114b9ad928SBarry Smith 
129371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi(PC pc) {
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));
285a7e1895SHong Zhang     PetscFunctionReturn(0);
295a7e1895SHong Zhang   }
305a7e1895SHong Zhang 
31*f1580f4eSBarry 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));
644b9ad928SBarry Smith         /* loop over blocks determing 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 
109*f1580f4eSBarry 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 
124*f1580f4eSBarry 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   }
1324b9ad928SBarry Smith   PetscFunctionReturn(0);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith /* Default destroy, if it has never been setup */
1369371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi(PC pc) {
1374b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
1419566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
1422e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL));
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL));
1452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL));
1462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL));
1479566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1484b9ad928SBarry Smith   PetscFunctionReturn(0);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith 
1519371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject) {
1524b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
153248bfaf0SJed Brown   PetscInt    blocks, i;
154ace3abfcSBarry Smith   PetscBool   flg;
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith   PetscFunctionBegin;
157d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
1589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
1599566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
1619566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
162248bfaf0SJed Brown   if (jac->ksp) {
163248bfaf0SJed Brown     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
164248bfaf0SJed Brown      * unless we had already been called. */
16548a46eb9SPierre Jolivet     for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
166248bfaf0SJed Brown   }
167d0609cedSBarry Smith   PetscOptionsHeadEnd();
1684b9ad928SBarry Smith   PetscFunctionReturn(0);
1694b9ad928SBarry Smith }
1704b9ad928SBarry Smith 
1719804daf3SBarry Smith #include <petscdraw.h>
1729371c9d4SSatish Balay static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer) {
1734b9ad928SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
174cbe18068SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
17569a612a9SBarry Smith   PetscMPIInt           rank;
17669a612a9SBarry Smith   PetscInt              i;
177d9884438SBarry Smith   PetscBool             iascii, isstring, isdraw;
1784b9ad928SBarry Smith   PetscViewer           sviewer;
179020d6619SPierre Jolivet   PetscViewerFormat     format;
180020d6619SPierre Jolivet   const char           *prefix;
1814b9ad928SBarry Smith 
1824b9ad928SBarry Smith   PetscFunctionBegin;
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
18632077d6dSBarry Smith   if (iascii) {
18748a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
18863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of blocks = %" PetscInt_FMT "\n", jac->n));
1899566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
1909566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
191020d6619SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1929566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
1939566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
1949566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
195933995ecSHong Zhang       if (jac->ksp && !jac->psubcomm) {
1969566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
197dd400576SPatrick Sanan         if (rank == 0) {
1989566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
1999566063dSJacob Faibussowitsch           PetscCall(KSPView(jac->ksp[0], sviewer));
2009566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2014b9ad928SBarry Smith         }
2029566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2039566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2049566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
205e4de9e1dSBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
207e4fa1014SBarry Smith       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
2089566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
209e4fa1014SBarry Smith         if (!mpjac->psubcomm->color) {
2109566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2119566063dSJacob Faibussowitsch           PetscCall(KSPView(*(jac->ksp), sviewer));
2129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
213cbe18068SHong Zhang         }
2149566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2159566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2169566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2179530cbd7SBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2189566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2194b9ad928SBarry Smith       } else {
2209566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
221e4de9e1dSBarry Smith       }
222e4de9e1dSBarry Smith     } else {
2234814766eSBarry Smith       PetscInt n_global;
2241c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2259566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2269566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
2279371c9d4SSatish 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));
2289566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
230dd2fa690SBarry Smith       for (i = 0; i < jac->n_local; i++) {
23163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2329566063dSJacob Faibussowitsch         PetscCall(KSPView(jac->ksp[i], sviewer));
2339566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2344b9ad928SBarry Smith       }
2359566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
2379566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2394b9ad928SBarry Smith     }
2404b9ad928SBarry Smith   } else if (isstring) {
24163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2439566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
245d9884438SBarry Smith   } else if (isdraw) {
246d9884438SBarry Smith     PetscDraw draw;
247d9884438SBarry Smith     char      str[25];
248d9884438SBarry Smith     PetscReal x, y, bottom, h;
249d9884438SBarry Smith 
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2519566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
25263a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2539566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
254d9884438SBarry Smith     bottom = y - h;
2559566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
256d9884438SBarry Smith     /* warning the communicator on viewer is different then on ksp in parallel */
2579566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2589566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2594b9ad928SBarry Smith   }
2604b9ad928SBarry Smith   PetscFunctionReturn(0);
2614b9ad928SBarry Smith }
2624b9ad928SBarry Smith 
2639371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) {
264feb237baSPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2654b9ad928SBarry Smith 
2664b9ad928SBarry Smith   PetscFunctionBegin;
26728b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith   if (n_local) *n_local = jac->n_local;
2704b9ad928SBarry Smith   if (first_local) *first_local = jac->first_local;
271020d6619SPierre Jolivet   if (ksp) *ksp = jac->ksp;
2724b9ad928SBarry Smith   PetscFunctionReturn(0);
2734b9ad928SBarry Smith }
2744b9ad928SBarry Smith 
2759371c9d4SSatish Balay static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, PetscInt *lens) {
2764b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2774b9ad928SBarry Smith 
2784b9ad928SBarry Smith   PetscFunctionBegin;
2792472a847SBarry 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");
2804b9ad928SBarry Smith   jac->n = blocks;
2810a545947SLisandro Dalcin   if (!lens) jac->g_lens = NULL;
2822fa5cd67SKarl Rupp   else {
2839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2849566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)pc, blocks * sizeof(PetscInt)));
2859566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2864b9ad928SBarry Smith   }
2874b9ad928SBarry Smith   PetscFunctionReturn(0);
2884b9ad928SBarry Smith }
2894b9ad928SBarry Smith 
2909371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
2914b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2924b9ad928SBarry Smith 
2934b9ad928SBarry Smith   PetscFunctionBegin;
2944b9ad928SBarry Smith   *blocks = jac->n;
2954b9ad928SBarry Smith   if (lens) *lens = jac->g_lens;
2964b9ad928SBarry Smith   PetscFunctionReturn(0);
2974b9ad928SBarry Smith }
2984b9ad928SBarry Smith 
2999371c9d4SSatish Balay static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[]) {
3004b9ad928SBarry Smith   PC_BJacobi *jac;
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith   PetscFunctionBegin;
3034b9ad928SBarry Smith   jac = (PC_BJacobi *)pc->data;
3044b9ad928SBarry Smith 
3054b9ad928SBarry Smith   jac->n_local = blocks;
3060a545947SLisandro Dalcin   if (!lens) jac->l_lens = NULL;
3072fa5cd67SKarl Rupp   else {
3089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3099566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)pc, blocks * sizeof(PetscInt)));
3109566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3114b9ad928SBarry Smith   }
3124b9ad928SBarry Smith   PetscFunctionReturn(0);
3134b9ad928SBarry Smith }
3144b9ad928SBarry Smith 
3159371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
3164b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3174b9ad928SBarry Smith 
3184b9ad928SBarry Smith   PetscFunctionBegin;
3194b9ad928SBarry Smith   *blocks = jac->n_local;
3204b9ad928SBarry Smith   if (lens) *lens = jac->l_lens;
3214b9ad928SBarry Smith   PetscFunctionReturn(0);
3224b9ad928SBarry Smith }
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith /*@C
325*f1580f4eSBarry Smith    PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3264b9ad928SBarry Smith    this processor.
3274b9ad928SBarry Smith 
3286da92b7fSHong Zhang    Not Collective
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Input Parameter:
3314b9ad928SBarry Smith .  pc - the preconditioner context
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith    Output Parameters:
3340298fd71SBarry Smith +  n_local - the number of blocks on this processor, or NULL
3350298fd71SBarry Smith .  first_local - the global number of the first block on this processor, or NULL
3364b9ad928SBarry Smith -  ksp - the array of KSP contexts
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith    Notes:
339*f1580f4eSBarry Smith    After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3404b9ad928SBarry Smith 
3414b9ad928SBarry Smith    Currently for some matrix implementations only 1 block per processor
3424b9ad928SBarry Smith    is supported.
3434b9ad928SBarry Smith 
344*f1580f4eSBarry Smith    You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3454b9ad928SBarry Smith 
346*f1580f4eSBarry Smith    Fortran Usage:
347*f1580f4eSBarry Smith    You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
348*f1580f4eSBarry Smith 
349*f1580f4eSBarry Smith    You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
350*f1580f4eSBarry Smith    `KSP` array must be.
351196cc216SBarry Smith 
3524b9ad928SBarry Smith    Level: advanced
3534b9ad928SBarry Smith 
354*f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3554b9ad928SBarry Smith @*/
3569371c9d4SSatish Balay PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) {
3574b9ad928SBarry Smith   PetscFunctionBegin;
3580700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
359cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3604b9ad928SBarry Smith   PetscFunctionReturn(0);
3614b9ad928SBarry Smith }
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith /*@
3644b9ad928SBarry Smith    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3654b9ad928SBarry Smith    Jacobi preconditioner.
3664b9ad928SBarry Smith 
367*f1580f4eSBarry Smith    Collective on pc
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith    Input Parameters:
3704b9ad928SBarry Smith +  pc - the preconditioner context
3714b9ad928SBarry Smith .  blocks - the number of blocks
3724b9ad928SBarry Smith -  lens - [optional] integer array containing the size of each block
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith    Options Database Key:
3754b9ad928SBarry Smith .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3764b9ad928SBarry Smith 
377*f1580f4eSBarry Smith    Note:
3784b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
379*f1580f4eSBarry Smith    All processors sharing the `PC` must call this routine with the same data.
3804b9ad928SBarry Smith 
3814b9ad928SBarry Smith    Level: intermediate
3824b9ad928SBarry Smith 
383*f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3844b9ad928SBarry Smith @*/
3859371c9d4SSatish Balay PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) {
3864b9ad928SBarry Smith   PetscFunctionBegin;
3870700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
38808401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
389cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3904b9ad928SBarry Smith   PetscFunctionReturn(0);
3914b9ad928SBarry Smith }
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith /*@C
3944b9ad928SBarry Smith    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
395*f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`, preconditioner.
3964b9ad928SBarry Smith 
397ad4df100SBarry Smith    Not Collective
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith    Input Parameter:
4004b9ad928SBarry Smith .  pc - the preconditioner context
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith    Output parameters:
4034b9ad928SBarry Smith +  blocks - the number of blocks
4044b9ad928SBarry Smith -  lens - integer array containing the size of each block
4054b9ad928SBarry Smith 
4064b9ad928SBarry Smith    Level: intermediate
4074b9ad928SBarry Smith 
408*f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4094b9ad928SBarry Smith @*/
4109371c9d4SSatish Balay PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
4114b9ad928SBarry Smith   PetscFunctionBegin;
4120700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4134482741eSBarry Smith   PetscValidIntPointer(blocks, 2);
414cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4154b9ad928SBarry Smith   PetscFunctionReturn(0);
4164b9ad928SBarry Smith }
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith /*@
4194b9ad928SBarry Smith    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
420*f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`,  preconditioner.
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith    Not Collective
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith    Input Parameters:
4254b9ad928SBarry Smith +  pc - the preconditioner context
4264b9ad928SBarry Smith .  blocks - the number of blocks
4274b9ad928SBarry Smith -  lens - [optional] integer array containing size of each block
4284b9ad928SBarry Smith 
429342c94f9SMatthew G. Knepley    Options Database Key:
430342c94f9SMatthew G. Knepley .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
431342c94f9SMatthew G. Knepley 
4324b9ad928SBarry Smith    Note:
4334b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith    Level: intermediate
4364b9ad928SBarry Smith 
437*f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4384b9ad928SBarry Smith @*/
4399371c9d4SSatish Balay PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) {
4404b9ad928SBarry Smith   PetscFunctionBegin;
4410700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
44208401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
443cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4444b9ad928SBarry Smith   PetscFunctionReturn(0);
4454b9ad928SBarry Smith }
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith /*@C
4484b9ad928SBarry Smith    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
449*f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`, preconditioner.
4504b9ad928SBarry Smith 
4514b9ad928SBarry Smith    Not Collective
4524b9ad928SBarry Smith 
4534b9ad928SBarry Smith    Input Parameters:
4544b9ad928SBarry Smith +  pc - the preconditioner context
4554b9ad928SBarry Smith .  blocks - the number of blocks
4564b9ad928SBarry Smith -  lens - [optional] integer array containing size of each block
4574b9ad928SBarry Smith 
4584b9ad928SBarry Smith    Note:
4594b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith    Level: intermediate
4624b9ad928SBarry Smith 
463*f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4644b9ad928SBarry Smith @*/
4659371c9d4SSatish Balay PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) {
4664b9ad928SBarry Smith   PetscFunctionBegin;
4670700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4684482741eSBarry Smith   PetscValidIntPointer(blocks, 2);
469cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4704b9ad928SBarry Smith   PetscFunctionReturn(0);
4714b9ad928SBarry Smith }
4724b9ad928SBarry Smith 
4734b9ad928SBarry Smith /*MC
4744b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
475*f1580f4eSBarry Smith            its own `KSP` object.
4764b9ad928SBarry Smith 
4774b9ad928SBarry Smith    Options Database Keys:
478011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
479011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4804b9ad928SBarry Smith 
48195452b02SPatrick Sanan    Notes:
482*f1580f4eSBarry 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 
486*f1580f4eSBarry 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 
489*f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
490*f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
491*f1580f4eSBarry Smith          `KSPGetPC())`
4924b9ad928SBarry Smith 
493*f1580f4eSBarry 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 
4994b9ad928SBarry Smith    Level: beginner
5004b9ad928SBarry Smith 
501*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
502db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
503db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5044b9ad928SBarry Smith M*/
5054b9ad928SBarry Smith 
5069371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc) {
50769a612a9SBarry Smith   PetscMPIInt rank;
5084b9ad928SBarry Smith   PC_BJacobi *jac;
5094b9ad928SBarry Smith 
5104b9ad928SBarry Smith   PetscFunctionBegin;
5119566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &jac));
5129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5132fa5cd67SKarl Rupp 
5140a545947SLisandro Dalcin   pc->ops->apply           = NULL;
5157b6e2003SPierre Jolivet   pc->ops->matapply        = NULL;
5160a545947SLisandro Dalcin   pc->ops->applytranspose  = 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));
5374b9ad928SBarry Smith   PetscFunctionReturn(0);
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 */
5439371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc) {
5444b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5454b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5464b9ad928SBarry Smith 
5474b9ad928SBarry Smith   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
551e91c6855SBarry Smith   PetscFunctionReturn(0);
552e91c6855SBarry Smith }
553e91c6855SBarry Smith 
5549371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) {
555e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
556e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
557e91c6855SBarry Smith 
558e91c6855SBarry Smith   PetscFunctionBegin;
5599566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5609566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5619566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5629566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5632e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5644b9ad928SBarry Smith   PetscFunctionReturn(0);
5654b9ad928SBarry Smith }
5664b9ad928SBarry Smith 
5679371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) {
5684b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5692295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
570539c167fSBarry Smith   KSPConvergedReason reason;
5714b9ad928SBarry Smith 
5724b9ad928SBarry Smith   PetscFunctionBegin;
5739566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5749566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
575ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5764b9ad928SBarry Smith   PetscFunctionReturn(0);
5774b9ad928SBarry Smith }
5784b9ad928SBarry Smith 
5799371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
5804b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5814b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5826e111a19SKarl Rupp 
5834b9ad928SBarry Smith   PetscFunctionBegin;
5849566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
5859566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
586bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
587910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
588910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
5899566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
5909566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
5919566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
5929566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
5939566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
5944b9ad928SBarry Smith   PetscFunctionReturn(0);
5954b9ad928SBarry Smith }
5964b9ad928SBarry Smith 
5979371c9d4SSatish Balay static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y) {
5987b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
5997b6e2003SPierre Jolivet   Mat         sX, sY;
6007b6e2003SPierre Jolivet 
6017b6e2003SPierre Jolivet   PetscFunctionBegin;
6027b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6037b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6047b6e2003SPierre Jolivet      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(MatDenseGetLocalMatrix(X, &sX));
6079566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
6089566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
6097b6e2003SPierre Jolivet   PetscFunctionReturn(0);
6107b6e2003SPierre Jolivet }
6117b6e2003SPierre Jolivet 
6129371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
6134b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6144b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
615d9ca1df4SBarry Smith   PetscScalar            *y_array;
616d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6174b9ad928SBarry Smith   PC                      subpc;
6184b9ad928SBarry Smith 
6194b9ad928SBarry Smith   PetscFunctionBegin;
6204b9ad928SBarry Smith   /*
6214b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6224b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6234b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6244b9ad928SBarry Smith     for the sequential solve.
6254b9ad928SBarry Smith   */
6269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6289566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6299566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6304b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
6314b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6329566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6339566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6349566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6359566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6384b9ad928SBarry Smith   PetscFunctionReturn(0);
6394b9ad928SBarry Smith }
6404b9ad928SBarry Smith 
6419371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
6424b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6434b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
644d9ca1df4SBarry Smith   PetscScalar            *y_array;
645d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6464b9ad928SBarry Smith   PC                      subpc;
6474b9ad928SBarry Smith 
6484b9ad928SBarry Smith   PetscFunctionBegin;
6494b9ad928SBarry Smith   /*
6504b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6514b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6524b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6534b9ad928SBarry Smith     for the sequential solve.
6544b9ad928SBarry Smith   */
6559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6579566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6589566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6594b9ad928SBarry Smith 
6604b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
6614b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6624b9ad928SBarry Smith 
6639566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6649566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6654b9ad928SBarry Smith 
6669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6684b9ad928SBarry Smith   PetscFunctionReturn(0);
6694b9ad928SBarry Smith }
6704b9ad928SBarry Smith 
6719371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y) {
6724b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6734b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
674d9ca1df4SBarry Smith   PetscScalar            *y_array;
675d9ca1df4SBarry Smith   const PetscScalar      *x_array;
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));
6889566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
6899566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
6909566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6919566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6944b9ad928SBarry Smith   PetscFunctionReturn(0);
6954b9ad928SBarry Smith }
6964b9ad928SBarry Smith 
6979371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat) {
6984b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
69969a612a9SBarry Smith   PetscInt                m;
7004b9ad928SBarry Smith   KSP                     ksp;
7014b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
702de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7033f6dc190SJunchao Zhang   VecType                 vectype;
704ea41da7aSPierre Jolivet   const char             *prefix;
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith   PetscFunctionBegin;
7074b9ad928SBarry Smith   if (!pc->setupcalled) {
708c2efdce3SBarry Smith     if (!jac->ksp) {
709302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7102fa5cd67SKarl Rupp 
7119566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7129566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7139566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7149566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp));
7159566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7169566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7179566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7189566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7194b9ad928SBarry Smith 
720e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7214b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7224b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7237b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7244b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7254b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7264b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7274b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7284b9ad928SBarry Smith 
7299566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7304b9ad928SBarry Smith       jac->ksp[0] = ksp;
731c2efdce3SBarry Smith 
7329566063dSJacob Faibussowitsch       PetscCall(PetscNewLog(pc, &bjac));
7334b9ad928SBarry Smith       jac->data = (void *)bjac;
7344b9ad928SBarry Smith     } else {
735c2efdce3SBarry Smith       ksp  = jac->ksp[0];
736c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
737c2efdce3SBarry Smith     }
738c2efdce3SBarry Smith 
739c2efdce3SBarry Smith     /*
740c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
741c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
742c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
743c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
744c2efdce3SBarry Smith       KSPSolve() on the block.
745c2efdce3SBarry Smith     */
7469566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7479566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7489566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7499566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7509566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7519566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
7529566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->x));
7539566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->y));
754c2efdce3SBarry Smith   } else {
7554b9ad928SBarry Smith     ksp  = jac->ksp[0];
7564b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7574b9ad928SBarry Smith   }
7589566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
75949517cdeSBarry Smith   if (pc->useAmat) {
7609566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7619566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7624b9ad928SBarry Smith   } else {
7639566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7644b9ad928SBarry Smith   }
7659566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
76690f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
767248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
7689566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
769302a9ddcSMatthew Knepley   }
7704b9ad928SBarry Smith   PetscFunctionReturn(0);
7714b9ad928SBarry Smith }
7724b9ad928SBarry Smith 
7739371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc) {
7744b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
7754b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
77669a612a9SBarry Smith   PetscInt               i;
7774b9ad928SBarry Smith 
7784b9ad928SBarry Smith   PetscFunctionBegin;
779e91c6855SBarry Smith   if (bjac && bjac->pmat) {
7809566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
78148a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
782e91c6855SBarry Smith   }
7834b9ad928SBarry Smith 
7844b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
7859566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
786e91c6855SBarry Smith     if (bjac && bjac->x) {
7879566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
7889566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
7899566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
7904b9ad928SBarry Smith     }
791e91c6855SBarry Smith   }
7929566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
7939566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
794e91c6855SBarry Smith   PetscFunctionReturn(0);
795e91c6855SBarry Smith }
796e91c6855SBarry Smith 
7979371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) {
798e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
799c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
800e91c6855SBarry Smith   PetscInt               i;
801e91c6855SBarry Smith 
802e91c6855SBarry Smith   PetscFunctionBegin;
8039566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
804c2efdce3SBarry Smith   if (bjac) {
8059566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8069566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8079566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
808c2efdce3SBarry Smith   }
8099566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
81048a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8119566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8122e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8134b9ad928SBarry Smith   PetscFunctionReturn(0);
8144b9ad928SBarry Smith }
8154b9ad928SBarry Smith 
8169371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) {
8174b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
81869a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
819539c167fSBarry Smith   KSPConvergedReason reason;
8204b9ad928SBarry Smith 
8214b9ad928SBarry Smith   PetscFunctionBegin;
8224b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8239566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8249566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
825ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8264b9ad928SBarry Smith   }
8274b9ad928SBarry Smith   PetscFunctionReturn(0);
8284b9ad928SBarry Smith }
8294b9ad928SBarry Smith 
8304b9ad928SBarry Smith /*
8314b9ad928SBarry Smith       Preconditioner for block Jacobi
8324b9ad928SBarry Smith */
8339371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y) {
8344b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
83569a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8364b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
837d9ca1df4SBarry Smith   PetscScalar           *yin;
838d9ca1df4SBarry Smith   const PetscScalar     *xin;
83958ebbce7SBarry Smith 
8404b9ad928SBarry Smith   PetscFunctionBegin;
8419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8434b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8444b9ad928SBarry Smith     /*
8454b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8464b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8474b9ad928SBarry Smith        the global vector.
8484b9ad928SBarry Smith     */
8499566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8509566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8514b9ad928SBarry Smith 
8529566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8539566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8549566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8559566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
856d11f3a42SBarry Smith 
8579566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8589566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8594b9ad928SBarry Smith   }
8609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8624b9ad928SBarry Smith   PetscFunctionReturn(0);
8634b9ad928SBarry Smith }
8644b9ad928SBarry Smith 
8654b9ad928SBarry Smith /*
8664b9ad928SBarry Smith       Preconditioner for block Jacobi
8674b9ad928SBarry Smith */
8689371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y) {
8694b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
87069a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8714b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
872d9ca1df4SBarry Smith   PetscScalar           *yin;
873d9ca1df4SBarry Smith   const PetscScalar     *xin;
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith   PetscFunctionBegin;
8769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8784b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8794b9ad928SBarry Smith     /*
8804b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8814b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8824b9ad928SBarry Smith        the global vector.
8834b9ad928SBarry Smith     */
8849566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8859566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8864b9ad928SBarry Smith 
8879566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8889566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
8899566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8909566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
891a7ff49e8SJed Brown 
8929566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8939566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8944b9ad928SBarry Smith   }
8959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8974b9ad928SBarry Smith   PetscFunctionReturn(0);
8984b9ad928SBarry Smith }
8994b9ad928SBarry Smith 
9009371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat) {
9014b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
90269a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
903ea41da7aSPierre Jolivet   const char            *prefix;
9044b9ad928SBarry Smith   KSP                    ksp;
9054b9ad928SBarry Smith   Vec                    x, y;
9064b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
9074b9ad928SBarry Smith   PC                     subpc;
9084b9ad928SBarry Smith   IS                     is;
909434a97beSBrad Aagaard   MatReuse               scall;
9103f6dc190SJunchao Zhang   VecType                vectype;
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith   PetscFunctionBegin;
9139566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith   n_local = jac->n_local;
9164b9ad928SBarry Smith 
91749517cdeSBarry Smith   if (pc->useAmat) {
918ace3abfcSBarry Smith     PetscBool same;
9199566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
92028b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
9214b9ad928SBarry Smith   }
9224b9ad928SBarry Smith 
9234b9ad928SBarry Smith   if (!pc->setupcalled) {
9244b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
925c2efdce3SBarry Smith 
926c2efdce3SBarry Smith     if (!jac->ksp) {
927e91c6855SBarry Smith       pc->ops->reset          = PCReset_BJacobi_Multiblock;
9284b9ad928SBarry Smith       pc->ops->destroy        = PCDestroy_BJacobi_Multiblock;
9294b9ad928SBarry Smith       pc->ops->apply          = PCApply_BJacobi_Multiblock;
9307b6e2003SPierre Jolivet       pc->ops->matapply       = NULL;
9314b9ad928SBarry Smith       pc->ops->applytranspose = PCApplyTranspose_BJacobi_Multiblock;
9324b9ad928SBarry Smith       pc->ops->setuponblocks  = PCSetUpOnBlocks_BJacobi_Multiblock;
9334b9ad928SBarry Smith 
9349566063dSJacob Faibussowitsch       PetscCall(PetscNewLog(pc, &bjac));
9359566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
9369566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(KSP))));
9379566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
9389566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
9399566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(PetscScalar))));
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith       jac->data = (void *)bjac;
9429566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
9439566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(n_local * sizeof(IS))));
9444b9ad928SBarry Smith 
9454b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
9469566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
9479566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
9489566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
9499566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp));
9509566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
9519566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
9529566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
9539566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
9549566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
9552fa5cd67SKarl Rupp 
956c2efdce3SBarry Smith         jac->ksp[i] = ksp;
957c2efdce3SBarry Smith       }
958c2efdce3SBarry Smith     } else {
959c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
960c2efdce3SBarry Smith     }
9614b9ad928SBarry Smith 
962c2efdce3SBarry Smith     start = 0;
9639566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
964c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
9654b9ad928SBarry Smith       m = jac->l_lens[i];
9664b9ad928SBarry Smith       /*
9674b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
9684b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
9694b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
9704b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
9714b9ad928SBarry Smith       KSPSolve() on the block.
9724b9ad928SBarry Smith 
9734b9ad928SBarry Smith       */
9749566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
9759566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
9769566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
9779566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
9789566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)x));
9799566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)y));
9802fa5cd67SKarl Rupp 
9814b9ad928SBarry Smith       bjac->x[i]      = x;
9824b9ad928SBarry Smith       bjac->y[i]      = y;
9834b9ad928SBarry Smith       bjac->starts[i] = start;
9844b9ad928SBarry Smith 
9859566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
9864b9ad928SBarry Smith       bjac->is[i] = is;
9879566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)is));
9884b9ad928SBarry Smith 
9894b9ad928SBarry Smith       start += m;
9904b9ad928SBarry Smith     }
9914b9ad928SBarry Smith   } else {
9924b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
9934b9ad928SBarry Smith     /*
9944b9ad928SBarry Smith        Destroy the blocks from the previous iteration
9954b9ad928SBarry Smith     */
9964b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
9979566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
99848a46eb9SPierre Jolivet       if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
9994b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
10002fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
10014b9ad928SBarry Smith   }
10024b9ad928SBarry Smith 
10039566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
100448a46eb9SPierre Jolivet   if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
10054b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
10064b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
10079566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
10109566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->pmat[i]));
10119566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
101249517cdeSBarry Smith     if (pc->useAmat) {
10139566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)bjac->mat[i]));
10149566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
10159566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
10164b9ad928SBarry Smith     } else {
10179566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
10184b9ad928SBarry Smith     }
10199566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
102048a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
102190f1c854SHong Zhang   }
10224b9ad928SBarry Smith   PetscFunctionReturn(0);
10234b9ad928SBarry Smith }
10245a7e1895SHong Zhang 
10255a7e1895SHong Zhang /*
1026fd0b8932SBarry Smith       These are for a single block with multiple processes
10275a7e1895SHong Zhang */
10289371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc) {
1029fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1030fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1031fd0b8932SBarry Smith   KSPConvergedReason reason;
1032fd0b8932SBarry Smith 
1033fd0b8932SBarry Smith   PetscFunctionBegin;
10349566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
10359566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1036ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1037fd0b8932SBarry Smith   PetscFunctionReturn(0);
1038fd0b8932SBarry Smith }
1039fd0b8932SBarry Smith 
10409371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc) {
10415a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
10425a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
10435a7e1895SHong Zhang 
10445a7e1895SHong Zhang   PetscFunctionBegin;
10459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
10469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
10479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
10489566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
1049e91c6855SBarry Smith   PetscFunctionReturn(0);
1050e91c6855SBarry Smith }
1051e91c6855SBarry Smith 
10529371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc) {
1053e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1054e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1055e91c6855SBarry Smith 
1056e91c6855SBarry Smith   PetscFunctionBegin;
10579566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
10589566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
10599566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
10609566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
10615a7e1895SHong Zhang 
10629566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
10632e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
10645a7e1895SHong Zhang   PetscFunctionReturn(0);
10655a7e1895SHong Zhang }
10665a7e1895SHong Zhang 
10679371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y) {
10685a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
10695a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1070d9ca1df4SBarry Smith   PetscScalar          *yarray;
1071d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1072539c167fSBarry Smith   KSPConvergedReason    reason;
10735a7e1895SHong Zhang 
10745a7e1895SHong Zhang   PetscFunctionBegin;
10755a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
10769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
10779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
10789566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
10799566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
10805a7e1895SHong Zhang 
10815a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
10829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
10839566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
10849566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
10859566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
10869566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1087ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1088250cf88bSHong Zhang 
10899566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
10909566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
10919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
10929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
10935a7e1895SHong Zhang   PetscFunctionReturn(0);
10945a7e1895SHong Zhang }
10955a7e1895SHong Zhang 
10969371c9d4SSatish Balay static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y) {
10977b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
10987b6e2003SPierre Jolivet   KSPConvergedReason reason;
10997b6e2003SPierre Jolivet   Mat                sX, sY;
11007b6e2003SPierre Jolivet   const PetscScalar *x;
11017b6e2003SPierre Jolivet   PetscScalar       *y;
11027b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
11037b6e2003SPierre Jolivet 
11047b6e2003SPierre Jolivet   PetscFunctionBegin;
11057b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
11069566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
11079566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
11089566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
11099566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
11109566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
11119566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
11129566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
11139566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
11149566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
11159566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
11169566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
11179566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
11189566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
11199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
11209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
11219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
11229566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
11239566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
11249566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1125ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11267b6e2003SPierre Jolivet   PetscFunctionReturn(0);
11277b6e2003SPierre Jolivet }
11287b6e2003SPierre Jolivet 
11299371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc) {
11305a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11315a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11325a7e1895SHong Zhang   PetscInt              m, n;
1133ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
11345a7e1895SHong Zhang   const char           *prefix;
1135de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
11363f6dc190SJunchao Zhang   VecType               vectype;
11371f62890fSKarl Rupp 
11385a7e1895SHong Zhang   PetscFunctionBegin;
11399566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
114008401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
11415a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
11425a7e1895SHong Zhang   if (!pc->setupcalled) {
1143de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
11449566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(pc, &mpjac));
11455a7e1895SHong Zhang     jac->data = (void *)mpjac;
11465a7e1895SHong Zhang 
11475a7e1895SHong Zhang     /* initialize datastructure mpjac */
11485a7e1895SHong Zhang     if (!jac->psubcomm) {
11495a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
11509566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
11519566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
11529566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
11539566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(PetscSubcomm)));
11545a7e1895SHong Zhang     }
11555a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1156306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
11575a7e1895SHong Zhang 
11585a7e1895SHong Zhang     /* Get matrix blocks of pmat */
11599566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
11605a7e1895SHong Zhang 
11615a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
11629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
11639566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
11649566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
11659566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
11669566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)jac->ksp[0]));
11679566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
11689566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
11695a7e1895SHong Zhang 
11709566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
11719566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
11729566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
11739566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
11749566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
11755a7e1895SHong Zhang 
11765a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
11779566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
11789566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
11799566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
11809566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
11819566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
11829566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
11839566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mpjac->xsub));
11849566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)mpjac->ysub));
11855a7e1895SHong Zhang 
1186fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1187e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
11885a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
11895a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
11907b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1191fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1192306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
11939e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
11945a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
11959566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
11969566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1197fc08c53fSHong Zhang     } else {
11989566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
11995a7e1895SHong Zhang     }
12009566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12015a7e1895SHong Zhang   }
12025a7e1895SHong Zhang 
120348a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
12045a7e1895SHong Zhang   PetscFunctionReturn(0);
12055a7e1895SHong Zhang }
1206