14b9ad928SBarry Smith /*
24b9ad928SBarry Smith Defines a block Jacobi preconditioner.
34b9ad928SBarry Smith */
400e125f8SBarry Smith
500e125f8SBarry Smith #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/
64b9ad928SBarry Smith
76849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat);
86849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat);
95a7e1895SHong Zhang static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC);
104b9ad928SBarry Smith
PCSetUp_BJacobi(PC pc)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi(PC pc)
12d71ae5a4SJacob Faibussowitsch {
134b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
144b9ad928SBarry Smith Mat mat = pc->mat, pmat = pc->pmat;
15976e8c5aSLisandro Dalcin PetscBool hasop;
1669a612a9SBarry Smith PetscInt N, M, start, i, sum, end;
1769a612a9SBarry Smith PetscInt bs, i_start = -1, i_end = -1;
1869a612a9SBarry Smith PetscMPIInt rank, size;
194b9ad928SBarry Smith
204b9ad928SBarry Smith PetscFunctionBegin;
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
239566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
249566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(pc->pmat, &bs));
254b9ad928SBarry Smith
265a7e1895SHong Zhang if (jac->n > 0 && jac->n < size) {
279566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Multiproc(pc));
283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
295a7e1895SHong Zhang }
305a7e1895SHong Zhang
31f1580f4eSBarry Smith /* Determines the number of blocks assigned to each processor */
324b9ad928SBarry Smith /* local block count given */
334b9ad928SBarry Smith if (jac->n_local > 0 && jac->n < 0) {
34462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
354b9ad928SBarry Smith if (jac->l_lens) { /* check that user set these correctly */
364b9ad928SBarry Smith sum = 0;
374b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) {
3808401ef6SPierre Jolivet PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
394b9ad928SBarry Smith sum += jac->l_lens[i];
404b9ad928SBarry Smith }
4108401ef6SPierre Jolivet PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly");
424b9ad928SBarry Smith } else {
439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
442fa5cd67SKarl Rupp for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
454b9ad928SBarry Smith }
464b9ad928SBarry Smith } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */
474b9ad928SBarry Smith /* global blocks given: determine which ones are local */
484b9ad928SBarry Smith if (jac->g_lens) {
494b9ad928SBarry Smith /* check if the g_lens is has valid entries */
504b9ad928SBarry Smith for (i = 0; i < jac->n; i++) {
517827d75bSBarry Smith PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed");
5208401ef6SPierre Jolivet PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout");
534b9ad928SBarry Smith }
544b9ad928SBarry Smith if (size == 1) {
554b9ad928SBarry Smith jac->n_local = jac->n;
569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
579566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local));
584b9ad928SBarry Smith /* check that user set these correctly */
594b9ad928SBarry Smith sum = 0;
604b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i];
6108401ef6SPierre Jolivet PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly");
624b9ad928SBarry Smith } else {
639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end));
64aaa8cc7dSPierre Jolivet /* loop over blocks determining first one owned by me */
654b9ad928SBarry Smith sum = 0;
664b9ad928SBarry Smith for (i = 0; i < jac->n + 1; i++) {
679371c9d4SSatish Balay if (sum == start) {
689371c9d4SSatish Balay i_start = i;
699371c9d4SSatish Balay goto start_1;
709371c9d4SSatish Balay }
714b9ad928SBarry Smith if (i < jac->n) sum += jac->g_lens[i];
724b9ad928SBarry Smith }
73e7e72b3dSBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
744b9ad928SBarry Smith start_1:
754b9ad928SBarry Smith for (i = i_start; i < jac->n + 1; i++) {
769371c9d4SSatish Balay if (sum == end) {
779371c9d4SSatish Balay i_end = i;
789371c9d4SSatish Balay goto end_1;
799371c9d4SSatish Balay }
804b9ad928SBarry Smith if (i < jac->n) sum += jac->g_lens[i];
814b9ad928SBarry Smith }
82e7e72b3dSBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout");
834b9ad928SBarry Smith end_1:
844b9ad928SBarry Smith jac->n_local = i_end - i_start;
859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
869566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local));
874b9ad928SBarry Smith }
884b9ad928SBarry Smith } else { /* no global blocks given, determine then using default layout */
894b9ad928SBarry Smith jac->n_local = jac->n / size + ((jac->n % size) > rank);
909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
914b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) {
924b9ad928SBarry Smith jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs;
937827d75bSBarry Smith PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given");
944b9ad928SBarry Smith }
954b9ad928SBarry Smith }
964b9ad928SBarry Smith } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */
974b9ad928SBarry Smith jac->n = size;
984b9ad928SBarry Smith jac->n_local = 1;
999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->l_lens));
1004b9ad928SBarry Smith jac->l_lens[0] = M;
1017271b979SBarry Smith } else { /* jac->n > 0 && jac->n_local > 0 */
1027271b979SBarry Smith if (!jac->l_lens) {
1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens));
1042fa5cd67SKarl Rupp for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i));
1057271b979SBarry Smith }
1064b9ad928SBarry Smith }
10708401ef6SPierre Jolivet PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors");
1084b9ad928SBarry Smith
109f1580f4eSBarry Smith /* Determines mat and pmat */
1109566063dSJacob Faibussowitsch PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop));
111976e8c5aSLisandro Dalcin if (!hasop && size == 1) {
1124b9ad928SBarry Smith mat = pc->mat;
1134b9ad928SBarry Smith pmat = pc->pmat;
1144b9ad928SBarry Smith } else {
11549517cdeSBarry Smith if (pc->useAmat) {
11649517cdeSBarry Smith /* use block from Amat matrix, not Pmat for local MatMult() */
1179566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(pc->mat, &mat));
1184b9ad928SBarry Smith }
119*ac530a7eSPierre Jolivet if (pc->pmat != pc->mat || !pc->useAmat) PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat));
120*ac530a7eSPierre Jolivet else pmat = mat;
1214b9ad928SBarry Smith }
1224b9ad928SBarry Smith
123f1580f4eSBarry Smith /*
1244b9ad928SBarry Smith Setup code depends on the number of blocks
1254b9ad928SBarry Smith */
126cc1d6079SHong Zhang if (jac->n_local == 1) {
1279566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat));
1284b9ad928SBarry Smith } else {
1299566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat));
1304b9ad928SBarry Smith }
1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1324b9ad928SBarry Smith }
1334b9ad928SBarry Smith
1344b9ad928SBarry Smith /* Default destroy, if it has never been setup */
PCDestroy_BJacobi(PC pc)135d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi(PC pc)
136d71ae5a4SJacob Faibussowitsch {
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));
1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith
PCSetFromOptions_BJacobi(PC pc,PetscOptionItems PetscOptionsObject)151ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems PetscOptionsObject)
152d71ae5a4SJacob Faibussowitsch {
1534b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
154248bfaf0SJed Brown PetscInt blocks, i;
155ace3abfcSBarry Smith PetscBool flg;
1564b9ad928SBarry Smith
1574b9ad928SBarry Smith PetscFunctionBegin;
158d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options");
1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg));
1609566063dSJacob Faibussowitsch if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL));
1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg));
1629566063dSJacob Faibussowitsch if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL));
163248bfaf0SJed Brown if (jac->ksp) {
164248bfaf0SJed Brown /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called
165248bfaf0SJed Brown * unless we had already been called. */
16648a46eb9SPierre Jolivet for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i]));
167248bfaf0SJed Brown }
168d0609cedSBarry Smith PetscOptionsHeadEnd();
1693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1704b9ad928SBarry Smith }
1714b9ad928SBarry Smith
1729804daf3SBarry Smith #include <petscdraw.h>
PCView_BJacobi(PC pc,PetscViewer viewer)173d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer)
174d71ae5a4SJacob Faibussowitsch {
1754b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
176cbe18068SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
17769a612a9SBarry Smith PetscMPIInt rank;
17869a612a9SBarry Smith PetscInt i;
1799f196a02SMartin Diehl PetscBool isascii, isstring, isdraw;
1804b9ad928SBarry Smith PetscViewer sviewer;
181020d6619SPierre Jolivet PetscViewerFormat format;
182020d6619SPierre Jolivet const char *prefix;
1834b9ad928SBarry Smith
1844b9ad928SBarry Smith PetscFunctionBegin;
1859f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1869566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
1879566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
1889f196a02SMartin Diehl if (isascii) {
18948a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n));
19063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " number of blocks = %" PetscInt_FMT "\n", jac->n));
1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
1929566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format));
193020d6619SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
1949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
1959566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix));
1969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
197933995ecSHong Zhang if (jac->ksp && !jac->psubcomm) {
1989566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
199dd400576SPatrick Sanan if (rank == 0) {
200b4025f61SBarry Smith PetscCall(PetscViewerASCIIPushTab(sviewer));
2019566063dSJacob Faibussowitsch PetscCall(KSPView(jac->ksp[0], sviewer));
202b4025f61SBarry Smith PetscCall(PetscViewerASCIIPopTab(sviewer));
2034b9ad928SBarry Smith }
2049566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
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) {
210b4025f61SBarry Smith PetscCall(PetscViewerASCIIPushTab(sviewer));
211f4f49eeaSPierre Jolivet PetscCall(KSPView(*jac->ksp, sviewer));
212b4025f61SBarry Smith PetscCall(PetscViewerASCIIPopTab(sviewer));
213cbe18068SHong Zhang }
2149566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer));
2159530cbd7SBarry Smith /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */
2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer));
217e4de9e1dSBarry Smith }
218e4de9e1dSBarry Smith } else {
2194814766eSBarry Smith PetscInt n_global;
220462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc)));
2219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following KSP and PC objects:\n"));
2239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
2249566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
225b4025f61SBarry Smith PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local));
226dd2fa690SBarry Smith for (i = 0; i < jac->n_local; i++) {
227b4025f61SBarry Smith PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i));
2289566063dSJacob Faibussowitsch PetscCall(KSPView(jac->ksp[i], sviewer));
229b4025f61SBarry Smith PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
2304b9ad928SBarry Smith }
2319566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2344b9ad928SBarry Smith }
2354b9ad928SBarry Smith } else if (isstring) {
23663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n));
2379566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
2389566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer));
2399566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
240d9884438SBarry Smith } else if (isdraw) {
241d9884438SBarry Smith PetscDraw draw;
242d9884438SBarry Smith char str[25];
243d9884438SBarry Smith PetscReal x, y, bottom, h;
244d9884438SBarry Smith
2459566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
2469566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
24763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n));
2489566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h));
249d9884438SBarry Smith bottom = y - h;
2509566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
251d9884438SBarry Smith /* warning the communicator on viewer is different then on ksp in parallel */
2529566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer));
2539566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw));
2544b9ad928SBarry Smith }
2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2564b9ad928SBarry Smith }
2574b9ad928SBarry Smith
PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt * n_local,PetscInt * first_local,KSP ** ksp)258d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
259d71ae5a4SJacob Faibussowitsch {
260feb237baSPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2614b9ad928SBarry Smith
2624b9ad928SBarry Smith PetscFunctionBegin;
26328b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first");
2644b9ad928SBarry Smith
2654b9ad928SBarry Smith if (n_local) *n_local = jac->n_local;
2664b9ad928SBarry Smith if (first_local) *first_local = jac->first_local;
267020d6619SPierre Jolivet if (ksp) *ksp = jac->ksp;
2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2694b9ad928SBarry Smith }
2704b9ad928SBarry Smith
PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt * lens)271f9663b93SPierre Jolivet static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt *lens)
272d71ae5a4SJacob Faibussowitsch {
2734b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2744b9ad928SBarry Smith
2754b9ad928SBarry Smith PetscFunctionBegin;
276371d2eb7SMartin Diehl PetscCheck(!pc->setupcalled || jac->n == blocks, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called");
2774b9ad928SBarry Smith jac->n = blocks;
2780a545947SLisandro Dalcin if (!lens) jac->g_lens = NULL;
2792fa5cd67SKarl Rupp else {
2809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(blocks, &jac->g_lens));
2819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->g_lens, lens, blocks));
2824b9ad928SBarry Smith }
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2844b9ad928SBarry Smith }
2854b9ad928SBarry Smith
PCBJacobiGetTotalBlocks_BJacobi(PC pc,PetscInt * blocks,const PetscInt * lens[])286d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
287d71ae5a4SJacob Faibussowitsch {
2884b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
2894b9ad928SBarry Smith
2904b9ad928SBarry Smith PetscFunctionBegin;
2914b9ad928SBarry Smith *blocks = jac->n;
2924b9ad928SBarry Smith if (lens) *lens = jac->g_lens;
2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2944b9ad928SBarry Smith }
2954b9ad928SBarry Smith
PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[])296d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[])
297d71ae5a4SJacob Faibussowitsch {
2984b9ad928SBarry Smith PC_BJacobi *jac;
2994b9ad928SBarry Smith
3004b9ad928SBarry Smith PetscFunctionBegin;
3014b9ad928SBarry Smith jac = (PC_BJacobi *)pc->data;
3024b9ad928SBarry Smith
3034b9ad928SBarry Smith jac->n_local = blocks;
3040a545947SLisandro Dalcin if (!lens) jac->l_lens = NULL;
3052fa5cd67SKarl Rupp else {
3069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(blocks, &jac->l_lens));
3079566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, lens, blocks));
3084b9ad928SBarry Smith }
3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3104b9ad928SBarry Smith }
3114b9ad928SBarry Smith
PCBJacobiGetLocalBlocks_BJacobi(PC pc,PetscInt * blocks,const PetscInt * lens[])312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[])
313d71ae5a4SJacob Faibussowitsch {
3144b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
3154b9ad928SBarry Smith
3164b9ad928SBarry Smith PetscFunctionBegin;
3174b9ad928SBarry Smith *blocks = jac->n_local;
3184b9ad928SBarry Smith if (lens) *lens = jac->l_lens;
3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3204b9ad928SBarry Smith }
3214b9ad928SBarry Smith
3224b9ad928SBarry Smith /*@C
323f1580f4eSBarry Smith PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on
3244b9ad928SBarry Smith this processor.
3254b9ad928SBarry Smith
3266da92b7fSHong Zhang Not Collective
3274b9ad928SBarry Smith
3284b9ad928SBarry Smith Input Parameter:
3294b9ad928SBarry Smith . pc - the preconditioner context
3304b9ad928SBarry Smith
3314b9ad928SBarry Smith Output Parameters:
3320298fd71SBarry Smith + n_local - the number of blocks on this processor, or NULL
3330298fd71SBarry Smith . first_local - the global number of the first block on this processor, or NULL
3344b9ad928SBarry Smith - ksp - the array of KSP contexts
3354b9ad928SBarry Smith
336ce78bad3SBarry Smith Level: advanced
337ce78bad3SBarry Smith
3384b9ad928SBarry Smith Notes:
339f1580f4eSBarry 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
344f1580f4eSBarry Smith You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`.
3454b9ad928SBarry Smith
34636083efbSBarry Smith Fortran Note:
34736083efbSBarry Smith Call `PCBJacobiRestoreSubKSP()` when you no longer need access to the array of `KSP`
34836083efbSBarry Smith
349562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()`
3504b9ad928SBarry Smith @*/
PCBJacobiGetSubKSP(PC pc,PetscInt * n_local,PetscInt * first_local,KSP * ksp[])351d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
352d71ae5a4SJacob Faibussowitsch {
3534b9ad928SBarry Smith PetscFunctionBegin;
3540700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
355cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3574b9ad928SBarry Smith }
3584b9ad928SBarry Smith
3594b9ad928SBarry Smith /*@
3604b9ad928SBarry Smith PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block
3614b9ad928SBarry Smith Jacobi preconditioner.
3624b9ad928SBarry Smith
363c3339decSBarry Smith Collective
3644b9ad928SBarry Smith
3654b9ad928SBarry Smith Input Parameters:
3664b9ad928SBarry Smith + pc - the preconditioner context
3674b9ad928SBarry Smith . blocks - the number of blocks
3684b9ad928SBarry Smith - lens - [optional] integer array containing the size of each block
3694b9ad928SBarry Smith
3704b9ad928SBarry Smith Options Database Key:
3714b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks
3724b9ad928SBarry Smith
373ce78bad3SBarry Smith Level: intermediate
374ce78bad3SBarry Smith
375f1580f4eSBarry Smith Note:
3764b9ad928SBarry Smith Currently only a limited number of blocking configurations are supported.
377f1580f4eSBarry Smith All processors sharing the `PC` must call this routine with the same data.
3784b9ad928SBarry Smith
379562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()`
3804b9ad928SBarry Smith @*/
PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])381d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
382d71ae5a4SJacob Faibussowitsch {
3834b9ad928SBarry Smith PetscFunctionBegin;
3840700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
38508401ef6SPierre Jolivet PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks");
386cac4c232SBarry Smith PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3884b9ad928SBarry Smith }
3894b9ad928SBarry Smith
3904b9ad928SBarry Smith /*@C
3914b9ad928SBarry Smith PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block
392f1580f4eSBarry Smith Jacobi, `PCBJACOBI`, preconditioner.
3934b9ad928SBarry Smith
394ad4df100SBarry Smith Not Collective
3954b9ad928SBarry Smith
3964b9ad928SBarry Smith Input Parameter:
3974b9ad928SBarry Smith . pc - the preconditioner context
3984b9ad928SBarry Smith
399feefa0e1SJacob Faibussowitsch Output Parameters:
4004b9ad928SBarry Smith + blocks - the number of blocks
4014b9ad928SBarry Smith - lens - integer array containing the size of each block
4024b9ad928SBarry Smith
4034b9ad928SBarry Smith Level: intermediate
4044b9ad928SBarry Smith
405562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()`
4064b9ad928SBarry Smith @*/
PCBJacobiGetTotalBlocks(PC pc,PetscInt * blocks,const PetscInt * lens[])407d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
408d71ae5a4SJacob Faibussowitsch {
4094b9ad928SBarry Smith PetscFunctionBegin;
4100700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4114f572ea9SToby Isaac PetscAssertPointer(blocks, 2);
412cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4144b9ad928SBarry Smith }
4154b9ad928SBarry Smith
4164b9ad928SBarry Smith /*@
4174b9ad928SBarry Smith PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block
418f1580f4eSBarry Smith Jacobi, `PCBJACOBI`, preconditioner.
4194b9ad928SBarry Smith
4204b9ad928SBarry Smith Not Collective
4214b9ad928SBarry Smith
4224b9ad928SBarry Smith Input Parameters:
4234b9ad928SBarry Smith + pc - the preconditioner context
4244b9ad928SBarry Smith . blocks - the number of blocks
4254b9ad928SBarry Smith - lens - [optional] integer array containing size of each block
4264b9ad928SBarry Smith
427342c94f9SMatthew G. Knepley Options Database Key:
428342c94f9SMatthew G. Knepley . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks
429342c94f9SMatthew G. Knepley
430ce78bad3SBarry Smith Level: intermediate
431ce78bad3SBarry Smith
4324b9ad928SBarry Smith Note:
4334b9ad928SBarry Smith Currently only a limited number of blocking configurations are supported.
4344b9ad928SBarry Smith
435562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()`
4364b9ad928SBarry Smith @*/
PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[])437d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[])
438d71ae5a4SJacob Faibussowitsch {
4394b9ad928SBarry Smith PetscFunctionBegin;
4400700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
44108401ef6SPierre Jolivet PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks");
442cac4c232SBarry Smith PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens));
4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4444b9ad928SBarry Smith }
4454b9ad928SBarry Smith
4464b9ad928SBarry Smith /*@C
4474b9ad928SBarry Smith PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block
448f1580f4eSBarry Smith Jacobi, `PCBJACOBI`, preconditioner.
4494b9ad928SBarry Smith
4504b9ad928SBarry Smith Not Collective
4514b9ad928SBarry Smith
4524b9ad928SBarry Smith Input Parameters:
4534b9ad928SBarry Smith + pc - the preconditioner context
4544b9ad928SBarry Smith . blocks - the number of blocks
4554b9ad928SBarry Smith - lens - [optional] integer array containing size of each block
4564b9ad928SBarry Smith
457ce78bad3SBarry Smith Level: intermediate
458ce78bad3SBarry Smith
4594b9ad928SBarry Smith Note:
4604b9ad928SBarry Smith Currently only a limited number of blocking configurations are supported.
4614b9ad928SBarry Smith
462562efe2eSBarry Smith .seealso: [](ch_ksp), `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()`
4634b9ad928SBarry Smith @*/
PCBJacobiGetLocalBlocks(PC pc,PetscInt * blocks,const PetscInt * lens[])464d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[])
465d71ae5a4SJacob Faibussowitsch {
4664b9ad928SBarry Smith PetscFunctionBegin;
4670700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4684f572ea9SToby Isaac PetscAssertPointer(blocks, 2);
469cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens));
4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4714b9ad928SBarry Smith }
4724b9ad928SBarry Smith
4734b9ad928SBarry Smith /*MC
4744b9ad928SBarry Smith PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with
475f1580f4eSBarry 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
481ce78bad3SBarry Smith Level: beginner
482ce78bad3SBarry Smith
48395452b02SPatrick Sanan Notes:
484f1580f4eSBarry Smith See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block
485468ae2e8SBarry Smith
48695452b02SPatrick 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.
4874b9ad928SBarry Smith
488f1580f4eSBarry Smith To set options on the solvers for each block append -sub_ to all the `KSP` and `PC`
489d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
4904b9ad928SBarry Smith
491f1580f4eSBarry Smith To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()`
492f1580f4eSBarry Smith and set the options directly on the resulting `KSP` object (you can access its `PC`
4930462cc06SPierre Jolivet `KSPGetPC()`)
4944b9ad928SBarry Smith
495f1580f4eSBarry Smith For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best
4962210c163SDominic Meiser performance. Different block partitioning may lead to additional data transfers
4972210c163SDominic Meiser between host and GPU that lead to degraded performance.
4982210c163SDominic Meiser
499011ea8aeSBarry Smith When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes.
500011ea8aeSBarry Smith
501562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`,
502db781477SPatrick Sanan `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`,
503db781477SPatrick Sanan `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI`
5044b9ad928SBarry Smith M*/
5054b9ad928SBarry Smith
PCCreate_BJacobi(PC pc)506d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc)
507d71ae5a4SJacob Faibussowitsch {
50869a612a9SBarry Smith PetscMPIInt rank;
5094b9ad928SBarry Smith PC_BJacobi *jac;
5104b9ad928SBarry Smith
5114b9ad928SBarry Smith PetscFunctionBegin;
5124dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac));
5139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
5142fa5cd67SKarl Rupp
5150a545947SLisandro Dalcin pc->ops->apply = NULL;
5167b6e2003SPierre Jolivet pc->ops->matapply = NULL;
5170a545947SLisandro Dalcin pc->ops->applytranspose = NULL;
5184dbf25a8SPierre Jolivet pc->ops->matapplytranspose = NULL;
5194b9ad928SBarry Smith pc->ops->setup = PCSetUp_BJacobi;
5204b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi;
5214b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_BJacobi;
5224b9ad928SBarry Smith pc->ops->view = PCView_BJacobi;
5230a545947SLisandro Dalcin pc->ops->applyrichardson = NULL;
5244b9ad928SBarry Smith
5254b9ad928SBarry Smith pc->data = (void *)jac;
5264b9ad928SBarry Smith jac->n = -1;
5274b9ad928SBarry Smith jac->n_local = -1;
5284b9ad928SBarry Smith jac->first_local = rank;
5290a545947SLisandro Dalcin jac->ksp = NULL;
5300a545947SLisandro Dalcin jac->g_lens = NULL;
5310a545947SLisandro Dalcin jac->l_lens = NULL;
5320a545947SLisandro Dalcin jac->psubcomm = NULL;
5334b9ad928SBarry Smith
5349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi));
5359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi));
5369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi));
5379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi));
5389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi));
5393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5404b9ad928SBarry Smith }
5414b9ad928SBarry Smith
5424b9ad928SBarry Smith /*
5434b9ad928SBarry Smith These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI
5444b9ad928SBarry Smith */
PCReset_BJacobi_Singleblock(PC pc)545d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc)
546d71ae5a4SJacob Faibussowitsch {
5474b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
5484b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5494b9ad928SBarry Smith
5504b9ad928SBarry Smith PetscFunctionBegin;
5519566063dSJacob Faibussowitsch PetscCall(KSPReset(jac->ksp[0]));
5529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->x));
5539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->y));
5543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
555e91c6855SBarry Smith }
556e91c6855SBarry Smith
PCDestroy_BJacobi_Singleblock(PC pc)557d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc)
558d71ae5a4SJacob Faibussowitsch {
559e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
560e91c6855SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
561e91c6855SBarry Smith
562e91c6855SBarry Smith PetscFunctionBegin;
5639566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Singleblock(pc));
5649566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->ksp[0]));
5659566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp));
5669566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac));
5672e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc));
5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5694b9ad928SBarry Smith }
5704b9ad928SBarry Smith
PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)571d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc)
572d71ae5a4SJacob Faibussowitsch {
5734b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
5742295b7c8SHong Zhang KSP subksp = jac->ksp[0];
575539c167fSBarry Smith KSPConvergedReason reason;
5764b9ad928SBarry Smith
5774b9ad928SBarry Smith PetscFunctionBegin;
5789566063dSJacob Faibussowitsch PetscCall(KSPSetUp(subksp));
5799566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(subksp, &reason));
580ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
5813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5824b9ad928SBarry Smith }
5834b9ad928SBarry Smith
PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y)584d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y)
585d71ae5a4SJacob Faibussowitsch {
5864b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
5874b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
5886e111a19SKarl Rupp
5894b9ad928SBarry Smith PetscFunctionBegin;
5909566063dSJacob Faibussowitsch PetscCall(VecGetLocalVectorRead(x, bjac->x));
5919566063dSJacob Faibussowitsch PetscCall(VecGetLocalVector(y, bjac->y));
592bba28a21SBarry Smith /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
593910cf402Sprj- matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
594910cf402Sprj- of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
5959566063dSJacob Faibussowitsch PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
596a9db3dc1SPierre Jolivet PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], bjac->x, bjac->y, 0));
5979566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y));
5989566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
599a9db3dc1SPierre Jolivet PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], bjac->x, bjac->y, 0));
6009566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVectorRead(x, bjac->x));
6019566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVector(y, bjac->y));
6023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6034b9ad928SBarry Smith }
6044b9ad928SBarry Smith
PCMatApply_BJacobi_Singleblock_Private(PC pc,Mat X,Mat Y,PetscBool transpose)6054dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_BJacobi_Singleblock_Private(PC pc, Mat X, Mat Y, PetscBool transpose)
606d71ae5a4SJacob Faibussowitsch {
6077b6e2003SPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6087b6e2003SPierre Jolivet Mat sX, sY;
6097b6e2003SPierre Jolivet
6107b6e2003SPierre Jolivet PetscFunctionBegin;
6117b6e2003SPierre Jolivet /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner
6127b6e2003SPierre Jolivet matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild
6137b6e2003SPierre Jolivet of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/
6149566063dSJacob Faibussowitsch PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner));
6159566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(X, &sX));
6169566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Y, &sY));
61778d05565SPierre Jolivet if (!transpose) {
61878d05565SPierre Jolivet PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], sX, sY, 0));
61978d05565SPierre Jolivet PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
62078d05565SPierre Jolivet PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], sX, sY, 0));
62178d05565SPierre Jolivet } else {
62278d05565SPierre Jolivet PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[0], sX, sY, 0));
62378d05565SPierre Jolivet PetscCall(KSPMatSolveTranspose(jac->ksp[0], sX, sY));
62478d05565SPierre Jolivet PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[0], sX, sY, 0));
62578d05565SPierre Jolivet }
6264dbf25a8SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS);
6274dbf25a8SPierre Jolivet }
6284dbf25a8SPierre Jolivet
PCMatApply_BJacobi_Singleblock(PC pc,Mat X,Mat Y)6294dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
6304dbf25a8SPierre Jolivet {
6314dbf25a8SPierre Jolivet PetscFunctionBegin;
6324dbf25a8SPierre Jolivet PetscCall(PCMatApply_BJacobi_Singleblock_Private(pc, X, Y, PETSC_FALSE));
6334dbf25a8SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS);
6344dbf25a8SPierre Jolivet }
6354dbf25a8SPierre Jolivet
PCMatApplyTranspose_BJacobi_Singleblock(PC pc,Mat X,Mat Y)6364dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_BJacobi_Singleblock(PC pc, Mat X, Mat Y)
6374dbf25a8SPierre Jolivet {
6384dbf25a8SPierre Jolivet PetscFunctionBegin;
6394dbf25a8SPierre Jolivet PetscCall(PCMatApply_BJacobi_Singleblock_Private(pc, X, Y, PETSC_TRUE));
6403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6417b6e2003SPierre Jolivet }
6427b6e2003SPierre Jolivet
PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y)643d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y)
644d71ae5a4SJacob Faibussowitsch {
6454b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6464b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
647d9ca1df4SBarry Smith PetscScalar *y_array;
648d9ca1df4SBarry Smith const PetscScalar *x_array;
6494b9ad928SBarry Smith PC subpc;
6504b9ad928SBarry Smith
6514b9ad928SBarry Smith PetscFunctionBegin;
6524b9ad928SBarry Smith /*
6534b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the
6544b9ad928SBarry Smith y vector into the bjac->x vector. The reason for
6554b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector
6564b9ad928SBarry Smith for the sequential solve.
6574b9ad928SBarry Smith */
6589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array));
6599566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array));
6609566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array));
6619566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array));
6624b9ad928SBarry Smith /* apply the symmetric left portion of the inner PC operator */
663c3f9dca2SPierre Jolivet /* note this bypasses the inner KSP and its options completely */
6649566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6659566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y));
6669566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x));
6679566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y));
6689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array));
6699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array));
6703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6714b9ad928SBarry Smith }
6724b9ad928SBarry Smith
PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y)673d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y)
674d71ae5a4SJacob Faibussowitsch {
6754b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
6764b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
677d9ca1df4SBarry Smith PetscScalar *y_array;
678d9ca1df4SBarry Smith const PetscScalar *x_array;
6794b9ad928SBarry Smith PC subpc;
6804b9ad928SBarry Smith
6814b9ad928SBarry Smith PetscFunctionBegin;
6824b9ad928SBarry Smith /*
6834b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the
6844b9ad928SBarry Smith y vector into the bjac->x vector. The reason for
6854b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector
6864b9ad928SBarry Smith for the sequential solve.
6874b9ad928SBarry Smith */
6889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array));
6899566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array));
6909566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array));
6919566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array));
6924b9ad928SBarry Smith
6934b9ad928SBarry Smith /* apply the symmetric right portion of the inner PC operator */
694c3f9dca2SPierre Jolivet /* note this bypasses the inner KSP and its options completely */
6954b9ad928SBarry Smith
6969566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &subpc));
6979566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y));
6984b9ad928SBarry Smith
69911e4f71cSBarry Smith PetscCall(VecResetArray(bjac->x));
70011e4f71cSBarry Smith PetscCall(VecResetArray(bjac->y));
7019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array));
7029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array));
7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7044b9ad928SBarry Smith }
7054b9ad928SBarry Smith
PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y)706d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y)
707d71ae5a4SJacob Faibussowitsch {
7084b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
7094b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data;
710d9ca1df4SBarry Smith PetscScalar *y_array;
711d9ca1df4SBarry Smith const PetscScalar *x_array;
7124b9ad928SBarry Smith
7134b9ad928SBarry Smith PetscFunctionBegin;
7144b9ad928SBarry Smith /*
7154b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the
7164b9ad928SBarry Smith y vector into the bjac->x vector. The reason for
7174b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector
7184b9ad928SBarry Smith for the sequential solve.
7194b9ad928SBarry Smith */
7209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array));
7219566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array));
7229566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array));
7239566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array));
724a9db3dc1SPierre Jolivet PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[0], bjac->x, bjac->y, 0));
7259566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y));
7269566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y));
727a9db3dc1SPierre Jolivet PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[0], bjac->x, bjac->y, 0));
7289566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x));
7299566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y));
7309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array));
7319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array));
7323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7334b9ad928SBarry Smith }
7344b9ad928SBarry Smith
PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat)735d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat)
736d71ae5a4SJacob Faibussowitsch {
7374b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
73869a612a9SBarry Smith PetscInt m;
7394b9ad928SBarry Smith KSP ksp;
7404b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac;
741de0953f6SBarry Smith PetscBool wasSetup = PETSC_TRUE;
7423f6dc190SJunchao Zhang VecType vectype;
743ea41da7aSPierre Jolivet const char *prefix;
7444b9ad928SBarry Smith
7454b9ad928SBarry Smith PetscFunctionBegin;
7464b9ad928SBarry Smith if (!pc->setupcalled) {
747c2efdce3SBarry Smith if (!jac->ksp) {
7483821be0aSBarry Smith PetscInt nestlevel;
7493821be0aSBarry Smith
750302a9ddcSMatthew Knepley wasSetup = PETSC_FALSE;
7512fa5cd67SKarl Rupp
7529566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
7533821be0aSBarry Smith PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
7543821be0aSBarry Smith PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
7559566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
7569566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
7579566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY));
7589566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix));
7599566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix));
7609566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
7614b9ad928SBarry Smith
762e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Singleblock;
7634b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi_Singleblock;
7644b9ad928SBarry Smith pc->ops->apply = PCApply_BJacobi_Singleblock;
7657b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_BJacobi_Singleblock;
7664dbf25a8SPierre Jolivet pc->ops->matapplytranspose = PCMatApplyTranspose_BJacobi_Singleblock;
7674b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock;
7684b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock;
7694b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock;
7704b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock;
7714b9ad928SBarry Smith
7729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->ksp));
7734b9ad928SBarry Smith jac->ksp[0] = ksp;
774c2efdce3SBarry Smith
7754dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bjac));
7764b9ad928SBarry Smith jac->data = (void *)bjac;
7774b9ad928SBarry Smith } else {
778c2efdce3SBarry Smith ksp = jac->ksp[0];
779c2efdce3SBarry Smith bjac = (PC_BJacobi_Singleblock *)jac->data;
780c2efdce3SBarry Smith }
781c2efdce3SBarry Smith
782c2efdce3SBarry Smith /*
783c2efdce3SBarry Smith The reason we need to generate these vectors is to serve
784c2efdce3SBarry Smith as the right-hand side and solution vector for the solve on the
785c2efdce3SBarry Smith block. We do not need to allocate space for the vectors since
786c2efdce3SBarry Smith that is provided via VecPlaceArray() just before the call to
787c2efdce3SBarry Smith KSPSolve() on the block.
788c2efdce3SBarry Smith */
7899566063dSJacob Faibussowitsch PetscCall(MatGetSize(pmat, &m, &m));
7909566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x));
7919566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y));
7929566063dSJacob Faibussowitsch PetscCall(MatGetVecType(pmat, &vectype));
7939566063dSJacob Faibussowitsch PetscCall(VecSetType(bjac->x, vectype));
7949566063dSJacob Faibussowitsch PetscCall(VecSetType(bjac->y, vectype));
795c2efdce3SBarry Smith } else {
7964b9ad928SBarry Smith ksp = jac->ksp[0];
7974b9ad928SBarry Smith bjac = (PC_BJacobi_Singleblock *)jac->data;
7984b9ad928SBarry Smith }
7999566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(ksp, &prefix));
80049517cdeSBarry Smith if (pc->useAmat) {
8019566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, mat, pmat));
8029566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(mat, prefix));
8034b9ad928SBarry Smith } else {
8049566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, pmat, pmat));
8054b9ad928SBarry Smith }
8069566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(pmat, prefix));
80790f1c854SHong Zhang if (!wasSetup && pc->setfromoptionscalled) {
808248bfaf0SJed Brown /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */
8099566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp));
810302a9ddcSMatthew Knepley }
8113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8124b9ad928SBarry Smith }
8134b9ad928SBarry Smith
PCReset_BJacobi_Multiblock(PC pc)814d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc)
815d71ae5a4SJacob Faibussowitsch {
8164b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
8174b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
81869a612a9SBarry Smith PetscInt i;
8194b9ad928SBarry Smith
8204b9ad928SBarry Smith PetscFunctionBegin;
821e91c6855SBarry Smith if (bjac && bjac->pmat) {
8229566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat));
82348a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat));
824e91c6855SBarry Smith }
8254b9ad928SBarry Smith
8264b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) {
8279566063dSJacob Faibussowitsch PetscCall(KSPReset(jac->ksp[i]));
828e91c6855SBarry Smith if (bjac && bjac->x) {
8299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->x[i]));
8309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->y[i]));
8319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bjac->is[i]));
8324b9ad928SBarry Smith }
833e91c6855SBarry Smith }
8349566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->l_lens));
8359566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->g_lens));
8363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
837e91c6855SBarry Smith }
838e91c6855SBarry Smith
PCDestroy_BJacobi_Multiblock(PC pc)839d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc)
840d71ae5a4SJacob Faibussowitsch {
841e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
842c2efdce3SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
843e91c6855SBarry Smith PetscInt i;
844e91c6855SBarry Smith
845e91c6855SBarry Smith PetscFunctionBegin;
8469566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Multiblock(pc));
847c2efdce3SBarry Smith if (bjac) {
8489566063dSJacob Faibussowitsch PetscCall(PetscFree2(bjac->x, bjac->y));
8499566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac->starts));
8509566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac->is));
851c2efdce3SBarry Smith }
8529566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->data));
85348a46eb9SPierre Jolivet for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i]));
8549566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp));
8552e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc));
8563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8574b9ad928SBarry Smith }
8584b9ad928SBarry Smith
PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)859d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc)
860d71ae5a4SJacob Faibussowitsch {
8614b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
86269a612a9SBarry Smith PetscInt i, n_local = jac->n_local;
863539c167fSBarry Smith KSPConvergedReason reason;
8644b9ad928SBarry Smith
8654b9ad928SBarry Smith PetscFunctionBegin;
8664b9ad928SBarry Smith for (i = 0; i < n_local; i++) {
8679566063dSJacob Faibussowitsch PetscCall(KSPSetUp(jac->ksp[i]));
8689566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason));
869ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
8704b9ad928SBarry Smith }
8713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8724b9ad928SBarry Smith }
8734b9ad928SBarry Smith
PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y)874d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y)
875d71ae5a4SJacob Faibussowitsch {
8764b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
87769a612a9SBarry Smith PetscInt i, n_local = jac->n_local;
8784b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
879d9ca1df4SBarry Smith PetscScalar *yin;
880d9ca1df4SBarry Smith const PetscScalar *xin;
88158ebbce7SBarry Smith
8824b9ad928SBarry Smith PetscFunctionBegin;
8839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xin));
8849566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yin));
8854b9ad928SBarry Smith for (i = 0; i < n_local; i++) {
8864b9ad928SBarry Smith /*
8874b9ad928SBarry Smith To avoid copying the subvector from x into a workspace we instead
8884b9ad928SBarry Smith make the workspace vector array point to the subpart of the array of
8894b9ad928SBarry Smith the global vector.
8904b9ad928SBarry Smith */
8919566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
8929566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
8934b9ad928SBarry Smith
8949566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
8959566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i]));
8969566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
8979566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
898d11f3a42SBarry Smith
8999566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x[i]));
9009566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y[i]));
9014b9ad928SBarry Smith }
9029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xin));
9039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yin));
9043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9054b9ad928SBarry Smith }
9064b9ad928SBarry Smith
PCApplySymmetricLeft_BJacobi_Multiblock(PC pc,Vec x,Vec y)907f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricLeft_BJacobi_Multiblock(PC pc, Vec x, Vec y)
908f4d694ddSBarry Smith {
909f4d694ddSBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
910f4d694ddSBarry Smith PetscInt i, n_local = jac->n_local;
911f4d694ddSBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
912f4d694ddSBarry Smith PetscScalar *yin;
913f4d694ddSBarry Smith const PetscScalar *xin;
914f4d694ddSBarry Smith PC subpc;
915f4d694ddSBarry Smith
916f4d694ddSBarry Smith PetscFunctionBegin;
917f4d694ddSBarry Smith PetscCall(VecGetArrayRead(x, &xin));
918f4d694ddSBarry Smith PetscCall(VecGetArray(y, &yin));
919f4d694ddSBarry Smith for (i = 0; i < n_local; i++) {
9204b9ad928SBarry Smith /*
921f4d694ddSBarry Smith To avoid copying the subvector from x into a workspace we instead
922f4d694ddSBarry Smith make the workspace vector array point to the subpart of the array of
923f4d694ddSBarry Smith the global vector.
9244b9ad928SBarry Smith */
925f4d694ddSBarry Smith PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
926f4d694ddSBarry Smith PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
927f4d694ddSBarry Smith
928f4d694ddSBarry Smith PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
929f4d694ddSBarry Smith /* apply the symmetric left portion of the inner PC operator */
930c3f9dca2SPierre Jolivet /* note this bypasses the inner KSP and its options completely */
931f4d694ddSBarry Smith PetscCall(KSPGetPC(jac->ksp[i], &subpc));
932f4d694ddSBarry Smith PetscCall(PCApplySymmetricLeft(subpc, bjac->x[i], bjac->y[i]));
933f4d694ddSBarry Smith PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
934f4d694ddSBarry Smith
935f4d694ddSBarry Smith PetscCall(VecResetArray(bjac->x[i]));
936f4d694ddSBarry Smith PetscCall(VecResetArray(bjac->y[i]));
937f4d694ddSBarry Smith }
938f4d694ddSBarry Smith PetscCall(VecRestoreArrayRead(x, &xin));
939f4d694ddSBarry Smith PetscCall(VecRestoreArray(y, &yin));
9403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
941f4d694ddSBarry Smith }
942f4d694ddSBarry Smith
PCApplySymmetricRight_BJacobi_Multiblock(PC pc,Vec x,Vec y)943f4d694ddSBarry Smith static PetscErrorCode PCApplySymmetricRight_BJacobi_Multiblock(PC pc, Vec x, Vec y)
944f4d694ddSBarry Smith {
945f4d694ddSBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
946f4d694ddSBarry Smith PetscInt i, n_local = jac->n_local;
947f4d694ddSBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
948f4d694ddSBarry Smith PetscScalar *yin;
949f4d694ddSBarry Smith const PetscScalar *xin;
950f4d694ddSBarry Smith PC subpc;
951f4d694ddSBarry Smith
952f4d694ddSBarry Smith PetscFunctionBegin;
953f4d694ddSBarry Smith PetscCall(VecGetArrayRead(x, &xin));
954f4d694ddSBarry Smith PetscCall(VecGetArray(y, &yin));
955f4d694ddSBarry Smith for (i = 0; i < n_local; i++) {
956f4d694ddSBarry Smith /*
957f4d694ddSBarry Smith To avoid copying the subvector from x into a workspace we instead
958f4d694ddSBarry Smith make the workspace vector array point to the subpart of the array of
959f4d694ddSBarry Smith the global vector.
960f4d694ddSBarry Smith */
961f4d694ddSBarry Smith PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
962f4d694ddSBarry Smith PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
963f4d694ddSBarry Smith
964f4d694ddSBarry Smith PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
965f4d694ddSBarry Smith /* apply the symmetric left portion of the inner PC operator */
966c3f9dca2SPierre Jolivet /* note this bypasses the inner KSP and its options completely */
967f4d694ddSBarry Smith PetscCall(KSPGetPC(jac->ksp[i], &subpc));
968f4d694ddSBarry Smith PetscCall(PCApplySymmetricRight(subpc, bjac->x[i], bjac->y[i]));
969f4d694ddSBarry Smith PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
970f4d694ddSBarry Smith
971f4d694ddSBarry Smith PetscCall(VecResetArray(bjac->x[i]));
972f4d694ddSBarry Smith PetscCall(VecResetArray(bjac->y[i]));
973f4d694ddSBarry Smith }
974f4d694ddSBarry Smith PetscCall(VecRestoreArrayRead(x, &xin));
975f4d694ddSBarry Smith PetscCall(VecRestoreArray(y, &yin));
9763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
977f4d694ddSBarry Smith }
978f4d694ddSBarry Smith
PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y)979d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y)
980d71ae5a4SJacob Faibussowitsch {
9814b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
98269a612a9SBarry Smith PetscInt i, n_local = jac->n_local;
9834b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
984d9ca1df4SBarry Smith PetscScalar *yin;
985d9ca1df4SBarry Smith const PetscScalar *xin;
9864b9ad928SBarry Smith
9874b9ad928SBarry Smith PetscFunctionBegin;
9889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xin));
9899566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yin));
9904b9ad928SBarry Smith for (i = 0; i < n_local; i++) {
9914b9ad928SBarry Smith /*
9924b9ad928SBarry Smith To avoid copying the subvector from x into a workspace we instead
9934b9ad928SBarry Smith make the workspace vector array point to the subpart of the array of
9944b9ad928SBarry Smith the global vector.
9954b9ad928SBarry Smith */
9969566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i]));
9979566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i]));
9984b9ad928SBarry Smith
9999566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
10009566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i]));
10019566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i]));
10029566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0));
1003a7ff49e8SJed Brown
10049566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x[i]));
10059566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y[i]));
10064b9ad928SBarry Smith }
10079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xin));
10089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yin));
10093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10104b9ad928SBarry Smith }
10114b9ad928SBarry Smith
PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat)1012d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat)
1013d71ae5a4SJacob Faibussowitsch {
10144b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
101569a612a9SBarry Smith PetscInt m, n_local, N, M, start, i;
1016ea41da7aSPierre Jolivet const char *prefix;
10174b9ad928SBarry Smith KSP ksp;
10184b9ad928SBarry Smith Vec x, y;
10194b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data;
10204b9ad928SBarry Smith PC subpc;
10214b9ad928SBarry Smith IS is;
1022434a97beSBrad Aagaard MatReuse scall;
10233f6dc190SJunchao Zhang VecType vectype;
10244849c82aSBarry Smith MatNullSpace *nullsp_mat = NULL, *nullsp_pmat = NULL;
10254b9ad928SBarry Smith
10264b9ad928SBarry Smith PetscFunctionBegin;
10279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &M, &N));
10284b9ad928SBarry Smith
10294b9ad928SBarry Smith n_local = jac->n_local;
10304b9ad928SBarry Smith
103149517cdeSBarry Smith if (pc->useAmat) {
1032ace3abfcSBarry Smith PetscBool same;
10339566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same));
103428b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type");
10354b9ad928SBarry Smith }
10364b9ad928SBarry Smith
10374b9ad928SBarry Smith if (!pc->setupcalled) {
10383821be0aSBarry Smith PetscInt nestlevel;
10393821be0aSBarry Smith
10404b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX;
1041c2efdce3SBarry Smith
1042c2efdce3SBarry Smith if (!jac->ksp) {
1043e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Multiblock;
10444b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi_Multiblock;
10454b9ad928SBarry Smith pc->ops->apply = PCApply_BJacobi_Multiblock;
10467b6e2003SPierre Jolivet pc->ops->matapply = NULL;
10474dbf25a8SPierre Jolivet pc->ops->matapplytranspose = NULL;
1048f4d694ddSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Multiblock;
1049f4d694ddSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Multiblock;
10504b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_BJacobi_Multiblock;
10514b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock;
10524b9ad928SBarry Smith
10534dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bjac));
10549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &jac->ksp));
10559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y));
10569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &bjac->starts));
10574b9ad928SBarry Smith
10584b9ad928SBarry Smith jac->data = (void *)bjac;
10599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &bjac->is));
10604b9ad928SBarry Smith
10614b9ad928SBarry Smith for (i = 0; i < n_local; i++) {
10629566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
10633821be0aSBarry Smith PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
10643821be0aSBarry Smith PetscCall(KSPSetNestLevel(ksp, nestlevel + 1));
10659566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
10669566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
10679566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY));
10689566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &subpc));
10699566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix));
10709566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix));
10719566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
10722fa5cd67SKarl Rupp
1073c2efdce3SBarry Smith jac->ksp[i] = ksp;
1074c2efdce3SBarry Smith }
1075c2efdce3SBarry Smith } else {
1076c2efdce3SBarry Smith bjac = (PC_BJacobi_Multiblock *)jac->data;
1077c2efdce3SBarry Smith }
10784b9ad928SBarry Smith
1079c2efdce3SBarry Smith start = 0;
10809566063dSJacob Faibussowitsch PetscCall(MatGetVecType(pmat, &vectype));
1081c2efdce3SBarry Smith for (i = 0; i < n_local; i++) {
10824b9ad928SBarry Smith m = jac->l_lens[i];
10834b9ad928SBarry Smith /*
10844b9ad928SBarry Smith The reason we need to generate these vectors is to serve
10854b9ad928SBarry Smith as the right-hand side and solution vector for the solve on the
10864b9ad928SBarry Smith block. We do not need to allocate space for the vectors since
10874b9ad928SBarry Smith that is provided via VecPlaceArray() just before the call to
10884b9ad928SBarry Smith KSPSolve() on the block.
10894b9ad928SBarry Smith
10904b9ad928SBarry Smith */
10919566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x));
10929566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y));
10939566063dSJacob Faibussowitsch PetscCall(VecSetType(x, vectype));
10949566063dSJacob Faibussowitsch PetscCall(VecSetType(y, vectype));
10952fa5cd67SKarl Rupp
10964b9ad928SBarry Smith bjac->x[i] = x;
10974b9ad928SBarry Smith bjac->y[i] = y;
10984b9ad928SBarry Smith bjac->starts[i] = start;
10994b9ad928SBarry Smith
11009566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is));
11014b9ad928SBarry Smith bjac->is[i] = is;
11024b9ad928SBarry Smith
11034b9ad928SBarry Smith start += m;
11044b9ad928SBarry Smith }
11054b9ad928SBarry Smith } else {
11064b9ad928SBarry Smith bjac = (PC_BJacobi_Multiblock *)jac->data;
11074b9ad928SBarry Smith /*
11084b9ad928SBarry Smith Destroy the blocks from the previous iteration
11094b9ad928SBarry Smith */
11104b9ad928SBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
11114849c82aSBarry Smith PetscCall(MatGetNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
11129566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(n_local, &bjac->pmat));
11134849c82aSBarry Smith if (pc->useAmat) {
11144849c82aSBarry Smith PetscCall(MatGetNullSpaces(n_local, bjac->mat, &nullsp_mat));
11154849c82aSBarry Smith PetscCall(MatDestroyMatrices(n_local, &bjac->mat));
11164849c82aSBarry Smith }
11174b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX;
11182fa5cd67SKarl Rupp } else scall = MAT_REUSE_MATRIX;
11194b9ad928SBarry Smith }
11204b9ad928SBarry Smith
11219566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat));
11224849c82aSBarry Smith if (nullsp_pmat) PetscCall(MatRestoreNullSpaces(n_local, bjac->pmat, &nullsp_pmat));
11234849c82aSBarry Smith if (pc->useAmat) {
11244849c82aSBarry Smith PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat));
11254849c82aSBarry Smith if (nullsp_mat) PetscCall(MatRestoreNullSpaces(n_local, bjac->mat, &nullsp_mat));
11264849c82aSBarry Smith }
11274b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply
11284b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */
11299566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP));
11304b9ad928SBarry Smith
11314b9ad928SBarry Smith for (i = 0; i < n_local; i++) {
11329566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix));
113349517cdeSBarry Smith if (pc->useAmat) {
11349566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i]));
11359566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix));
11364b9ad928SBarry Smith } else {
11379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i]));
11384b9ad928SBarry Smith }
11399566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix));
114048a46eb9SPierre Jolivet if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i]));
114190f1c854SHong Zhang }
11423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11434b9ad928SBarry Smith }
11445a7e1895SHong Zhang
11455a7e1895SHong Zhang /*
1146fd0b8932SBarry Smith These are for a single block with multiple processes
11475a7e1895SHong Zhang */
PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)1148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc)
1149d71ae5a4SJacob Faibussowitsch {
1150fd0b8932SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1151fd0b8932SBarry Smith KSP subksp = jac->ksp[0];
1152fd0b8932SBarry Smith KSPConvergedReason reason;
1153fd0b8932SBarry Smith
1154fd0b8932SBarry Smith PetscFunctionBegin;
11559566063dSJacob Faibussowitsch PetscCall(KSPSetUp(subksp));
11569566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(subksp, &reason));
1157ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1159fd0b8932SBarry Smith }
1160fd0b8932SBarry Smith
PCReset_BJacobi_Multiproc(PC pc)1161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc)
1162d71ae5a4SJacob Faibussowitsch {
11635a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data;
11645a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
11655a7e1895SHong Zhang
11665a7e1895SHong Zhang PetscFunctionBegin;
11679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mpjac->ysub));
11689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mpjac->xsub));
11699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mpjac->submats));
11709566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPReset(jac->ksp[0]));
11713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1172e91c6855SBarry Smith }
1173e91c6855SBarry Smith
PCDestroy_BJacobi_Multiproc(PC pc)1174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc)
1175d71ae5a4SJacob Faibussowitsch {
1176e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data;
1177e91c6855SBarry Smith PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1178e91c6855SBarry Smith
1179e91c6855SBarry Smith PetscFunctionBegin;
11809566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Multiproc(pc));
11819566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->ksp[0]));
11829566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp));
11839566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&mpjac->psubcomm));
11845a7e1895SHong Zhang
11859566063dSJacob Faibussowitsch PetscCall(PetscFree(mpjac));
11862e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc));
11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11885a7e1895SHong Zhang }
11895a7e1895SHong Zhang
PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y)1190d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y)
1191d71ae5a4SJacob Faibussowitsch {
11925a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data;
11935a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
1194d9ca1df4SBarry Smith PetscScalar *yarray;
1195d9ca1df4SBarry Smith const PetscScalar *xarray;
1196539c167fSBarry Smith KSPConvergedReason reason;
11975a7e1895SHong Zhang
11985a7e1895SHong Zhang PetscFunctionBegin;
11995a7e1895SHong Zhang /* place x's and y's local arrays into xsub and ysub */
12009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xarray));
12019566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yarray));
12029566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mpjac->xsub, xarray));
12039566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mpjac->ysub, yarray));
12045a7e1895SHong Zhang
12055a7e1895SHong Zhang /* apply preconditioner on each matrix block */
12069566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
12079566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub));
12089566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub));
12099566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0));
12109566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1211ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
1212250cf88bSHong Zhang
12139566063dSJacob Faibussowitsch PetscCall(VecResetArray(mpjac->xsub));
12149566063dSJacob Faibussowitsch PetscCall(VecResetArray(mpjac->ysub));
12159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xarray));
12169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yarray));
12173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12185a7e1895SHong Zhang }
12195a7e1895SHong Zhang
PCMatApply_BJacobi_Multiproc(PC pc,Mat X,Mat Y)1220d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y)
1221d71ae5a4SJacob Faibussowitsch {
12227b6e2003SPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data;
12237b6e2003SPierre Jolivet KSPConvergedReason reason;
12247b6e2003SPierre Jolivet Mat sX, sY;
12257b6e2003SPierre Jolivet const PetscScalar *x;
12267b6e2003SPierre Jolivet PetscScalar *y;
12277b6e2003SPierre Jolivet PetscInt m, N, lda, ldb;
12287b6e2003SPierre Jolivet
12297b6e2003SPierre Jolivet PetscFunctionBegin;
12307b6e2003SPierre Jolivet /* apply preconditioner on each matrix block */
12319566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m, NULL));
12329566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N));
12339566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda));
12349566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(Y, &ldb));
12359566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X, &x));
12369566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(Y, &y));
12379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX));
12389566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY));
12399566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(sX, lda));
12409566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(sY, ldb));
12419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12429566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(jac->ksp[0], sX, sY));
12439566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL));
12449566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0));
12459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sY));
12469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sX));
12479566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(Y, &y));
12489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X, &x));
12499566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason));
1250ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
12513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12527b6e2003SPierre Jolivet }
12537b6e2003SPierre Jolivet
PCSetUp_BJacobi_Multiproc(PC pc)1254d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc)
1255d71ae5a4SJacob Faibussowitsch {
12565a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data;
12575a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data;
12585a7e1895SHong Zhang PetscInt m, n;
1259ce94432eSBarry Smith MPI_Comm comm, subcomm = 0;
12605a7e1895SHong Zhang const char *prefix;
1261de0953f6SBarry Smith PetscBool wasSetup = PETSC_TRUE;
12623f6dc190SJunchao Zhang VecType vectype;
12631f62890fSKarl Rupp
12645a7e1895SHong Zhang PetscFunctionBegin;
12659566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
126608401ef6SPierre Jolivet PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported");
12675a7e1895SHong Zhang jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */
12685a7e1895SHong Zhang if (!pc->setupcalled) {
12693821be0aSBarry Smith PetscInt nestlevel;
12703821be0aSBarry Smith
1271de0953f6SBarry Smith wasSetup = PETSC_FALSE;
12724dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mpjac));
12735a7e1895SHong Zhang jac->data = (void *)mpjac;
12745a7e1895SHong Zhang
12755a7e1895SHong Zhang /* initialize datastructure mpjac */
12765a7e1895SHong Zhang if (!jac->psubcomm) {
12775a7e1895SHong Zhang /* Create default contiguous subcommunicatiors if user does not provide them */
12789566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(comm, &jac->psubcomm));
12799566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n));
12809566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
12815a7e1895SHong Zhang }
12825a7e1895SHong Zhang mpjac->psubcomm = jac->psubcomm;
1283306c2d5bSBarry Smith subcomm = PetscSubcommChild(mpjac->psubcomm);
12845a7e1895SHong Zhang
12855a7e1895SHong Zhang /* Get matrix blocks of pmat */
12869566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
12875a7e1895SHong Zhang
12885a7e1895SHong Zhang /* create a new PC that processors in each subcomm have copy of */
12899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->ksp));
12909566063dSJacob Faibussowitsch PetscCall(KSPCreate(subcomm, &jac->ksp[0]));
12913821be0aSBarry Smith PetscCall(PCGetKSPNestLevel(pc, &nestlevel));
12923821be0aSBarry Smith PetscCall(KSPSetNestLevel(jac->ksp[0], nestlevel + 1));
12939566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure));
12949566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1));
12959566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
12969566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc));
12975a7e1895SHong Zhang
12989566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix));
12999566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix));
13009566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_"));
13019566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix));
13029566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix));
13035a7e1895SHong Zhang
13049566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mpjac->submats, &m, &n));
13059566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub));
13069566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub));
13079566063dSJacob Faibussowitsch PetscCall(MatGetVecType(mpjac->submats, &vectype));
13089566063dSJacob Faibussowitsch PetscCall(VecSetType(mpjac->xsub, vectype));
13099566063dSJacob Faibussowitsch PetscCall(VecSetType(mpjac->ysub, vectype));
13105a7e1895SHong Zhang
1311fd0b8932SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc;
1312e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Multiproc;
13135a7e1895SHong Zhang pc->ops->destroy = PCDestroy_BJacobi_Multiproc;
13145a7e1895SHong Zhang pc->ops->apply = PCApply_BJacobi_Multiproc;
13157b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_BJacobi_Multiproc;
1316fc08c53fSHong Zhang } else { /* pc->setupcalled */
1317306c2d5bSBarry Smith subcomm = PetscSubcommChild(mpjac->psubcomm);
13189e0ae222SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
13195a7e1895SHong Zhang /* destroy old matrix blocks, then get new matrix blocks */
13209566063dSJacob Faibussowitsch if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats));
13219566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats));
1322fc08c53fSHong Zhang } else {
13239566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats));
13245a7e1895SHong Zhang }
13259566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats));
13265a7e1895SHong Zhang }
13275a7e1895SHong Zhang
132848a46eb9SPierre Jolivet if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0]));
13293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13305a7e1895SHong Zhang }
1331