1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a block Jacobi preconditioner. 44b9ad928SBarry Smith */ 500e125f8SBarry Smith 600e125f8SBarry Smith #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/ 74b9ad928SBarry Smith 86849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat); 96849ba73SBarry Smith static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat); 105a7e1895SHong Zhang static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC); 114b9ad928SBarry Smith 129371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi(PC pc) { 134b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 144b9ad928SBarry Smith Mat mat = pc->mat, pmat = pc->pmat; 15976e8c5aSLisandro Dalcin PetscBool hasop; 1669a612a9SBarry Smith PetscInt N, M, start, i, sum, end; 1769a612a9SBarry Smith PetscInt bs, i_start = -1, i_end = -1; 1869a612a9SBarry Smith PetscMPIInt rank, size; 194b9ad928SBarry Smith 204b9ad928SBarry Smith PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 239566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &M, &N)); 249566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(pc->pmat, &bs)); 254b9ad928SBarry Smith 265a7e1895SHong Zhang if (jac->n > 0 && jac->n < size) { 279566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Multiproc(pc)); 285a7e1895SHong Zhang PetscFunctionReturn(0); 295a7e1895SHong Zhang } 305a7e1895SHong Zhang 31f1580f4eSBarry Smith /* Determines the number of blocks assigned to each processor */ 324b9ad928SBarry Smith /* local block count given */ 334b9ad928SBarry Smith if (jac->n_local > 0 && jac->n < 0) { 341c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 354b9ad928SBarry Smith if (jac->l_lens) { /* check that user set these correctly */ 364b9ad928SBarry Smith sum = 0; 374b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) { 3808401ef6SPierre Jolivet PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout"); 394b9ad928SBarry Smith sum += jac->l_lens[i]; 404b9ad928SBarry Smith } 4108401ef6SPierre Jolivet PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly"); 424b9ad928SBarry Smith } else { 439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 442fa5cd67SKarl Rupp for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)); 454b9ad928SBarry Smith } 464b9ad928SBarry Smith } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */ 474b9ad928SBarry Smith /* global blocks given: determine which ones are local */ 484b9ad928SBarry Smith if (jac->g_lens) { 494b9ad928SBarry Smith /* check if the g_lens is has valid entries */ 504b9ad928SBarry Smith for (i = 0; i < jac->n; i++) { 517827d75bSBarry Smith PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed"); 5208401ef6SPierre Jolivet PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout"); 534b9ad928SBarry Smith } 544b9ad928SBarry Smith if (size == 1) { 554b9ad928SBarry Smith jac->n_local = jac->n; 569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 579566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local)); 584b9ad928SBarry Smith /* check that user set these correctly */ 594b9ad928SBarry Smith sum = 0; 604b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i]; 6108401ef6SPierre Jolivet PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly"); 624b9ad928SBarry Smith } else { 639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end)); 644b9ad928SBarry Smith /* loop over blocks determing first one owned by me */ 654b9ad928SBarry Smith sum = 0; 664b9ad928SBarry Smith for (i = 0; i < jac->n + 1; i++) { 679371c9d4SSatish Balay if (sum == start) { 689371c9d4SSatish Balay i_start = i; 699371c9d4SSatish Balay goto start_1; 709371c9d4SSatish Balay } 714b9ad928SBarry Smith if (i < jac->n) sum += jac->g_lens[i]; 724b9ad928SBarry Smith } 73e7e72b3dSBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout"); 744b9ad928SBarry Smith start_1: 754b9ad928SBarry Smith for (i = i_start; i < jac->n + 1; i++) { 769371c9d4SSatish Balay if (sum == end) { 779371c9d4SSatish Balay i_end = i; 789371c9d4SSatish Balay goto end_1; 799371c9d4SSatish Balay } 804b9ad928SBarry Smith if (i < jac->n) sum += jac->g_lens[i]; 814b9ad928SBarry Smith } 82e7e72b3dSBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout"); 834b9ad928SBarry Smith end_1: 844b9ad928SBarry Smith jac->n_local = i_end - i_start; 859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 869566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local)); 874b9ad928SBarry Smith } 884b9ad928SBarry Smith } else { /* no global blocks given, determine then using default layout */ 894b9ad928SBarry Smith jac->n_local = jac->n / size + ((jac->n % size) > rank); 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 914b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) { 924b9ad928SBarry Smith jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs; 937827d75bSBarry Smith PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given"); 944b9ad928SBarry Smith } 954b9ad928SBarry Smith } 964b9ad928SBarry Smith } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */ 974b9ad928SBarry Smith jac->n = size; 984b9ad928SBarry Smith jac->n_local = 1; 999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->l_lens)); 1004b9ad928SBarry Smith jac->l_lens[0] = M; 1017271b979SBarry Smith } else { /* jac->n > 0 && jac->n_local > 0 */ 1027271b979SBarry Smith if (!jac->l_lens) { 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 1042fa5cd67SKarl Rupp for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)); 1057271b979SBarry Smith } 1064b9ad928SBarry Smith } 10708401ef6SPierre Jolivet PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors"); 1084b9ad928SBarry Smith 109f1580f4eSBarry Smith /* Determines mat and pmat */ 1109566063dSJacob Faibussowitsch PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 111976e8c5aSLisandro Dalcin if (!hasop && size == 1) { 1124b9ad928SBarry Smith mat = pc->mat; 1134b9ad928SBarry Smith pmat = pc->pmat; 1144b9ad928SBarry Smith } else { 11549517cdeSBarry Smith if (pc->useAmat) { 11649517cdeSBarry Smith /* use block from Amat matrix, not Pmat for local MatMult() */ 1179566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(pc->mat, &mat)); 1184b9ad928SBarry Smith } 11949517cdeSBarry Smith if (pc->pmat != pc->mat || !pc->useAmat) { 1209566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat)); 1212fa5cd67SKarl Rupp } else pmat = mat; 1224b9ad928SBarry Smith } 1234b9ad928SBarry Smith 124f1580f4eSBarry Smith /* 1254b9ad928SBarry Smith Setup code depends on the number of blocks 1264b9ad928SBarry Smith */ 127cc1d6079SHong Zhang if (jac->n_local == 1) { 1289566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat)); 1294b9ad928SBarry Smith } else { 1309566063dSJacob Faibussowitsch PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat)); 1314b9ad928SBarry Smith } 1324b9ad928SBarry Smith PetscFunctionReturn(0); 1334b9ad928SBarry Smith } 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith /* Default destroy, if it has never been setup */ 1369371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi(PC pc) { 1374b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith PetscFunctionBegin; 1409566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->g_lens)); 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->l_lens)); 1422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL)); 1432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL)); 1442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL)); 1452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL)); 1462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1484b9ad928SBarry Smith PetscFunctionReturn(0); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith 1519371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject) { 1524b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 153248bfaf0SJed Brown PetscInt blocks, i; 154ace3abfcSBarry Smith PetscBool flg; 1554b9ad928SBarry Smith 1564b9ad928SBarry Smith PetscFunctionBegin; 157d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options"); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg)); 1599566063dSJacob Faibussowitsch if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL)); 1609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg)); 1619566063dSJacob Faibussowitsch if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL)); 162248bfaf0SJed Brown if (jac->ksp) { 163248bfaf0SJed Brown /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called 164248bfaf0SJed Brown * unless we had already been called. */ 16548a46eb9SPierre Jolivet for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i])); 166248bfaf0SJed Brown } 167d0609cedSBarry Smith PetscOptionsHeadEnd(); 1684b9ad928SBarry Smith PetscFunctionReturn(0); 1694b9ad928SBarry Smith } 1704b9ad928SBarry Smith 1719804daf3SBarry Smith #include <petscdraw.h> 1729371c9d4SSatish Balay static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer) { 1734b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 174cbe18068SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 17569a612a9SBarry Smith PetscMPIInt rank; 17669a612a9SBarry Smith PetscInt i; 177d9884438SBarry Smith PetscBool iascii, isstring, isdraw; 1784b9ad928SBarry Smith PetscViewer sviewer; 179020d6619SPierre Jolivet PetscViewerFormat format; 180020d6619SPierre Jolivet const char *prefix; 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 18632077d6dSBarry Smith if (iascii) { 18748a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n)); 18863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " number of blocks = %" PetscInt_FMT "\n", jac->n)); 1899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 1909566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 191020d6619SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 1929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 1939566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 195933995ecSHong Zhang if (jac->ksp && !jac->psubcomm) { 1969566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 197dd400576SPatrick Sanan if (rank == 0) { 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1999566063dSJacob Faibussowitsch PetscCall(KSPView(jac->ksp[0], sviewer)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2014b9ad928SBarry Smith } 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(sviewer)); 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 205e4de9e1dSBarry Smith /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 2069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 207e4fa1014SBarry Smith } else if (mpjac && jac->ksp && mpjac->psubcomm) { 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer)); 209e4fa1014SBarry Smith if (!mpjac->psubcomm->color) { 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2119566063dSJacob Faibussowitsch PetscCall(KSPView(*(jac->ksp), sviewer)); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 213cbe18068SHong Zhang } 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(sviewer)); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer)); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2179530cbd7SBarry Smith /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2194b9ad928SBarry Smith } else { 2209566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 221e4de9e1dSBarry Smith } 222e4de9e1dSBarry Smith } else { 2234814766eSBarry Smith PetscInt n_global; 2241c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following KSP and PC objects:\n")); 2279371c9d4SSatish Balay PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local)); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 230dd2fa690SBarry Smith for (i = 0; i < jac->n_local; i++) { 23163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i)); 2329566063dSJacob Faibussowitsch PetscCall(KSPView(jac->ksp[i], sviewer)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 2344b9ad928SBarry Smith } 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith } else if (isstring) { 24163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2439566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer)); 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 245d9884438SBarry Smith } else if (isdraw) { 246d9884438SBarry Smith PetscDraw draw; 247d9884438SBarry Smith char str[25]; 248d9884438SBarry Smith PetscReal x, y, bottom, h; 249d9884438SBarry Smith 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 2519566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 25263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n)); 2539566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 254d9884438SBarry Smith bottom = y - h; 2559566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 256d9884438SBarry Smith /* warning the communicator on viewer is different then on ksp in parallel */ 2579566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer)); 2589566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 2594b9ad928SBarry Smith } 2604b9ad928SBarry Smith PetscFunctionReturn(0); 2614b9ad928SBarry Smith } 2624b9ad928SBarry Smith 2639371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) { 264feb237baSPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data; 2654b9ad928SBarry Smith 2664b9ad928SBarry Smith PetscFunctionBegin; 26728b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first"); 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith if (n_local) *n_local = jac->n_local; 2704b9ad928SBarry Smith if (first_local) *first_local = jac->first_local; 271020d6619SPierre Jolivet if (ksp) *ksp = jac->ksp; 2724b9ad928SBarry Smith PetscFunctionReturn(0); 2734b9ad928SBarry Smith } 2744b9ad928SBarry Smith 2759371c9d4SSatish Balay static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, PetscInt *lens) { 2764b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith PetscFunctionBegin; 2792472a847SBarry Smith PetscCheck(pc->setupcalled <= 0 || jac->n == blocks, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called"); 2804b9ad928SBarry Smith jac->n = blocks; 2810a545947SLisandro Dalcin if (!lens) jac->g_lens = NULL; 2822fa5cd67SKarl Rupp else { 2839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(blocks, &jac->g_lens)); 2849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->g_lens, lens, blocks)); 2854b9ad928SBarry Smith } 2864b9ad928SBarry Smith PetscFunctionReturn(0); 2874b9ad928SBarry Smith } 2884b9ad928SBarry Smith 2899371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) { 2904b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith PetscFunctionBegin; 2934b9ad928SBarry Smith *blocks = jac->n; 2944b9ad928SBarry Smith if (lens) *lens = jac->g_lens; 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 2989371c9d4SSatish Balay static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[]) { 2994b9ad928SBarry Smith PC_BJacobi *jac; 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith PetscFunctionBegin; 3024b9ad928SBarry Smith jac = (PC_BJacobi *)pc->data; 3034b9ad928SBarry Smith 3044b9ad928SBarry Smith jac->n_local = blocks; 3050a545947SLisandro Dalcin if (!lens) jac->l_lens = NULL; 3062fa5cd67SKarl Rupp else { 3079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(blocks, &jac->l_lens)); 3089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(jac->l_lens, lens, blocks)); 3094b9ad928SBarry Smith } 3104b9ad928SBarry Smith PetscFunctionReturn(0); 3114b9ad928SBarry Smith } 3124b9ad928SBarry Smith 3139371c9d4SSatish Balay static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) { 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; 3194b9ad928SBarry Smith PetscFunctionReturn(0); 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 3364b9ad928SBarry Smith Notes: 337f1580f4eSBarry Smith After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed. 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Currently for some matrix implementations only 1 block per processor 3404b9ad928SBarry Smith is supported. 3414b9ad928SBarry Smith 342f1580f4eSBarry Smith You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`. 3434b9ad928SBarry Smith 344f1580f4eSBarry Smith Fortran Usage: 345f1580f4eSBarry Smith You must pass in a `KSP` array that is large enough to contain all the local `KSP`s. 346f1580f4eSBarry Smith 347f1580f4eSBarry Smith You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the 348f1580f4eSBarry Smith `KSP` array must be. 349196cc216SBarry Smith 3504b9ad928SBarry Smith Level: advanced 3514b9ad928SBarry Smith 352f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()` 3534b9ad928SBarry Smith @*/ 3549371c9d4SSatish Balay PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) { 3554b9ad928SBarry Smith PetscFunctionBegin; 3560700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 357cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 3584b9ad928SBarry Smith PetscFunctionReturn(0); 3594b9ad928SBarry Smith } 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith /*@ 3624b9ad928SBarry Smith PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block 3634b9ad928SBarry Smith Jacobi preconditioner. 3644b9ad928SBarry Smith 365f1580f4eSBarry Smith Collective on pc 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith Input Parameters: 3684b9ad928SBarry Smith + pc - the preconditioner context 3694b9ad928SBarry Smith . blocks - the number of blocks 3704b9ad928SBarry Smith - lens - [optional] integer array containing the size of each block 3714b9ad928SBarry Smith 3724b9ad928SBarry Smith Options Database Key: 3734b9ad928SBarry Smith . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks 3744b9ad928SBarry 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 3794b9ad928SBarry Smith Level: intermediate 3804b9ad928SBarry Smith 381f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()` 3824b9ad928SBarry Smith @*/ 3839371c9d4SSatish Balay PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) { 3844b9ad928SBarry Smith PetscFunctionBegin; 3850700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 38608401ef6SPierre Jolivet PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks"); 387cac4c232SBarry Smith PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens)); 3884b9ad928SBarry Smith PetscFunctionReturn(0); 3894b9ad928SBarry Smith } 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith /*@C 3924b9ad928SBarry Smith PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block 393f1580f4eSBarry Smith Jacobi, `PCBJACOBI`, preconditioner. 3944b9ad928SBarry Smith 395ad4df100SBarry Smith Not Collective 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith Input Parameter: 3984b9ad928SBarry Smith . pc - the preconditioner context 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith Output parameters: 4014b9ad928SBarry Smith + blocks - the number of blocks 4024b9ad928SBarry Smith - lens - integer array containing the size of each block 4034b9ad928SBarry Smith 4044b9ad928SBarry Smith Level: intermediate 4054b9ad928SBarry Smith 406f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()` 4074b9ad928SBarry Smith @*/ 4089371c9d4SSatish Balay PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) { 4094b9ad928SBarry Smith PetscFunctionBegin; 4100700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4114482741eSBarry Smith PetscValidIntPointer(blocks, 2); 412cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens)); 4134b9ad928SBarry Smith PetscFunctionReturn(0); 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 4304b9ad928SBarry Smith Note: 4314b9ad928SBarry Smith Currently only a limited number of blocking configurations are supported. 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith Level: intermediate 4344b9ad928SBarry Smith 435f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()` 4364b9ad928SBarry Smith @*/ 4379371c9d4SSatish Balay PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) { 4384b9ad928SBarry Smith PetscFunctionBegin; 4390700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 44008401ef6SPierre Jolivet PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks"); 441cac4c232SBarry Smith PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens)); 4424b9ad928SBarry Smith PetscFunctionReturn(0); 4434b9ad928SBarry Smith } 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith /*@C 4464b9ad928SBarry Smith PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block 447f1580f4eSBarry Smith Jacobi, `PCBJACOBI`, preconditioner. 4484b9ad928SBarry Smith 4494b9ad928SBarry Smith Not Collective 4504b9ad928SBarry Smith 4514b9ad928SBarry Smith Input Parameters: 4524b9ad928SBarry Smith + pc - the preconditioner context 4534b9ad928SBarry Smith . blocks - the number of blocks 4544b9ad928SBarry Smith - lens - [optional] integer array containing size of each block 4554b9ad928SBarry Smith 4564b9ad928SBarry Smith Note: 4574b9ad928SBarry Smith Currently only a limited number of blocking configurations are supported. 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith Level: intermediate 4604b9ad928SBarry Smith 461f1580f4eSBarry Smith .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()` 4624b9ad928SBarry Smith @*/ 4639371c9d4SSatish Balay PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) { 4644b9ad928SBarry Smith PetscFunctionBegin; 4650700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4664482741eSBarry Smith PetscValidIntPointer(blocks, 2); 467cac4c232SBarry Smith PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens)); 4684b9ad928SBarry Smith PetscFunctionReturn(0); 4694b9ad928SBarry Smith } 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith /*MC 4724b9ad928SBarry Smith PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with 473f1580f4eSBarry Smith its own `KSP` object. 4744b9ad928SBarry Smith 4754b9ad928SBarry Smith Options Database Keys: 476011ea8aeSBarry Smith + -pc_use_amat - use Amat to apply block of operator in inner Krylov method 477011ea8aeSBarry Smith - -pc_bjacobi_blocks <n> - use n total blocks 4784b9ad928SBarry Smith 47995452b02SPatrick Sanan Notes: 480f1580f4eSBarry Smith See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block 481468ae2e8SBarry Smith 48295452b02SPatrick 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. 4834b9ad928SBarry Smith 484f1580f4eSBarry Smith To set options on the solvers for each block append -sub_ to all the `KSP` and `PC` 485d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 4864b9ad928SBarry Smith 487f1580f4eSBarry Smith To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()` 488f1580f4eSBarry Smith and set the options directly on the resulting `KSP` object (you can access its `PC` 489f1580f4eSBarry Smith `KSPGetPC())` 4904b9ad928SBarry Smith 491f1580f4eSBarry Smith For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best 4922210c163SDominic Meiser performance. Different block partitioning may lead to additional data transfers 4932210c163SDominic Meiser between host and GPU that lead to degraded performance. 4942210c163SDominic Meiser 495011ea8aeSBarry Smith When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes. 496011ea8aeSBarry Smith 4974b9ad928SBarry Smith Level: beginner 4984b9ad928SBarry Smith 499f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`, 500db781477SPatrick Sanan `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`, 501db781477SPatrick Sanan `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI` 5024b9ad928SBarry Smith M*/ 5034b9ad928SBarry Smith 5049371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc) { 50569a612a9SBarry Smith PetscMPIInt rank; 5064b9ad928SBarry Smith PC_BJacobi *jac; 5074b9ad928SBarry Smith 5084b9ad928SBarry Smith PetscFunctionBegin; 509*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 5109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 5112fa5cd67SKarl Rupp 5120a545947SLisandro Dalcin pc->ops->apply = NULL; 5137b6e2003SPierre Jolivet pc->ops->matapply = NULL; 5140a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 5154b9ad928SBarry Smith pc->ops->setup = PCSetUp_BJacobi; 5164b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi; 5174b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_BJacobi; 5184b9ad928SBarry Smith pc->ops->view = PCView_BJacobi; 5190a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 5204b9ad928SBarry Smith 5214b9ad928SBarry Smith pc->data = (void *)jac; 5224b9ad928SBarry Smith jac->n = -1; 5234b9ad928SBarry Smith jac->n_local = -1; 5244b9ad928SBarry Smith jac->first_local = rank; 5250a545947SLisandro Dalcin jac->ksp = NULL; 5260a545947SLisandro Dalcin jac->g_lens = NULL; 5270a545947SLisandro Dalcin jac->l_lens = NULL; 5280a545947SLisandro Dalcin jac->psubcomm = NULL; 5294b9ad928SBarry Smith 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi)); 5319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi)); 5329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi)); 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi)); 5349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi)); 5354b9ad928SBarry Smith PetscFunctionReturn(0); 5364b9ad928SBarry Smith } 5374b9ad928SBarry Smith 5384b9ad928SBarry Smith /* 5394b9ad928SBarry Smith These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI 5404b9ad928SBarry Smith */ 5419371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc) { 5424b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 5434b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 5444b9ad928SBarry Smith 5454b9ad928SBarry Smith PetscFunctionBegin; 5469566063dSJacob Faibussowitsch PetscCall(KSPReset(jac->ksp[0])); 5479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->x)); 5489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->y)); 549e91c6855SBarry Smith PetscFunctionReturn(0); 550e91c6855SBarry Smith } 551e91c6855SBarry Smith 5529371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) { 553e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 554e91c6855SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 555e91c6855SBarry Smith 556e91c6855SBarry Smith PetscFunctionBegin; 5579566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Singleblock(pc)); 5589566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->ksp[0])); 5599566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp)); 5609566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac)); 5612e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc)); 5624b9ad928SBarry Smith PetscFunctionReturn(0); 5634b9ad928SBarry Smith } 5644b9ad928SBarry Smith 5659371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) { 5664b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 5672295b7c8SHong Zhang KSP subksp = jac->ksp[0]; 568539c167fSBarry Smith KSPConvergedReason reason; 5694b9ad928SBarry Smith 5704b9ad928SBarry Smith PetscFunctionBegin; 5719566063dSJacob Faibussowitsch PetscCall(KSPSetUp(subksp)); 5729566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(subksp, &reason)); 573ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 5744b9ad928SBarry Smith PetscFunctionReturn(0); 5754b9ad928SBarry Smith } 5764b9ad928SBarry Smith 5779371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y) { 5784b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 5794b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 5806e111a19SKarl Rupp 5814b9ad928SBarry Smith PetscFunctionBegin; 5829566063dSJacob Faibussowitsch PetscCall(VecGetLocalVectorRead(x, bjac->x)); 5839566063dSJacob Faibussowitsch PetscCall(VecGetLocalVector(y, bjac->y)); 584bba28a21SBarry Smith /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner 585910cf402Sprj- matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild 586910cf402Sprj- of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/ 5879566063dSJacob Faibussowitsch PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner)); 5889566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y)); 5899566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y)); 5909566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVectorRead(x, bjac->x)); 5919566063dSJacob Faibussowitsch PetscCall(VecRestoreLocalVector(y, bjac->y)); 5924b9ad928SBarry Smith PetscFunctionReturn(0); 5934b9ad928SBarry Smith } 5944b9ad928SBarry Smith 5959371c9d4SSatish Balay static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y) { 5967b6e2003SPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data; 5977b6e2003SPierre Jolivet Mat sX, sY; 5987b6e2003SPierre Jolivet 5997b6e2003SPierre Jolivet PetscFunctionBegin; 6007b6e2003SPierre Jolivet /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner 6017b6e2003SPierre Jolivet matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild 6027b6e2003SPierre Jolivet of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/ 6039566063dSJacob Faibussowitsch PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner)); 6049566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(X, &sX)); 6059566063dSJacob Faibussowitsch PetscCall(MatDenseGetLocalMatrix(Y, &sY)); 6069566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(jac->ksp[0], sX, sY)); 6077b6e2003SPierre Jolivet PetscFunctionReturn(0); 6087b6e2003SPierre Jolivet } 6097b6e2003SPierre Jolivet 6109371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y) { 6114b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 6124b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 613d9ca1df4SBarry Smith PetscScalar *y_array; 614d9ca1df4SBarry Smith const PetscScalar *x_array; 6154b9ad928SBarry Smith PC subpc; 6164b9ad928SBarry Smith 6174b9ad928SBarry Smith PetscFunctionBegin; 6184b9ad928SBarry Smith /* 6194b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the 6204b9ad928SBarry Smith y vector into the bjac->x vector. The reason for 6214b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector 6224b9ad928SBarry Smith for the sequential solve. 6234b9ad928SBarry Smith */ 6249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array)); 6259566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array)); 6269566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array)); 6279566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array)); 6284b9ad928SBarry Smith /* apply the symmetric left portion of the inner PC operator */ 6294b9ad928SBarry Smith /* note this by-passes the inner KSP and its options completely */ 6309566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &subpc)); 6319566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y)); 6329566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x)); 6339566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y)); 6349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array)); 6359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array)); 6364b9ad928SBarry Smith PetscFunctionReturn(0); 6374b9ad928SBarry Smith } 6384b9ad928SBarry Smith 6399371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y) { 6404b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 6414b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 642d9ca1df4SBarry Smith PetscScalar *y_array; 643d9ca1df4SBarry Smith const PetscScalar *x_array; 6444b9ad928SBarry Smith PC subpc; 6454b9ad928SBarry Smith 6464b9ad928SBarry Smith PetscFunctionBegin; 6474b9ad928SBarry Smith /* 6484b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the 6494b9ad928SBarry Smith y vector into the bjac->x vector. The reason for 6504b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector 6514b9ad928SBarry Smith for the sequential solve. 6524b9ad928SBarry Smith */ 6539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array)); 6549566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array)); 6559566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array)); 6569566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array)); 6574b9ad928SBarry Smith 6584b9ad928SBarry Smith /* apply the symmetric right portion of the inner PC operator */ 6594b9ad928SBarry Smith /* note this by-passes the inner KSP and its options completely */ 6604b9ad928SBarry Smith 6619566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &subpc)); 6629566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y)); 6634b9ad928SBarry Smith 6649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array)); 6659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array)); 6664b9ad928SBarry Smith PetscFunctionReturn(0); 6674b9ad928SBarry Smith } 6684b9ad928SBarry Smith 6699371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y) { 6704b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 6714b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 672d9ca1df4SBarry Smith PetscScalar *y_array; 673d9ca1df4SBarry Smith const PetscScalar *x_array; 6744b9ad928SBarry Smith 6754b9ad928SBarry Smith PetscFunctionBegin; 6764b9ad928SBarry Smith /* 6774b9ad928SBarry Smith The VecPlaceArray() is to avoid having to copy the 6784b9ad928SBarry Smith y vector into the bjac->x vector. The reason for 6794b9ad928SBarry Smith the bjac->x vector is that we need a sequential vector 6804b9ad928SBarry Smith for the sequential solve. 6814b9ad928SBarry Smith */ 6829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &x_array)); 6839566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array)); 6849566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x, x_array)); 6859566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y, y_array)); 6869566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y)); 6879566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y)); 6889566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x)); 6899566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y)); 6909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &x_array)); 6919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array)); 6924b9ad928SBarry Smith PetscFunctionReturn(0); 6934b9ad928SBarry Smith } 6944b9ad928SBarry Smith 6959371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat) { 6964b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 69769a612a9SBarry Smith PetscInt m; 6984b9ad928SBarry Smith KSP ksp; 6994b9ad928SBarry Smith PC_BJacobi_Singleblock *bjac; 700de0953f6SBarry Smith PetscBool wasSetup = PETSC_TRUE; 7013f6dc190SJunchao Zhang VecType vectype; 702ea41da7aSPierre Jolivet const char *prefix; 7034b9ad928SBarry Smith 7044b9ad928SBarry Smith PetscFunctionBegin; 7054b9ad928SBarry Smith if (!pc->setupcalled) { 706c2efdce3SBarry Smith if (!jac->ksp) { 707302a9ddcSMatthew Knepley wasSetup = PETSC_FALSE; 7082fa5cd67SKarl Rupp 7099566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 7109566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 7129566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY)); 7139566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 7149566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 7159566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 7164b9ad928SBarry Smith 717e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Singleblock; 7184b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi_Singleblock; 7194b9ad928SBarry Smith pc->ops->apply = PCApply_BJacobi_Singleblock; 7207b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_BJacobi_Singleblock; 7214b9ad928SBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock; 7224b9ad928SBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock; 7234b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock; 7244b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock; 7254b9ad928SBarry Smith 7269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->ksp)); 7274b9ad928SBarry Smith jac->ksp[0] = ksp; 728c2efdce3SBarry Smith 729*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bjac)); 7304b9ad928SBarry Smith jac->data = (void *)bjac; 7314b9ad928SBarry Smith } else { 732c2efdce3SBarry Smith ksp = jac->ksp[0]; 733c2efdce3SBarry Smith bjac = (PC_BJacobi_Singleblock *)jac->data; 734c2efdce3SBarry Smith } 735c2efdce3SBarry Smith 736c2efdce3SBarry Smith /* 737c2efdce3SBarry Smith The reason we need to generate these vectors is to serve 738c2efdce3SBarry Smith as the right-hand side and solution vector for the solve on the 739c2efdce3SBarry Smith block. We do not need to allocate space for the vectors since 740c2efdce3SBarry Smith that is provided via VecPlaceArray() just before the call to 741c2efdce3SBarry Smith KSPSolve() on the block. 742c2efdce3SBarry Smith */ 7439566063dSJacob Faibussowitsch PetscCall(MatGetSize(pmat, &m, &m)); 7449566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x)); 7459566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y)); 7469566063dSJacob Faibussowitsch PetscCall(MatGetVecType(pmat, &vectype)); 7479566063dSJacob Faibussowitsch PetscCall(VecSetType(bjac->x, vectype)); 7489566063dSJacob Faibussowitsch PetscCall(VecSetType(bjac->y, vectype)); 749c2efdce3SBarry Smith } else { 7504b9ad928SBarry Smith ksp = jac->ksp[0]; 7514b9ad928SBarry Smith bjac = (PC_BJacobi_Singleblock *)jac->data; 7524b9ad928SBarry Smith } 7539566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(ksp, &prefix)); 75449517cdeSBarry Smith if (pc->useAmat) { 7559566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, mat, pmat)); 7569566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(mat, prefix)); 7574b9ad928SBarry Smith } else { 7589566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, pmat, pmat)); 7594b9ad928SBarry Smith } 7609566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(pmat, prefix)); 76190f1c854SHong Zhang if (!wasSetup && pc->setfromoptionscalled) { 762248bfaf0SJed Brown /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */ 7639566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 764302a9ddcSMatthew Knepley } 7654b9ad928SBarry Smith PetscFunctionReturn(0); 7664b9ad928SBarry Smith } 7674b9ad928SBarry Smith 7689371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc) { 7694b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 7704b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 77169a612a9SBarry Smith PetscInt i; 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith PetscFunctionBegin; 774e91c6855SBarry Smith if (bjac && bjac->pmat) { 7759566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat)); 77648a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat)); 777e91c6855SBarry Smith } 7784b9ad928SBarry Smith 7794b9ad928SBarry Smith for (i = 0; i < jac->n_local; i++) { 7809566063dSJacob Faibussowitsch PetscCall(KSPReset(jac->ksp[i])); 781e91c6855SBarry Smith if (bjac && bjac->x) { 7829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->x[i])); 7839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bjac->y[i])); 7849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bjac->is[i])); 7854b9ad928SBarry Smith } 786e91c6855SBarry Smith } 7879566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->l_lens)); 7889566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->g_lens)); 789e91c6855SBarry Smith PetscFunctionReturn(0); 790e91c6855SBarry Smith } 791e91c6855SBarry Smith 7929371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) { 793e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 794c2efdce3SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 795e91c6855SBarry Smith PetscInt i; 796e91c6855SBarry Smith 797e91c6855SBarry Smith PetscFunctionBegin; 7989566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Multiblock(pc)); 799c2efdce3SBarry Smith if (bjac) { 8009566063dSJacob Faibussowitsch PetscCall(PetscFree2(bjac->x, bjac->y)); 8019566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac->starts)); 8029566063dSJacob Faibussowitsch PetscCall(PetscFree(bjac->is)); 803c2efdce3SBarry Smith } 8049566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->data)); 80548a46eb9SPierre Jolivet for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i])); 8069566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp)); 8072e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc)); 8084b9ad928SBarry Smith PetscFunctionReturn(0); 8094b9ad928SBarry Smith } 8104b9ad928SBarry Smith 8119371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) { 8124b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 81369a612a9SBarry Smith PetscInt i, n_local = jac->n_local; 814539c167fSBarry Smith KSPConvergedReason reason; 8154b9ad928SBarry Smith 8164b9ad928SBarry Smith PetscFunctionBegin; 8174b9ad928SBarry Smith for (i = 0; i < n_local; i++) { 8189566063dSJacob Faibussowitsch PetscCall(KSPSetUp(jac->ksp[i])); 8199566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason)); 820ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 8214b9ad928SBarry Smith } 8224b9ad928SBarry Smith PetscFunctionReturn(0); 8234b9ad928SBarry Smith } 8244b9ad928SBarry Smith 8254b9ad928SBarry Smith /* 8264b9ad928SBarry Smith Preconditioner for block Jacobi 8274b9ad928SBarry Smith */ 8289371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y) { 8294b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 83069a612a9SBarry Smith PetscInt i, n_local = jac->n_local; 8314b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 832d9ca1df4SBarry Smith PetscScalar *yin; 833d9ca1df4SBarry Smith const PetscScalar *xin; 83458ebbce7SBarry Smith 8354b9ad928SBarry Smith PetscFunctionBegin; 8369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xin)); 8379566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yin)); 8384b9ad928SBarry Smith for (i = 0; i < n_local; i++) { 8394b9ad928SBarry Smith /* 8404b9ad928SBarry Smith To avoid copying the subvector from x into a workspace we instead 8414b9ad928SBarry Smith make the workspace vector array point to the subpart of the array of 8424b9ad928SBarry Smith the global vector. 8434b9ad928SBarry Smith */ 8449566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i])); 8459566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i])); 8464b9ad928SBarry Smith 8479566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 8489566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i])); 8499566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i])); 8509566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 851d11f3a42SBarry Smith 8529566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x[i])); 8539566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y[i])); 8544b9ad928SBarry Smith } 8559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xin)); 8569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yin)); 8574b9ad928SBarry Smith PetscFunctionReturn(0); 8584b9ad928SBarry Smith } 8594b9ad928SBarry Smith 8604b9ad928SBarry Smith /* 8614b9ad928SBarry Smith Preconditioner for block Jacobi 8624b9ad928SBarry Smith */ 8639371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y) { 8644b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 86569a612a9SBarry Smith PetscInt i, n_local = jac->n_local; 8664b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 867d9ca1df4SBarry Smith PetscScalar *yin; 868d9ca1df4SBarry Smith const PetscScalar *xin; 8694b9ad928SBarry Smith 8704b9ad928SBarry Smith PetscFunctionBegin; 8719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xin)); 8729566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yin)); 8734b9ad928SBarry Smith for (i = 0; i < n_local; i++) { 8744b9ad928SBarry Smith /* 8754b9ad928SBarry Smith To avoid copying the subvector from x into a workspace we instead 8764b9ad928SBarry Smith make the workspace vector array point to the subpart of the array of 8774b9ad928SBarry Smith the global vector. 8784b9ad928SBarry Smith */ 8799566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i])); 8809566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i])); 8814b9ad928SBarry Smith 8829566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 8839566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i])); 8849566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i])); 8859566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 886a7ff49e8SJed Brown 8879566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->x[i])); 8889566063dSJacob Faibussowitsch PetscCall(VecResetArray(bjac->y[i])); 8894b9ad928SBarry Smith } 8909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xin)); 8919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yin)); 8924b9ad928SBarry Smith PetscFunctionReturn(0); 8934b9ad928SBarry Smith } 8944b9ad928SBarry Smith 8959371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat) { 8964b9ad928SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 89769a612a9SBarry Smith PetscInt m, n_local, N, M, start, i; 898ea41da7aSPierre Jolivet const char *prefix; 8994b9ad928SBarry Smith KSP ksp; 9004b9ad928SBarry Smith Vec x, y; 9014b9ad928SBarry Smith PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 9024b9ad928SBarry Smith PC subpc; 9034b9ad928SBarry Smith IS is; 904434a97beSBrad Aagaard MatReuse scall; 9053f6dc190SJunchao Zhang VecType vectype; 9064b9ad928SBarry Smith 9074b9ad928SBarry Smith PetscFunctionBegin; 9089566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &M, &N)); 9094b9ad928SBarry Smith 9104b9ad928SBarry Smith n_local = jac->n_local; 9114b9ad928SBarry Smith 91249517cdeSBarry Smith if (pc->useAmat) { 913ace3abfcSBarry Smith PetscBool same; 9149566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same)); 91528b400f6SJacob Faibussowitsch PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type"); 9164b9ad928SBarry Smith } 9174b9ad928SBarry Smith 9184b9ad928SBarry Smith if (!pc->setupcalled) { 9194b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 920c2efdce3SBarry Smith 921c2efdce3SBarry Smith if (!jac->ksp) { 922e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Multiblock; 9234b9ad928SBarry Smith pc->ops->destroy = PCDestroy_BJacobi_Multiblock; 9244b9ad928SBarry Smith pc->ops->apply = PCApply_BJacobi_Multiblock; 9257b6e2003SPierre Jolivet pc->ops->matapply = NULL; 9264b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_BJacobi_Multiblock; 9274b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock; 9284b9ad928SBarry Smith 929*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&bjac)); 9309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &jac->ksp)); 9319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y)); 9329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &bjac->starts)); 9334b9ad928SBarry Smith 9344b9ad928SBarry Smith jac->data = (void *)bjac; 9359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_local, &bjac->is)); 9364b9ad928SBarry Smith 9374b9ad928SBarry Smith for (i = 0; i < n_local; i++) { 9389566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 9399566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 9409566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 9419566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY)); 9429566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &subpc)); 9439566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 9449566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 9459566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 9462fa5cd67SKarl Rupp 947c2efdce3SBarry Smith jac->ksp[i] = ksp; 948c2efdce3SBarry Smith } 949c2efdce3SBarry Smith } else { 950c2efdce3SBarry Smith bjac = (PC_BJacobi_Multiblock *)jac->data; 951c2efdce3SBarry Smith } 9524b9ad928SBarry Smith 953c2efdce3SBarry Smith start = 0; 9549566063dSJacob Faibussowitsch PetscCall(MatGetVecType(pmat, &vectype)); 955c2efdce3SBarry Smith for (i = 0; i < n_local; i++) { 9564b9ad928SBarry Smith m = jac->l_lens[i]; 9574b9ad928SBarry Smith /* 9584b9ad928SBarry Smith The reason we need to generate these vectors is to serve 9594b9ad928SBarry Smith as the right-hand side and solution vector for the solve on the 9604b9ad928SBarry Smith block. We do not need to allocate space for the vectors since 9614b9ad928SBarry Smith that is provided via VecPlaceArray() just before the call to 9624b9ad928SBarry Smith KSPSolve() on the block. 9634b9ad928SBarry Smith 9644b9ad928SBarry Smith */ 9659566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x)); 9669566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y)); 9679566063dSJacob Faibussowitsch PetscCall(VecSetType(x, vectype)); 9689566063dSJacob Faibussowitsch PetscCall(VecSetType(y, vectype)); 9692fa5cd67SKarl Rupp 9704b9ad928SBarry Smith bjac->x[i] = x; 9714b9ad928SBarry Smith bjac->y[i] = y; 9724b9ad928SBarry Smith bjac->starts[i] = start; 9734b9ad928SBarry Smith 9749566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is)); 9754b9ad928SBarry Smith bjac->is[i] = is; 9764b9ad928SBarry Smith 9774b9ad928SBarry Smith start += m; 9784b9ad928SBarry Smith } 9794b9ad928SBarry Smith } else { 9804b9ad928SBarry Smith bjac = (PC_BJacobi_Multiblock *)jac->data; 9814b9ad928SBarry Smith /* 9824b9ad928SBarry Smith Destroy the blocks from the previous iteration 9834b9ad928SBarry Smith */ 9844b9ad928SBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 9859566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(n_local, &bjac->pmat)); 98648a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat)); 9874b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 9882fa5cd67SKarl Rupp } else scall = MAT_REUSE_MATRIX; 9894b9ad928SBarry Smith } 9904b9ad928SBarry Smith 9919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat)); 99248a46eb9SPierre Jolivet if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat)); 9934b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 9944b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 9959566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP)); 9964b9ad928SBarry Smith 9974b9ad928SBarry Smith for (i = 0; i < n_local; i++) { 9989566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix)); 99949517cdeSBarry Smith if (pc->useAmat) { 10009566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i])); 10019566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix)); 10024b9ad928SBarry Smith } else { 10039566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i])); 10044b9ad928SBarry Smith } 10059566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix)); 100648a46eb9SPierre Jolivet if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i])); 100790f1c854SHong Zhang } 10084b9ad928SBarry Smith PetscFunctionReturn(0); 10094b9ad928SBarry Smith } 10105a7e1895SHong Zhang 10115a7e1895SHong Zhang /* 1012fd0b8932SBarry Smith These are for a single block with multiple processes 10135a7e1895SHong Zhang */ 10149371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc) { 1015fd0b8932SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1016fd0b8932SBarry Smith KSP subksp = jac->ksp[0]; 1017fd0b8932SBarry Smith KSPConvergedReason reason; 1018fd0b8932SBarry Smith 1019fd0b8932SBarry Smith PetscFunctionBegin; 10209566063dSJacob Faibussowitsch PetscCall(KSPSetUp(subksp)); 10219566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(subksp, &reason)); 1022ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1023fd0b8932SBarry Smith PetscFunctionReturn(0); 1024fd0b8932SBarry Smith } 1025fd0b8932SBarry Smith 10269371c9d4SSatish Balay static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc) { 10275a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data; 10285a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 10295a7e1895SHong Zhang 10305a7e1895SHong Zhang PetscFunctionBegin; 10319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mpjac->ysub)); 10329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mpjac->xsub)); 10339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mpjac->submats)); 10349566063dSJacob Faibussowitsch if (jac->ksp) PetscCall(KSPReset(jac->ksp[0])); 1035e91c6855SBarry Smith PetscFunctionReturn(0); 1036e91c6855SBarry Smith } 1037e91c6855SBarry Smith 10389371c9d4SSatish Balay static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc) { 1039e91c6855SBarry Smith PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1040e91c6855SBarry Smith PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1041e91c6855SBarry Smith 1042e91c6855SBarry Smith PetscFunctionBegin; 10439566063dSJacob Faibussowitsch PetscCall(PCReset_BJacobi_Multiproc(pc)); 10449566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&jac->ksp[0])); 10459566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->ksp)); 10469566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&mpjac->psubcomm)); 10475a7e1895SHong Zhang 10489566063dSJacob Faibussowitsch PetscCall(PetscFree(mpjac)); 10492e956fe4SStefano Zampini PetscCall(PCDestroy_BJacobi(pc)); 10505a7e1895SHong Zhang PetscFunctionReturn(0); 10515a7e1895SHong Zhang } 10525a7e1895SHong Zhang 10539371c9d4SSatish Balay static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y) { 10545a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data; 10555a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1056d9ca1df4SBarry Smith PetscScalar *yarray; 1057d9ca1df4SBarry Smith const PetscScalar *xarray; 1058539c167fSBarry Smith KSPConvergedReason reason; 10595a7e1895SHong Zhang 10605a7e1895SHong Zhang PetscFunctionBegin; 10615a7e1895SHong Zhang /* place x's and y's local arrays into xsub and ysub */ 10629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xarray)); 10639566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yarray)); 10649566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mpjac->xsub, xarray)); 10659566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(mpjac->ysub, yarray)); 10665a7e1895SHong Zhang 10675a7e1895SHong Zhang /* apply preconditioner on each matrix block */ 10689566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0)); 10699566063dSJacob Faibussowitsch PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub)); 10709566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub)); 10719566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0)); 10729566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason)); 1073ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1074250cf88bSHong Zhang 10759566063dSJacob Faibussowitsch PetscCall(VecResetArray(mpjac->xsub)); 10769566063dSJacob Faibussowitsch PetscCall(VecResetArray(mpjac->ysub)); 10779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xarray)); 10789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yarray)); 10795a7e1895SHong Zhang PetscFunctionReturn(0); 10805a7e1895SHong Zhang } 10815a7e1895SHong Zhang 10829371c9d4SSatish Balay static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y) { 10837b6e2003SPierre Jolivet PC_BJacobi *jac = (PC_BJacobi *)pc->data; 10847b6e2003SPierre Jolivet KSPConvergedReason reason; 10857b6e2003SPierre Jolivet Mat sX, sY; 10867b6e2003SPierre Jolivet const PetscScalar *x; 10877b6e2003SPierre Jolivet PetscScalar *y; 10887b6e2003SPierre Jolivet PetscInt m, N, lda, ldb; 10897b6e2003SPierre Jolivet 10907b6e2003SPierre Jolivet PetscFunctionBegin; 10917b6e2003SPierre Jolivet /* apply preconditioner on each matrix block */ 10929566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m, NULL)); 10939566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N)); 10949566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda)); 10959566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(Y, &ldb)); 10969566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X, &x)); 10979566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(Y, &y)); 10989566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX)); 10999566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY)); 11009566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(sX, lda)); 11019566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(sY, ldb)); 11029566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0)); 11039566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(jac->ksp[0], sX, sY)); 11049566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL)); 11059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0)); 11069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sY)); 11079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sX)); 11089566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(Y, &y)); 11099566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X, &x)); 11109566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason)); 1111ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 11127b6e2003SPierre Jolivet PetscFunctionReturn(0); 11137b6e2003SPierre Jolivet } 11147b6e2003SPierre Jolivet 11159371c9d4SSatish Balay static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc) { 11165a7e1895SHong Zhang PC_BJacobi *jac = (PC_BJacobi *)pc->data; 11175a7e1895SHong Zhang PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 11185a7e1895SHong Zhang PetscInt m, n; 1119ce94432eSBarry Smith MPI_Comm comm, subcomm = 0; 11205a7e1895SHong Zhang const char *prefix; 1121de0953f6SBarry Smith PetscBool wasSetup = PETSC_TRUE; 11223f6dc190SJunchao Zhang VecType vectype; 11231f62890fSKarl Rupp 11245a7e1895SHong Zhang PetscFunctionBegin; 11259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 112608401ef6SPierre Jolivet PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported"); 11275a7e1895SHong Zhang jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */ 11285a7e1895SHong Zhang if (!pc->setupcalled) { 1129de0953f6SBarry Smith wasSetup = PETSC_FALSE; 1130*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&mpjac)); 11315a7e1895SHong Zhang jac->data = (void *)mpjac; 11325a7e1895SHong Zhang 11335a7e1895SHong Zhang /* initialize datastructure mpjac */ 11345a7e1895SHong Zhang if (!jac->psubcomm) { 11355a7e1895SHong Zhang /* Create default contiguous subcommunicatiors if user does not provide them */ 11369566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(comm, &jac->psubcomm)); 11379566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n)); 11389566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS)); 11395a7e1895SHong Zhang } 11405a7e1895SHong Zhang mpjac->psubcomm = jac->psubcomm; 1141306c2d5bSBarry Smith subcomm = PetscSubcommChild(mpjac->psubcomm); 11425a7e1895SHong Zhang 11435a7e1895SHong Zhang /* Get matrix blocks of pmat */ 11449566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats)); 11455a7e1895SHong Zhang 11465a7e1895SHong Zhang /* create a new PC that processors in each subcomm have copy of */ 11479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &jac->ksp)); 11489566063dSJacob Faibussowitsch PetscCall(KSPCreate(subcomm, &jac->ksp[0])); 11499566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure)); 11509566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1)); 11519566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats)); 11529566063dSJacob Faibussowitsch PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc)); 11535a7e1895SHong Zhang 11549566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 11559566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix)); 11569566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_")); 11579566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix)); 11589566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix)); 11595a7e1895SHong Zhang 11605a7e1895SHong Zhang /* create dummy vectors xsub and ysub */ 11619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(mpjac->submats, &m, &n)); 11629566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub)); 11639566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub)); 11649566063dSJacob Faibussowitsch PetscCall(MatGetVecType(mpjac->submats, &vectype)); 11659566063dSJacob Faibussowitsch PetscCall(VecSetType(mpjac->xsub, vectype)); 11669566063dSJacob Faibussowitsch PetscCall(VecSetType(mpjac->ysub, vectype)); 11675a7e1895SHong Zhang 1168fd0b8932SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc; 1169e91c6855SBarry Smith pc->ops->reset = PCReset_BJacobi_Multiproc; 11705a7e1895SHong Zhang pc->ops->destroy = PCDestroy_BJacobi_Multiproc; 11715a7e1895SHong Zhang pc->ops->apply = PCApply_BJacobi_Multiproc; 11727b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_BJacobi_Multiproc; 1173fc08c53fSHong Zhang } else { /* pc->setupcalled */ 1174306c2d5bSBarry Smith subcomm = PetscSubcommChild(mpjac->psubcomm); 11759e0ae222SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 11765a7e1895SHong Zhang /* destroy old matrix blocks, then get new matrix blocks */ 11779566063dSJacob Faibussowitsch if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats)); 11789566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats)); 1179fc08c53fSHong Zhang } else { 11809566063dSJacob Faibussowitsch PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats)); 11815a7e1895SHong Zhang } 11829566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats)); 11835a7e1895SHong Zhang } 11845a7e1895SHong Zhang 118548a46eb9SPierre Jolivet if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0])); 11865a7e1895SHong Zhang PetscFunctionReturn(0); 11875a7e1895SHong Zhang } 1188