xref: /petsc/src/ksp/pc/impls/bjacobi/bjacobi.c (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
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 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi(PC pc)
13d71ae5a4SJacob Faibussowitsch {
144b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
154b9ad928SBarry Smith   Mat         mat = pc->mat, pmat = pc->pmat;
16976e8c5aSLisandro Dalcin   PetscBool   hasop;
1769a612a9SBarry Smith   PetscInt    N, M, start, i, sum, end;
1869a612a9SBarry Smith   PetscInt    bs, i_start = -1, i_end = -1;
1969a612a9SBarry Smith   PetscMPIInt rank, size;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
249566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
259566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(pc->pmat, &bs));
264b9ad928SBarry Smith 
275a7e1895SHong Zhang   if (jac->n > 0 && jac->n < size) {
289566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiproc(pc));
293ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
305a7e1895SHong Zhang   }
315a7e1895SHong Zhang 
32f1580f4eSBarry Smith   /*    Determines the number of blocks assigned to each processor */
334b9ad928SBarry Smith   /*   local block count  given */
344b9ad928SBarry Smith   if (jac->n_local > 0 && jac->n < 0) {
351c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
364b9ad928SBarry Smith     if (jac->l_lens) { /* check that user set these correctly */
374b9ad928SBarry Smith       sum = 0;
384b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
3908401ef6SPierre 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");
404b9ad928SBarry Smith         sum += jac->l_lens[i];
414b9ad928SBarry Smith       }
4208401ef6SPierre Jolivet       PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly");
434b9ad928SBarry Smith     } else {
449566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
452fa5cd67SKarl 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));
464b9ad928SBarry Smith     }
474b9ad928SBarry Smith   } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
484b9ad928SBarry Smith     /* global blocks given: determine which ones are local */
494b9ad928SBarry Smith     if (jac->g_lens) {
504b9ad928SBarry Smith       /* check if the g_lens is has valid entries */
514b9ad928SBarry Smith       for (i = 0; i < jac->n; i++) {
527827d75bSBarry Smith         PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed");
5308401ef6SPierre 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");
544b9ad928SBarry Smith       }
554b9ad928SBarry Smith       if (size == 1) {
564b9ad928SBarry Smith         jac->n_local = jac->n;
579566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
589566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local));
594b9ad928SBarry Smith         /* check that user set these correctly */
604b9ad928SBarry Smith         sum = 0;
614b9ad928SBarry Smith         for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i];
6208401ef6SPierre Jolivet         PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly");
634b9ad928SBarry Smith       } else {
649566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end));
65*aaa8cc7dSPierre Jolivet         /* loop over blocks determining first one owned by me */
664b9ad928SBarry Smith         sum = 0;
674b9ad928SBarry Smith         for (i = 0; i < jac->n + 1; i++) {
689371c9d4SSatish Balay           if (sum == start) {
699371c9d4SSatish Balay             i_start = i;
709371c9d4SSatish Balay             goto start_1;
719371c9d4SSatish Balay           }
724b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
734b9ad928SBarry Smith         }
74e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
754b9ad928SBarry Smith       start_1:
764b9ad928SBarry Smith         for (i = i_start; i < jac->n + 1; i++) {
779371c9d4SSatish Balay           if (sum == end) {
789371c9d4SSatish Balay             i_end = i;
799371c9d4SSatish Balay             goto end_1;
809371c9d4SSatish Balay           }
814b9ad928SBarry Smith           if (i < jac->n) sum += jac->g_lens[i];
824b9ad928SBarry Smith         }
83e7e72b3dSBarry Smith         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
844b9ad928SBarry Smith       end_1:
854b9ad928SBarry Smith         jac->n_local = i_end - i_start;
869566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
879566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local));
884b9ad928SBarry Smith       }
894b9ad928SBarry Smith     } else { /* no global blocks given, determine then using default layout */
904b9ad928SBarry Smith       jac->n_local = jac->n / size + ((jac->n % size) > rank);
919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
924b9ad928SBarry Smith       for (i = 0; i < jac->n_local; i++) {
934b9ad928SBarry Smith         jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs;
947827d75bSBarry Smith         PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given");
954b9ad928SBarry Smith       }
964b9ad928SBarry Smith     }
974b9ad928SBarry Smith   } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
984b9ad928SBarry Smith     jac->n       = size;
994b9ad928SBarry Smith     jac->n_local = 1;
1009566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->l_lens));
1014b9ad928SBarry Smith     jac->l_lens[0] = M;
1027271b979SBarry Smith   } else { /* jac->n > 0 && jac->n_local > 0 */
1037271b979SBarry Smith     if (!jac->l_lens) {
1049566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
1052fa5cd67SKarl 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));
1067271b979SBarry Smith     }
1074b9ad928SBarry Smith   }
10808401ef6SPierre Jolivet   PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors");
1094b9ad928SBarry Smith 
110f1580f4eSBarry Smith   /*    Determines mat and pmat */
1119566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
112976e8c5aSLisandro Dalcin   if (!hasop && size == 1) {
1134b9ad928SBarry Smith     mat  = pc->mat;
1144b9ad928SBarry Smith     pmat = pc->pmat;
1154b9ad928SBarry Smith   } else {
11649517cdeSBarry Smith     if (pc->useAmat) {
11749517cdeSBarry Smith       /* use block from Amat matrix, not Pmat for local MatMult() */
1189566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->mat, &mat));
1194b9ad928SBarry Smith     }
12049517cdeSBarry Smith     if (pc->pmat != pc->mat || !pc->useAmat) {
1219566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat));
1222fa5cd67SKarl Rupp     } else pmat = mat;
1234b9ad928SBarry Smith   }
1244b9ad928SBarry Smith 
125f1580f4eSBarry Smith   /*
1264b9ad928SBarry Smith      Setup code depends on the number of blocks
1274b9ad928SBarry Smith   */
128cc1d6079SHong Zhang   if (jac->n_local == 1) {
1299566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat));
1304b9ad928SBarry Smith   } else {
1319566063dSJacob Faibussowitsch     PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat));
1324b9ad928SBarry Smith   }
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1344b9ad928SBarry Smith }
1354b9ad928SBarry Smith 
1364b9ad928SBarry Smith /* Default destroy, if it has never been setup */
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi(PC pc)
138d71ae5a4SJacob Faibussowitsch {
1394b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1404b9ad928SBarry Smith 
1414b9ad928SBarry Smith   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
1439566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL));
1452e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL));
1462e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL));
1472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL));
1482e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL));
1499566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1514b9ad928SBarry Smith }
1524b9ad928SBarry Smith 
153d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject)
154d71ae5a4SJacob Faibussowitsch {
1554b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
156248bfaf0SJed Brown   PetscInt    blocks, i;
157ace3abfcSBarry Smith   PetscBool   flg;
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith   PetscFunctionBegin;
160d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
1619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
1629566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
1639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
1649566063dSJacob Faibussowitsch   if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
165248bfaf0SJed Brown   if (jac->ksp) {
166248bfaf0SJed Brown     /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
167248bfaf0SJed Brown      * unless we had already been called. */
16848a46eb9SPierre Jolivet     for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
169248bfaf0SJed Brown   }
170d0609cedSBarry Smith   PetscOptionsHeadEnd();
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1724b9ad928SBarry Smith }
1734b9ad928SBarry Smith 
1749804daf3SBarry Smith #include <petscdraw.h>
175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer)
176d71ae5a4SJacob Faibussowitsch {
1774b9ad928SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
178cbe18068SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
17969a612a9SBarry Smith   PetscMPIInt           rank;
18069a612a9SBarry Smith   PetscInt              i;
181d9884438SBarry Smith   PetscBool             iascii, isstring, isdraw;
1824b9ad928SBarry Smith   PetscViewer           sviewer;
183020d6619SPierre Jolivet   PetscViewerFormat     format;
184020d6619SPierre Jolivet   const char           *prefix;
1854b9ad928SBarry Smith 
1864b9ad928SBarry Smith   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1899566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
19032077d6dSBarry Smith   if (iascii) {
19148a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, "  using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
19263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  number of blocks = %" PetscInt_FMT "\n", jac->n));
1939566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
1949566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
195020d6619SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1969566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
1979566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
1989566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
199933995ecSHong Zhang       if (jac->ksp && !jac->psubcomm) {
2009566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
201dd400576SPatrick Sanan         if (rank == 0) {
2029566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2039566063dSJacob Faibussowitsch           PetscCall(KSPView(jac->ksp[0], sviewer));
2049566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2054b9ad928SBarry Smith         }
2069566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2079566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2089566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
209e4de9e1dSBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
211e4fa1014SBarry Smith       } else if (mpjac && jac->ksp && mpjac->psubcomm) {
2129566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
213e4fa1014SBarry Smith         if (!mpjac->psubcomm->color) {
2149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2159566063dSJacob Faibussowitsch           PetscCall(KSPView(*(jac->ksp), sviewer));
2169566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
217cbe18068SHong Zhang         }
2189566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(sviewer));
2199566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2209566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
2219530cbd7SBarry Smith         /*  extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2234b9ad928SBarry Smith       } else {
2249566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
225e4de9e1dSBarry Smith       }
226e4de9e1dSBarry Smith     } else {
2274814766eSBarry Smith       PetscInt n_global;
2281c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
2319371c9d4SSatish 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));
2329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
2339566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
234dd2fa690SBarry Smith       for (i = 0; i < jac->n_local; i++) {
23563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2369566063dSJacob Faibussowitsch         PetscCall(KSPView(jac->ksp[i], sviewer));
2379566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
2384b9ad928SBarry Smith       }
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2434b9ad928SBarry Smith     }
2444b9ad928SBarry Smith   } else if (isstring) {
24563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2479566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
249d9884438SBarry Smith   } else if (isdraw) {
250d9884438SBarry Smith     PetscDraw draw;
251d9884438SBarry Smith     char      str[25];
252d9884438SBarry Smith     PetscReal x, y, bottom, h;
253d9884438SBarry Smith 
2549566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2559566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
25663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2579566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
258d9884438SBarry Smith     bottom = y - h;
2599566063dSJacob Faibussowitsch     PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
260d9884438SBarry Smith     /* warning the communicator on viewer is different then on ksp in parallel */
2619566063dSJacob Faibussowitsch     if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2629566063dSJacob Faibussowitsch     PetscCall(PetscDrawPopCurrentPoint(draw));
2634b9ad928SBarry Smith   }
2643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2654b9ad928SBarry Smith }
2664b9ad928SBarry Smith 
267d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
268d71ae5a4SJacob Faibussowitsch {
269feb237baSPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2704b9ad928SBarry Smith 
2714b9ad928SBarry Smith   PetscFunctionBegin;
27228b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2734b9ad928SBarry Smith 
2744b9ad928SBarry Smith   if (n_local) *n_local = jac->n_local;
2754b9ad928SBarry Smith   if (first_local) *first_local = jac->first_local;
276020d6619SPierre Jolivet   if (ksp) *ksp = jac->ksp;
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2784b9ad928SBarry Smith }
2794b9ad928SBarry Smith 
280d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, PetscInt *lens)
281d71ae5a4SJacob Faibussowitsch {
2824b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2834b9ad928SBarry Smith 
2844b9ad928SBarry Smith   PetscFunctionBegin;
2852472a847SBarry 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");
2864b9ad928SBarry Smith   jac->n = blocks;
2870a545947SLisandro Dalcin   if (!lens) jac->g_lens = NULL;
2882fa5cd67SKarl Rupp   else {
2899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2909566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2914b9ad928SBarry Smith   }
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2934b9ad928SBarry Smith }
2944b9ad928SBarry Smith 
295d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
296d71ae5a4SJacob Faibussowitsch {
2974b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2984b9ad928SBarry Smith 
2994b9ad928SBarry Smith   PetscFunctionBegin;
3004b9ad928SBarry Smith   *blocks = jac->n;
3014b9ad928SBarry Smith   if (lens) *lens = jac->g_lens;
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3034b9ad928SBarry Smith }
3044b9ad928SBarry Smith 
305d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[])
306d71ae5a4SJacob Faibussowitsch {
3074b9ad928SBarry Smith   PC_BJacobi *jac;
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith   PetscFunctionBegin;
3104b9ad928SBarry Smith   jac = (PC_BJacobi *)pc->data;
3114b9ad928SBarry Smith 
3124b9ad928SBarry Smith   jac->n_local = blocks;
3130a545947SLisandro Dalcin   if (!lens) jac->l_lens = NULL;
3142fa5cd67SKarl Rupp   else {
3159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3169566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3174b9ad928SBarry Smith   }
3183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3194b9ad928SBarry Smith }
3204b9ad928SBarry Smith 
321d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
322d71ae5a4SJacob Faibussowitsch {
3234b9ad928SBarry Smith   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith   PetscFunctionBegin;
3264b9ad928SBarry Smith   *blocks = jac->n_local;
3274b9ad928SBarry Smith   if (lens) *lens = jac->l_lens;
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3294b9ad928SBarry Smith }
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith /*@C
332f1580f4eSBarry Smith    PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3334b9ad928SBarry Smith    this processor.
3344b9ad928SBarry Smith 
3356da92b7fSHong Zhang    Not Collective
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith    Input Parameter:
3384b9ad928SBarry Smith .  pc - the preconditioner context
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith    Output Parameters:
3410298fd71SBarry Smith +  n_local - the number of blocks on this processor, or NULL
3420298fd71SBarry Smith .  first_local - the global number of the first block on this processor, or NULL
3434b9ad928SBarry Smith -  ksp - the array of KSP contexts
3444b9ad928SBarry Smith 
3454b9ad928SBarry Smith    Notes:
346f1580f4eSBarry Smith    After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed.
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith    Currently for some matrix implementations only 1 block per processor
3494b9ad928SBarry Smith    is supported.
3504b9ad928SBarry Smith 
351f1580f4eSBarry Smith    You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3524b9ad928SBarry Smith 
353f1580f4eSBarry Smith    Fortran Usage:
354f1580f4eSBarry Smith    You must pass in a `KSP` array that is large enough to contain all the local `KSP`s.
355f1580f4eSBarry Smith 
356f1580f4eSBarry Smith    You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the
357f1580f4eSBarry Smith    `KSP` array must be.
358196cc216SBarry Smith 
3594b9ad928SBarry Smith    Level: advanced
3604b9ad928SBarry Smith 
361f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3624b9ad928SBarry Smith @*/
363d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
364d71ae5a4SJacob Faibussowitsch {
3654b9ad928SBarry Smith   PetscFunctionBegin;
3660700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
367cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3694b9ad928SBarry Smith }
3704b9ad928SBarry Smith 
3714b9ad928SBarry Smith /*@
3724b9ad928SBarry Smith    PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3734b9ad928SBarry Smith    Jacobi preconditioner.
3744b9ad928SBarry Smith 
375c3339decSBarry Smith    Collective
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith    Input Parameters:
3784b9ad928SBarry Smith +  pc - the preconditioner context
3794b9ad928SBarry Smith .  blocks - the number of blocks
3804b9ad928SBarry Smith -  lens - [optional] integer array containing the size of each block
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith    Options Database Key:
3834b9ad928SBarry Smith .  -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3844b9ad928SBarry Smith 
385f1580f4eSBarry Smith    Note:
3864b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
387f1580f4eSBarry Smith    All processors sharing the `PC` must call this routine with the same data.
3884b9ad928SBarry Smith 
3894b9ad928SBarry Smith    Level: intermediate
3904b9ad928SBarry Smith 
391f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3924b9ad928SBarry Smith @*/
393d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
394d71ae5a4SJacob Faibussowitsch {
3954b9ad928SBarry Smith   PetscFunctionBegin;
3960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
39708401ef6SPierre Jolivet   PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
398cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4004b9ad928SBarry Smith }
4014b9ad928SBarry Smith 
4024b9ad928SBarry Smith /*@C
4034b9ad928SBarry Smith    PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
404f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`, preconditioner.
4054b9ad928SBarry Smith 
406ad4df100SBarry Smith    Not Collective
4074b9ad928SBarry Smith 
4084b9ad928SBarry Smith    Input Parameter:
4094b9ad928SBarry Smith .  pc - the preconditioner context
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith    Output parameters:
4124b9ad928SBarry Smith +  blocks - the number of blocks
4134b9ad928SBarry Smith -  lens - integer array containing the size of each block
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith    Level: intermediate
4164b9ad928SBarry Smith 
417f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4184b9ad928SBarry Smith @*/
419d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
420d71ae5a4SJacob Faibussowitsch {
4214b9ad928SBarry Smith   PetscFunctionBegin;
4220700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4234482741eSBarry Smith   PetscValidIntPointer(blocks, 2);
424cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4264b9ad928SBarry Smith }
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith /*@
4294b9ad928SBarry Smith    PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
430f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`,  preconditioner.
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith    Not Collective
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith    Input Parameters:
4354b9ad928SBarry Smith +  pc - the preconditioner context
4364b9ad928SBarry Smith .  blocks - the number of blocks
4374b9ad928SBarry Smith -  lens - [optional] integer array containing size of each block
4384b9ad928SBarry Smith 
439342c94f9SMatthew G. Knepley    Options Database Key:
440342c94f9SMatthew G. Knepley .  -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
441342c94f9SMatthew G. Knepley 
4424b9ad928SBarry Smith    Note:
4434b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
4444b9ad928SBarry Smith 
4454b9ad928SBarry Smith    Level: intermediate
4464b9ad928SBarry Smith 
447f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4484b9ad928SBarry Smith @*/
449d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
450d71ae5a4SJacob Faibussowitsch {
4514b9ad928SBarry Smith   PetscFunctionBegin;
4520700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
45308401ef6SPierre Jolivet   PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
454cac4c232SBarry Smith   PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4564b9ad928SBarry Smith }
4574b9ad928SBarry Smith 
4584b9ad928SBarry Smith /*@C
4594b9ad928SBarry Smith    PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
460f1580f4eSBarry Smith    Jacobi, `PCBJACOBI`, preconditioner.
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith    Not Collective
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith    Input Parameters:
4654b9ad928SBarry Smith +  pc - the preconditioner context
4664b9ad928SBarry Smith .  blocks - the number of blocks
4674b9ad928SBarry Smith -  lens - [optional] integer array containing size of each block
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith    Note:
4704b9ad928SBarry Smith    Currently only a limited number of blocking configurations are supported.
4714b9ad928SBarry Smith 
4724b9ad928SBarry Smith    Level: intermediate
4734b9ad928SBarry Smith 
474f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4754b9ad928SBarry Smith @*/
476d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
477d71ae5a4SJacob Faibussowitsch {
4784b9ad928SBarry Smith   PetscFunctionBegin;
4790700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4804482741eSBarry Smith   PetscValidIntPointer(blocks, 2);
481cac4c232SBarry Smith   PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4834b9ad928SBarry Smith }
4844b9ad928SBarry Smith 
4854b9ad928SBarry Smith /*MC
4864b9ad928SBarry Smith    PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
487f1580f4eSBarry Smith            its own `KSP` object.
4884b9ad928SBarry Smith 
4894b9ad928SBarry Smith    Options Database Keys:
490011ea8aeSBarry Smith +  -pc_use_amat - use Amat to apply block of operator in inner Krylov method
491011ea8aeSBarry Smith -  -pc_bjacobi_blocks <n> - use n total blocks
4924b9ad928SBarry Smith 
49395452b02SPatrick Sanan    Notes:
494f1580f4eSBarry Smith     See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
495468ae2e8SBarry Smith 
49695452b02SPatrick 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.
4974b9ad928SBarry Smith 
498f1580f4eSBarry Smith      To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
499d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
5004b9ad928SBarry Smith 
501f1580f4eSBarry Smith      To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
502f1580f4eSBarry Smith          and set the options directly on the resulting `KSP` object (you can access its `PC`
5030462cc06SPierre Jolivet          `KSPGetPC()`)
5044b9ad928SBarry Smith 
505f1580f4eSBarry Smith      For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
5062210c163SDominic Meiser          performance.  Different block partitioning may lead to additional data transfers
5072210c163SDominic Meiser          between host and GPU that lead to degraded performance.
5082210c163SDominic Meiser 
509011ea8aeSBarry Smith      When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
510011ea8aeSBarry Smith 
5114b9ad928SBarry Smith    Level: beginner
5124b9ad928SBarry Smith 
513f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
514db781477SPatrick Sanan           `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
515db781477SPatrick Sanan           `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5164b9ad928SBarry Smith M*/
5174b9ad928SBarry Smith 
518d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
519d71ae5a4SJacob Faibussowitsch {
52069a612a9SBarry Smith   PetscMPIInt rank;
5214b9ad928SBarry Smith   PC_BJacobi *jac;
5224b9ad928SBarry Smith 
5234b9ad928SBarry Smith   PetscFunctionBegin;
5244dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
5259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5262fa5cd67SKarl Rupp 
5270a545947SLisandro Dalcin   pc->ops->apply           = NULL;
5287b6e2003SPierre Jolivet   pc->ops->matapply        = NULL;
5290a545947SLisandro Dalcin   pc->ops->applytranspose  = NULL;
5304b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_BJacobi;
5314b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_BJacobi;
5324b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_BJacobi;
5334b9ad928SBarry Smith   pc->ops->view            = PCView_BJacobi;
5340a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
5354b9ad928SBarry Smith 
5364b9ad928SBarry Smith   pc->data         = (void *)jac;
5374b9ad928SBarry Smith   jac->n           = -1;
5384b9ad928SBarry Smith   jac->n_local     = -1;
5394b9ad928SBarry Smith   jac->first_local = rank;
5400a545947SLisandro Dalcin   jac->ksp         = NULL;
5410a545947SLisandro Dalcin   jac->g_lens      = NULL;
5420a545947SLisandro Dalcin   jac->l_lens      = NULL;
5430a545947SLisandro Dalcin   jac->psubcomm    = NULL;
5444b9ad928SBarry Smith 
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5514b9ad928SBarry Smith }
5524b9ad928SBarry Smith 
5534b9ad928SBarry Smith /*
5544b9ad928SBarry Smith         These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5554b9ad928SBarry Smith */
556d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
557d71ae5a4SJacob Faibussowitsch {
5584b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5594b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5604b9ad928SBarry Smith 
5614b9ad928SBarry Smith   PetscFunctionBegin;
5629566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp[0]));
5639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->x));
5649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bjac->y));
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
566e91c6855SBarry Smith }
567e91c6855SBarry Smith 
568d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
569d71ae5a4SJacob Faibussowitsch {
570e91c6855SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
571e91c6855SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
572e91c6855SBarry Smith 
573e91c6855SBarry Smith   PetscFunctionBegin;
5749566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Singleblock(pc));
5759566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
5769566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
5779566063dSJacob Faibussowitsch   PetscCall(PetscFree(bjac));
5782e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5804b9ad928SBarry Smith }
5814b9ad928SBarry Smith 
582d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
583d71ae5a4SJacob Faibussowitsch {
5844b9ad928SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
5852295b7c8SHong Zhang   KSP                subksp = jac->ksp[0];
586539c167fSBarry Smith   KSPConvergedReason reason;
5874b9ad928SBarry Smith 
5884b9ad928SBarry Smith   PetscFunctionBegin;
5899566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
5909566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
591ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5934b9ad928SBarry Smith }
5944b9ad928SBarry Smith 
595d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
596d71ae5a4SJacob Faibussowitsch {
5974b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
5984b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5996e111a19SKarl Rupp 
6004b9ad928SBarry Smith   PetscFunctionBegin;
6019566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVectorRead(x, bjac->x));
6029566063dSJacob Faibussowitsch   PetscCall(VecGetLocalVector(y, bjac->y));
603bba28a21SBarry Smith   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
604910cf402Sprj-      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
605910cf402Sprj-      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6069566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6079566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
6089566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
6099566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
6109566063dSJacob Faibussowitsch   PetscCall(VecRestoreLocalVector(y, bjac->y));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6124b9ad928SBarry Smith }
6134b9ad928SBarry Smith 
614d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
615d71ae5a4SJacob Faibussowitsch {
6167b6e2003SPierre Jolivet   PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6177b6e2003SPierre Jolivet   Mat         sX, sY;
6187b6e2003SPierre Jolivet 
6197b6e2003SPierre Jolivet   PetscFunctionBegin;
6207b6e2003SPierre Jolivet   /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6217b6e2003SPierre Jolivet      matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6227b6e2003SPierre Jolivet      of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6239566063dSJacob Faibussowitsch   PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6249566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(X, &sX));
6259566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLocalMatrix(Y, &sY));
6269566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6287b6e2003SPierre Jolivet }
6297b6e2003SPierre Jolivet 
630d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
631d71ae5a4SJacob Faibussowitsch {
6324b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6334b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
634d9ca1df4SBarry Smith   PetscScalar            *y_array;
635d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6364b9ad928SBarry Smith   PC                      subpc;
6374b9ad928SBarry Smith 
6384b9ad928SBarry Smith   PetscFunctionBegin;
6394b9ad928SBarry Smith   /*
6404b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6414b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6424b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6434b9ad928SBarry Smith     for the sequential solve.
6444b9ad928SBarry Smith   */
6459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6479566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6489566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6494b9ad928SBarry Smith   /* apply the symmetric left portion of the inner PC operator */
6504b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6519566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6529566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6539566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
6549566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
6559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6584b9ad928SBarry Smith }
6594b9ad928SBarry Smith 
660d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
661d71ae5a4SJacob Faibussowitsch {
6624b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6634b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
664d9ca1df4SBarry Smith   PetscScalar            *y_array;
665d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6664b9ad928SBarry Smith   PC                      subpc;
6674b9ad928SBarry Smith 
6684b9ad928SBarry Smith   PetscFunctionBegin;
6694b9ad928SBarry Smith   /*
6704b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
6714b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
6724b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
6734b9ad928SBarry Smith     for the sequential solve.
6744b9ad928SBarry Smith   */
6759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
6769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
6779566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
6789566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
6794b9ad928SBarry Smith 
6804b9ad928SBarry Smith   /* apply the symmetric right portion of the inner PC operator */
6814b9ad928SBarry Smith   /* note this by-passes the inner KSP and its options completely */
6824b9ad928SBarry Smith 
6839566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6849566063dSJacob Faibussowitsch   PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6854b9ad928SBarry Smith 
68611e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->x));
68711e4f71cSBarry Smith   PetscCall(VecResetArray(bjac->y));
6889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
6899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6914b9ad928SBarry Smith }
6924b9ad928SBarry Smith 
693d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
694d71ae5a4SJacob Faibussowitsch {
6954b9ad928SBarry Smith   PC_BJacobi             *jac  = (PC_BJacobi *)pc->data;
6964b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
697d9ca1df4SBarry Smith   PetscScalar            *y_array;
698d9ca1df4SBarry Smith   const PetscScalar      *x_array;
6994b9ad928SBarry Smith 
7004b9ad928SBarry Smith   PetscFunctionBegin;
7014b9ad928SBarry Smith   /*
7024b9ad928SBarry Smith       The VecPlaceArray() is to avoid having to copy the
7034b9ad928SBarry Smith     y vector into the bjac->x vector. The reason for
7044b9ad928SBarry Smith     the bjac->x vector is that we need a sequential vector
7054b9ad928SBarry Smith     for the sequential solve.
7064b9ad928SBarry Smith   */
7079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &x_array));
7089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_array));
7099566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->x, x_array));
7109566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(bjac->y, y_array));
7119566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7129566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
7139566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->x));
7149566063dSJacob Faibussowitsch   PetscCall(VecResetArray(bjac->y));
7159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &x_array));
7169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_array));
7173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7184b9ad928SBarry Smith }
7194b9ad928SBarry Smith 
720d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
721d71ae5a4SJacob Faibussowitsch {
7224b9ad928SBarry Smith   PC_BJacobi             *jac = (PC_BJacobi *)pc->data;
72369a612a9SBarry Smith   PetscInt                m;
7244b9ad928SBarry Smith   KSP                     ksp;
7254b9ad928SBarry Smith   PC_BJacobi_Singleblock *bjac;
726de0953f6SBarry Smith   PetscBool               wasSetup = PETSC_TRUE;
7273f6dc190SJunchao Zhang   VecType                 vectype;
728ea41da7aSPierre Jolivet   const char             *prefix;
7294b9ad928SBarry Smith 
7304b9ad928SBarry Smith   PetscFunctionBegin;
7314b9ad928SBarry Smith   if (!pc->setupcalled) {
732c2efdce3SBarry Smith     if (!jac->ksp) {
733302a9ddcSMatthew Knepley       wasSetup = PETSC_FALSE;
7342fa5cd67SKarl Rupp 
7359566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7369566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7379566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7389566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPPREONLY));
7399566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
7409566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7419566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7424b9ad928SBarry Smith 
743e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Singleblock;
7444b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Singleblock;
7454b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Singleblock;
7467b6e2003SPierre Jolivet       pc->ops->matapply            = PCMatApply_BJacobi_Singleblock;
7474b9ad928SBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Singleblock;
7484b9ad928SBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7494b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Singleblock;
7504b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Singleblock;
7514b9ad928SBarry Smith 
7529566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(1, &jac->ksp));
7534b9ad928SBarry Smith       jac->ksp[0] = ksp;
754c2efdce3SBarry Smith 
7554dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
7564b9ad928SBarry Smith       jac->data = (void *)bjac;
7574b9ad928SBarry Smith     } else {
758c2efdce3SBarry Smith       ksp  = jac->ksp[0];
759c2efdce3SBarry Smith       bjac = (PC_BJacobi_Singleblock *)jac->data;
760c2efdce3SBarry Smith     }
761c2efdce3SBarry Smith 
762c2efdce3SBarry Smith     /*
763c2efdce3SBarry Smith       The reason we need to generate these vectors is to serve
764c2efdce3SBarry Smith       as the right-hand side and solution vector for the solve on the
765c2efdce3SBarry Smith       block. We do not need to allocate space for the vectors since
766c2efdce3SBarry Smith       that is provided via VecPlaceArray() just before the call to
767c2efdce3SBarry Smith       KSPSolve() on the block.
768c2efdce3SBarry Smith     */
7699566063dSJacob Faibussowitsch     PetscCall(MatGetSize(pmat, &m, &m));
7709566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7719566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7729566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
7739566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->x, vectype));
7749566063dSJacob Faibussowitsch     PetscCall(VecSetType(bjac->y, vectype));
775c2efdce3SBarry Smith   } else {
7764b9ad928SBarry Smith     ksp  = jac->ksp[0];
7774b9ad928SBarry Smith     bjac = (PC_BJacobi_Singleblock *)jac->data;
7784b9ad928SBarry Smith   }
7799566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
78049517cdeSBarry Smith   if (pc->useAmat) {
7819566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, mat, pmat));
7829566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mat, prefix));
7834b9ad928SBarry Smith   } else {
7849566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ksp, pmat, pmat));
7854b9ad928SBarry Smith   }
7869566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(pmat, prefix));
78790f1c854SHong Zhang   if (!wasSetup && pc->setfromoptionscalled) {
788248bfaf0SJed Brown     /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
7899566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ksp));
790302a9ddcSMatthew Knepley   }
7913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7924b9ad928SBarry Smith }
7934b9ad928SBarry Smith 
794d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
795d71ae5a4SJacob Faibussowitsch {
7964b9ad928SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
7974b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
79869a612a9SBarry Smith   PetscInt               i;
7994b9ad928SBarry Smith 
8004b9ad928SBarry Smith   PetscFunctionBegin;
801e91c6855SBarry Smith   if (bjac && bjac->pmat) {
8029566063dSJacob Faibussowitsch     PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
80348a46eb9SPierre Jolivet     if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
804e91c6855SBarry Smith   }
8054b9ad928SBarry Smith 
8064b9ad928SBarry Smith   for (i = 0; i < jac->n_local; i++) {
8079566063dSJacob Faibussowitsch     PetscCall(KSPReset(jac->ksp[i]));
808e91c6855SBarry Smith     if (bjac && bjac->x) {
8099566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->x[i]));
8109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&bjac->y[i]));
8119566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&bjac->is[i]));
8124b9ad928SBarry Smith     }
813e91c6855SBarry Smith   }
8149566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->l_lens));
8159566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->g_lens));
8163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
817e91c6855SBarry Smith }
818e91c6855SBarry Smith 
819d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
820d71ae5a4SJacob Faibussowitsch {
821e91c6855SBarry Smith   PC_BJacobi            *jac  = (PC_BJacobi *)pc->data;
822c2efdce3SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
823e91c6855SBarry Smith   PetscInt               i;
824e91c6855SBarry Smith 
825e91c6855SBarry Smith   PetscFunctionBegin;
8269566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiblock(pc));
827c2efdce3SBarry Smith   if (bjac) {
8289566063dSJacob Faibussowitsch     PetscCall(PetscFree2(bjac->x, bjac->y));
8299566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->starts));
8309566063dSJacob Faibussowitsch     PetscCall(PetscFree(bjac->is));
831c2efdce3SBarry Smith   }
8329566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->data));
83348a46eb9SPierre Jolivet   for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8349566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
8352e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
8363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8374b9ad928SBarry Smith }
8384b9ad928SBarry Smith 
839d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
840d71ae5a4SJacob Faibussowitsch {
8414b9ad928SBarry Smith   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
84269a612a9SBarry Smith   PetscInt           i, n_local = jac->n_local;
843539c167fSBarry Smith   KSPConvergedReason reason;
8444b9ad928SBarry Smith 
8454b9ad928SBarry Smith   PetscFunctionBegin;
8464b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8479566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(jac->ksp[i]));
8489566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
849ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8504b9ad928SBarry Smith   }
8513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8524b9ad928SBarry Smith }
8534b9ad928SBarry Smith 
854d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
855d71ae5a4SJacob Faibussowitsch {
8564b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
85769a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
8584b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
859d9ca1df4SBarry Smith   PetscScalar           *yin;
860d9ca1df4SBarry Smith   const PetscScalar     *xin;
86158ebbce7SBarry Smith 
8624b9ad928SBarry Smith   PetscFunctionBegin;
8639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
8649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
8654b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
8664b9ad928SBarry Smith     /*
8674b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
8684b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
8694b9ad928SBarry Smith        the global vector.
8704b9ad928SBarry Smith     */
8719566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8729566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8734b9ad928SBarry Smith 
8749566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8759566063dSJacob Faibussowitsch     PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8769566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8779566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
878d11f3a42SBarry Smith 
8799566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
8809566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
8814b9ad928SBarry Smith   }
8829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
8839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
8843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8854b9ad928SBarry Smith }
8864b9ad928SBarry Smith 
887f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
888f4d694ddSBarry Smith {
889f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
890f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
891f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
892f4d694ddSBarry Smith   PetscScalar           *yin;
893f4d694ddSBarry Smith   const PetscScalar     *xin;
894f4d694ddSBarry Smith   PC                     subpc;
895f4d694ddSBarry Smith 
896f4d694ddSBarry Smith   PetscFunctionBegin;
897f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
898f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
899f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
9004b9ad928SBarry Smith     /*
901f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
902f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
903f4d694ddSBarry Smith        the global vector.
9044b9ad928SBarry Smith     */
905f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
906f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
907f4d694ddSBarry Smith 
908f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
909f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
910f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
911f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
912f4d694ddSBarry Smith     PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
913f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
914f4d694ddSBarry Smith 
915f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
916f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
917f4d694ddSBarry Smith   }
918f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
919f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
921f4d694ddSBarry Smith }
922f4d694ddSBarry Smith 
923f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
924f4d694ddSBarry Smith {
925f4d694ddSBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
926f4d694ddSBarry Smith   PetscInt               i, n_local = jac->n_local;
927f4d694ddSBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
928f4d694ddSBarry Smith   PetscScalar           *yin;
929f4d694ddSBarry Smith   const PetscScalar     *xin;
930f4d694ddSBarry Smith   PC                     subpc;
931f4d694ddSBarry Smith 
932f4d694ddSBarry Smith   PetscFunctionBegin;
933f4d694ddSBarry Smith   PetscCall(VecGetArrayRead(x, &xin));
934f4d694ddSBarry Smith   PetscCall(VecGetArray(y, &yin));
935f4d694ddSBarry Smith   for (i = 0; i < n_local; i++) {
936f4d694ddSBarry Smith     /*
937f4d694ddSBarry Smith        To avoid copying the subvector from x into a workspace we instead
938f4d694ddSBarry Smith        make the workspace vector array point to the subpart of the array of
939f4d694ddSBarry Smith        the global vector.
940f4d694ddSBarry Smith     */
941f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
942f4d694ddSBarry Smith     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
943f4d694ddSBarry Smith 
944f4d694ddSBarry Smith     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
945f4d694ddSBarry Smith     /* apply the symmetric left portion of the inner PC operator */
946f4d694ddSBarry Smith     /* note this by-passes the inner KSP and its options completely */
947f4d694ddSBarry Smith     PetscCall(KSPGetPC(jac->ksp[i], &subpc));
948f4d694ddSBarry Smith     PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
949f4d694ddSBarry Smith     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
950f4d694ddSBarry Smith 
951f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->x[i]));
952f4d694ddSBarry Smith     PetscCall(VecResetArray(bjac->y[i]));
953f4d694ddSBarry Smith   }
954f4d694ddSBarry Smith   PetscCall(VecRestoreArrayRead(x, &xin));
955f4d694ddSBarry Smith   PetscCall(VecRestoreArray(y, &yin));
9563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
957f4d694ddSBarry Smith }
958f4d694ddSBarry Smith 
959d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
960d71ae5a4SJacob Faibussowitsch {
9614b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
96269a612a9SBarry Smith   PetscInt               i, n_local = jac->n_local;
9634b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
964d9ca1df4SBarry Smith   PetscScalar           *yin;
965d9ca1df4SBarry Smith   const PetscScalar     *xin;
9664b9ad928SBarry Smith 
9674b9ad928SBarry Smith   PetscFunctionBegin;
9689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xin));
9699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yin));
9704b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
9714b9ad928SBarry Smith     /*
9724b9ad928SBarry Smith        To avoid copying the subvector from x into a workspace we instead
9734b9ad928SBarry Smith        make the workspace vector array point to the subpart of the array of
9744b9ad928SBarry Smith        the global vector.
9754b9ad928SBarry Smith     */
9769566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9779566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9784b9ad928SBarry Smith 
9799566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
9809566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
9819566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
9829566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
983a7ff49e8SJed Brown 
9849566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->x[i]));
9859566063dSJacob Faibussowitsch     PetscCall(VecResetArray(bjac->y[i]));
9864b9ad928SBarry Smith   }
9879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xin));
9889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yin));
9893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9904b9ad928SBarry Smith }
9914b9ad928SBarry Smith 
992d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
993d71ae5a4SJacob Faibussowitsch {
9944b9ad928SBarry Smith   PC_BJacobi            *jac = (PC_BJacobi *)pc->data;
99569a612a9SBarry Smith   PetscInt               m, n_local, N, M, start, i;
996ea41da7aSPierre Jolivet   const char            *prefix;
9974b9ad928SBarry Smith   KSP                    ksp;
9984b9ad928SBarry Smith   Vec                    x, y;
9994b9ad928SBarry Smith   PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
10004b9ad928SBarry Smith   PC                     subpc;
10014b9ad928SBarry Smith   IS                     is;
1002434a97beSBrad Aagaard   MatReuse               scall;
10033f6dc190SJunchao Zhang   VecType                vectype;
10044b9ad928SBarry Smith 
10054b9ad928SBarry Smith   PetscFunctionBegin;
10069566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10074b9ad928SBarry Smith 
10084b9ad928SBarry Smith   n_local = jac->n_local;
10094b9ad928SBarry Smith 
101049517cdeSBarry Smith   if (pc->useAmat) {
1011ace3abfcSBarry Smith     PetscBool same;
10129566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
101328b400f6SJacob Faibussowitsch     PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10144b9ad928SBarry Smith   }
10154b9ad928SBarry Smith 
10164b9ad928SBarry Smith   if (!pc->setupcalled) {
10174b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
1018c2efdce3SBarry Smith 
1019c2efdce3SBarry Smith     if (!jac->ksp) {
1020e91c6855SBarry Smith       pc->ops->reset               = PCReset_BJacobi_Multiblock;
10214b9ad928SBarry Smith       pc->ops->destroy             = PCDestroy_BJacobi_Multiblock;
10224b9ad928SBarry Smith       pc->ops->apply               = PCApply_BJacobi_Multiblock;
10237b6e2003SPierre Jolivet       pc->ops->matapply            = NULL;
1024f4d694ddSBarry Smith       pc->ops->applysymmetricleft  = PCApplySymmetricLeft_BJacobi_Multiblock;
1025f4d694ddSBarry Smith       pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10264b9ad928SBarry Smith       pc->ops->applytranspose      = PCApplyTranspose_BJacobi_Multiblock;
10274b9ad928SBarry Smith       pc->ops->setuponblocks       = PCSetUpOnBlocks_BJacobi_Multiblock;
10284b9ad928SBarry Smith 
10294dfa11a4SJacob Faibussowitsch       PetscCall(PetscNew(&bjac));
10309566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &jac->ksp));
10319566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10329566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->starts));
10334b9ad928SBarry Smith 
10344b9ad928SBarry Smith       jac->data = (void *)bjac;
10359566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n_local, &bjac->is));
10364b9ad928SBarry Smith 
10374b9ad928SBarry Smith       for (i = 0; i < n_local; i++) {
10389566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10399566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10409566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10419566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
10429566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
10439566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
10449566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10459566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10462fa5cd67SKarl Rupp 
1047c2efdce3SBarry Smith         jac->ksp[i] = ksp;
1048c2efdce3SBarry Smith       }
1049c2efdce3SBarry Smith     } else {
1050c2efdce3SBarry Smith       bjac = (PC_BJacobi_Multiblock *)jac->data;
1051c2efdce3SBarry Smith     }
10524b9ad928SBarry Smith 
1053c2efdce3SBarry Smith     start = 0;
10549566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(pmat, &vectype));
1055c2efdce3SBarry Smith     for (i = 0; i < n_local; i++) {
10564b9ad928SBarry Smith       m = jac->l_lens[i];
10574b9ad928SBarry Smith       /*
10584b9ad928SBarry Smith       The reason we need to generate these vectors is to serve
10594b9ad928SBarry Smith       as the right-hand side and solution vector for the solve on the
10604b9ad928SBarry Smith       block. We do not need to allocate space for the vectors since
10614b9ad928SBarry Smith       that is provided via VecPlaceArray() just before the call to
10624b9ad928SBarry Smith       KSPSolve() on the block.
10634b9ad928SBarry Smith 
10644b9ad928SBarry Smith       */
10659566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10669566063dSJacob Faibussowitsch       PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10679566063dSJacob Faibussowitsch       PetscCall(VecSetType(x, vectype));
10689566063dSJacob Faibussowitsch       PetscCall(VecSetType(y, vectype));
10692fa5cd67SKarl Rupp 
10704b9ad928SBarry Smith       bjac->x[i]      = x;
10714b9ad928SBarry Smith       bjac->y[i]      = y;
10724b9ad928SBarry Smith       bjac->starts[i] = start;
10734b9ad928SBarry Smith 
10749566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
10754b9ad928SBarry Smith       bjac->is[i] = is;
10764b9ad928SBarry Smith 
10774b9ad928SBarry Smith       start += m;
10784b9ad928SBarry Smith     }
10794b9ad928SBarry Smith   } else {
10804b9ad928SBarry Smith     bjac = (PC_BJacobi_Multiblock *)jac->data;
10814b9ad928SBarry Smith     /*
10824b9ad928SBarry Smith        Destroy the blocks from the previous iteration
10834b9ad928SBarry Smith     */
10844b9ad928SBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
10859566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
108648a46eb9SPierre Jolivet       if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
10874b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
10882fa5cd67SKarl Rupp     } else scall = MAT_REUSE_MATRIX;
10894b9ad928SBarry Smith   }
10904b9ad928SBarry Smith 
10919566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
109248a46eb9SPierre Jolivet   if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
10934b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
10944b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
10959566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
10964b9ad928SBarry Smith 
10974b9ad928SBarry Smith   for (i = 0; i < n_local; i++) {
10989566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
109949517cdeSBarry Smith     if (pc->useAmat) {
11009566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11019566063dSJacob Faibussowitsch       PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11024b9ad928SBarry Smith     } else {
11039566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11044b9ad928SBarry Smith     }
11059566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
110648a46eb9SPierre Jolivet     if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
110790f1c854SHong Zhang   }
11083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11094b9ad928SBarry Smith }
11105a7e1895SHong Zhang 
11115a7e1895SHong Zhang /*
1112fd0b8932SBarry Smith       These are for a single block with multiple processes
11135a7e1895SHong Zhang */
1114d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1115d71ae5a4SJacob Faibussowitsch {
1116fd0b8932SBarry Smith   PC_BJacobi        *jac    = (PC_BJacobi *)pc->data;
1117fd0b8932SBarry Smith   KSP                subksp = jac->ksp[0];
1118fd0b8932SBarry Smith   KSPConvergedReason reason;
1119fd0b8932SBarry Smith 
1120fd0b8932SBarry Smith   PetscFunctionBegin;
11219566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(subksp));
11229566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(subksp, &reason));
1123ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1125fd0b8932SBarry Smith }
1126fd0b8932SBarry Smith 
1127d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1128d71ae5a4SJacob Faibussowitsch {
11295a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11305a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11315a7e1895SHong Zhang 
11325a7e1895SHong Zhang   PetscFunctionBegin;
11339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->ysub));
11349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mpjac->xsub));
11359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpjac->submats));
11369566063dSJacob Faibussowitsch   if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1138e91c6855SBarry Smith }
1139e91c6855SBarry Smith 
1140d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1141d71ae5a4SJacob Faibussowitsch {
1142e91c6855SBarry Smith   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
1143e91c6855SBarry Smith   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1144e91c6855SBarry Smith 
1145e91c6855SBarry Smith   PetscFunctionBegin;
11469566063dSJacob Faibussowitsch   PetscCall(PCReset_BJacobi_Multiproc(pc));
11479566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp[0]));
11489566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->ksp));
11499566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11505a7e1895SHong Zhang 
11519566063dSJacob Faibussowitsch   PetscCall(PetscFree(mpjac));
11522e956fe4SStefano Zampini   PetscCall(PCDestroy_BJacobi(pc));
11533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11545a7e1895SHong Zhang }
11555a7e1895SHong Zhang 
1156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1157d71ae5a4SJacob Faibussowitsch {
11585a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
11595a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1160d9ca1df4SBarry Smith   PetscScalar          *yarray;
1161d9ca1df4SBarry Smith   const PetscScalar    *xarray;
1162539c167fSBarry Smith   KSPConvergedReason    reason;
11635a7e1895SHong Zhang 
11645a7e1895SHong Zhang   PetscFunctionBegin;
11655a7e1895SHong Zhang   /* place x's and y's local arrays into xsub and ysub */
11669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xarray));
11679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yarray));
11689566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->xsub, xarray));
11699566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(mpjac->ysub, yarray));
11705a7e1895SHong Zhang 
11715a7e1895SHong Zhang   /* apply preconditioner on each matrix block */
11729566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11739566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
11749566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
11759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
11769566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1177ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1178250cf88bSHong Zhang 
11799566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->xsub));
11809566063dSJacob Faibussowitsch   PetscCall(VecResetArray(mpjac->ysub));
11819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xarray));
11829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yarray));
11833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11845a7e1895SHong Zhang }
11855a7e1895SHong Zhang 
1186d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1187d71ae5a4SJacob Faibussowitsch {
11887b6e2003SPierre Jolivet   PC_BJacobi        *jac = (PC_BJacobi *)pc->data;
11897b6e2003SPierre Jolivet   KSPConvergedReason reason;
11907b6e2003SPierre Jolivet   Mat                sX, sY;
11917b6e2003SPierre Jolivet   const PetscScalar *x;
11927b6e2003SPierre Jolivet   PetscScalar       *y;
11937b6e2003SPierre Jolivet   PetscInt           m, N, lda, ldb;
11947b6e2003SPierre Jolivet 
11957b6e2003SPierre Jolivet   PetscFunctionBegin;
11967b6e2003SPierre Jolivet   /* apply preconditioner on each matrix block */
11979566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, NULL));
11989566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
11999566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
12009566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(Y, &ldb));
12019566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(X, &x));
12029566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(Y, &y));
12039566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12049566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12059566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sX, lda));
12069566063dSJacob Faibussowitsch   PetscCall(MatDenseSetLDA(sY, ldb));
12079566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12089566063dSJacob Faibussowitsch   PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12099566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sY));
12129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sX));
12139566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12149566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(X, &x));
12159566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1216ad540459SPierre Jolivet   if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12187b6e2003SPierre Jolivet }
12197b6e2003SPierre Jolivet 
1220d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1221d71ae5a4SJacob Faibussowitsch {
12225a7e1895SHong Zhang   PC_BJacobi           *jac   = (PC_BJacobi *)pc->data;
12235a7e1895SHong Zhang   PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12245a7e1895SHong Zhang   PetscInt              m, n;
1225ce94432eSBarry Smith   MPI_Comm              comm, subcomm = 0;
12265a7e1895SHong Zhang   const char           *prefix;
1227de0953f6SBarry Smith   PetscBool             wasSetup = PETSC_TRUE;
12283f6dc190SJunchao Zhang   VecType               vectype;
12291f62890fSKarl Rupp 
12305a7e1895SHong Zhang   PetscFunctionBegin;
12319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
123208401ef6SPierre Jolivet   PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12335a7e1895SHong Zhang   jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12345a7e1895SHong Zhang   if (!pc->setupcalled) {
1235de0953f6SBarry Smith     wasSetup = PETSC_FALSE;
12364dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mpjac));
12375a7e1895SHong Zhang     jac->data = (void *)mpjac;
12385a7e1895SHong Zhang 
12395a7e1895SHong Zhang     /* initialize datastructure mpjac */
12405a7e1895SHong Zhang     if (!jac->psubcomm) {
12415a7e1895SHong Zhang       /* Create default contiguous subcommunicatiors if user does not provide them */
12429566063dSJacob Faibussowitsch       PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12439566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12449566063dSJacob Faibussowitsch       PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12455a7e1895SHong Zhang     }
12465a7e1895SHong Zhang     mpjac->psubcomm = jac->psubcomm;
1247306c2d5bSBarry Smith     subcomm         = PetscSubcommChild(mpjac->psubcomm);
12485a7e1895SHong Zhang 
12495a7e1895SHong Zhang     /* Get matrix blocks of pmat */
12509566063dSJacob Faibussowitsch     PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12515a7e1895SHong Zhang 
12525a7e1895SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
12539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &jac->ksp));
12549566063dSJacob Faibussowitsch     PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12559566063dSJacob Faibussowitsch     PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12569566063dSJacob Faibussowitsch     PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12579566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12589566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12595a7e1895SHong Zhang 
12609566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
12619566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
12629566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
12639566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
12649566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
12655a7e1895SHong Zhang 
12665a7e1895SHong Zhang     /* create dummy vectors xsub and ysub */
12679566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
12689566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
12699566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
12709566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(mpjac->submats, &vectype));
12719566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->xsub, vectype));
12729566063dSJacob Faibussowitsch     PetscCall(VecSetType(mpjac->ysub, vectype));
12735a7e1895SHong Zhang 
1274fd0b8932SBarry Smith     pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1275e91c6855SBarry Smith     pc->ops->reset         = PCReset_BJacobi_Multiproc;
12765a7e1895SHong Zhang     pc->ops->destroy       = PCDestroy_BJacobi_Multiproc;
12775a7e1895SHong Zhang     pc->ops->apply         = PCApply_BJacobi_Multiproc;
12787b6e2003SPierre Jolivet     pc->ops->matapply      = PCMatApply_BJacobi_Multiproc;
1279fc08c53fSHong Zhang   } else { /* pc->setupcalled */
1280306c2d5bSBarry Smith     subcomm = PetscSubcommChild(mpjac->psubcomm);
12819e0ae222SHong Zhang     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
12825a7e1895SHong Zhang       /* destroy old matrix blocks, then get new matrix blocks */
12839566063dSJacob Faibussowitsch       if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
12849566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1285fc08c53fSHong Zhang     } else {
12869566063dSJacob Faibussowitsch       PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
12875a7e1895SHong Zhang     }
12889566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12895a7e1895SHong Zhang   }
12905a7e1895SHong Zhang 
129148a46eb9SPierre Jolivet   if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
12923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12935a7e1895SHong Zhang }
1294