1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Include files needed for the PBJacobi preconditioner: 44b9ad928SBarry Smith pcimpl.h - private include file intended for use by all preconditioners 54b9ad928SBarry Smith */ 64b9ad928SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 84b9ad928SBarry Smith 94b9ad928SBarry Smith /* 104b9ad928SBarry Smith Private context (data structure) for the PBJacobi preconditioner. 114b9ad928SBarry Smith */ 124b9ad928SBarry Smith typedef struct { 13713ccfa9SJed Brown const MatScalar *diag; 14c1ac3661SBarry Smith PetscInt bs, mbs; 154b9ad928SBarry Smith } PC_PBJacobi; 164b9ad928SBarry Smith 179371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_1(PC pc, Vec x, Vec y) { 18bbead8a2SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 19bbead8a2SBarry Smith PetscInt i, m = jac->mbs; 20bbead8a2SBarry Smith const MatScalar *diag = jac->diag; 21bbead8a2SBarry Smith const PetscScalar *xx; 22bbead8a2SBarry Smith PetscScalar *yy; 23bbead8a2SBarry Smith 24bbead8a2SBarry Smith PetscFunctionBegin; 259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 269566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 272fa5cd67SKarl Rupp for (i = 0; i < m; i++) yy[i] = diag[i] * xx[i]; 289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 309566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(m)); 31bbead8a2SBarry Smith PetscFunctionReturn(0); 32bbead8a2SBarry Smith } 33bbead8a2SBarry Smith 349371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_2(PC pc, Vec x, Vec y) { 354b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 36c1ac3661SBarry Smith PetscInt i, m = jac->mbs; 3785f4f44aSBarry Smith const MatScalar *diag = jac->diag; 380e666e78SToby Isaac PetscScalar x0, x1, *yy; 390e666e78SToby Isaac const PetscScalar *xx; 404b9ad928SBarry Smith 414b9ad928SBarry Smith PetscFunctionBegin; 429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 439566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 444b9ad928SBarry Smith for (i = 0; i < m; i++) { 459371c9d4SSatish Balay x0 = xx[2 * i]; 469371c9d4SSatish Balay x1 = xx[2 * i + 1]; 474b9ad928SBarry Smith yy[2 * i] = diag[0] * x0 + diag[2] * x1; 484b9ad928SBarry Smith yy[2 * i + 1] = diag[1] * x0 + diag[3] * x1; 494b9ad928SBarry Smith diag += 4; 504b9ad928SBarry Smith } 519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 539566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(6.0 * m)); 544b9ad928SBarry Smith PetscFunctionReturn(0); 554b9ad928SBarry Smith } 569371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_3(PC pc, Vec x, Vec y) { 574b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 58c1ac3661SBarry Smith PetscInt i, m = jac->mbs; 5985f4f44aSBarry Smith const MatScalar *diag = jac->diag; 600e666e78SToby Isaac PetscScalar x0, x1, x2, *yy; 610e666e78SToby Isaac const PetscScalar *xx; 624b9ad928SBarry Smith 634b9ad928SBarry Smith PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 659566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 664b9ad928SBarry Smith for (i = 0; i < m; i++) { 679371c9d4SSatish Balay x0 = xx[3 * i]; 689371c9d4SSatish Balay x1 = xx[3 * i + 1]; 699371c9d4SSatish Balay x2 = xx[3 * i + 2]; 702fa5cd67SKarl Rupp 714b9ad928SBarry Smith yy[3 * i] = diag[0] * x0 + diag[3] * x1 + diag[6] * x2; 724b9ad928SBarry Smith yy[3 * i + 1] = diag[1] * x0 + diag[4] * x1 + diag[7] * x2; 734b9ad928SBarry Smith yy[3 * i + 2] = diag[2] * x0 + diag[5] * x1 + diag[8] * x2; 744b9ad928SBarry Smith diag += 9; 754b9ad928SBarry Smith } 769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 789566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(15.0 * m)); 794b9ad928SBarry Smith PetscFunctionReturn(0); 804b9ad928SBarry Smith } 819371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_4(PC pc, Vec x, Vec y) { 824b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 83c1ac3661SBarry Smith PetscInt i, m = jac->mbs; 8485f4f44aSBarry Smith const MatScalar *diag = jac->diag; 850e666e78SToby Isaac PetscScalar x0, x1, x2, x3, *yy; 860e666e78SToby Isaac const PetscScalar *xx; 874b9ad928SBarry Smith 884b9ad928SBarry Smith PetscFunctionBegin; 899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 909566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 914b9ad928SBarry Smith for (i = 0; i < m; i++) { 929371c9d4SSatish Balay x0 = xx[4 * i]; 939371c9d4SSatish Balay x1 = xx[4 * i + 1]; 949371c9d4SSatish Balay x2 = xx[4 * i + 2]; 959371c9d4SSatish Balay x3 = xx[4 * i + 3]; 962fa5cd67SKarl Rupp 974b9ad928SBarry Smith yy[4 * i] = diag[0] * x0 + diag[4] * x1 + diag[8] * x2 + diag[12] * x3; 984b9ad928SBarry Smith yy[4 * i + 1] = diag[1] * x0 + diag[5] * x1 + diag[9] * x2 + diag[13] * x3; 994b9ad928SBarry Smith yy[4 * i + 2] = diag[2] * x0 + diag[6] * x1 + diag[10] * x2 + diag[14] * x3; 1004b9ad928SBarry Smith yy[4 * i + 3] = diag[3] * x0 + diag[7] * x1 + diag[11] * x2 + diag[15] * x3; 1014b9ad928SBarry Smith diag += 16; 1024b9ad928SBarry Smith } 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 1059566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(28.0 * m)); 1064b9ad928SBarry Smith PetscFunctionReturn(0); 1074b9ad928SBarry Smith } 1089371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_5(PC pc, Vec x, Vec y) { 1094b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 110c1ac3661SBarry Smith PetscInt i, m = jac->mbs; 11185f4f44aSBarry Smith const MatScalar *diag = jac->diag; 1120e666e78SToby Isaac PetscScalar x0, x1, x2, x3, x4, *yy; 1130e666e78SToby Isaac const PetscScalar *xx; 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith PetscFunctionBegin; 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 1184b9ad928SBarry Smith for (i = 0; i < m; i++) { 1199371c9d4SSatish Balay x0 = xx[5 * i]; 1209371c9d4SSatish Balay x1 = xx[5 * i + 1]; 1219371c9d4SSatish Balay x2 = xx[5 * i + 2]; 1229371c9d4SSatish Balay x3 = xx[5 * i + 3]; 1239371c9d4SSatish Balay x4 = xx[5 * i + 4]; 1242fa5cd67SKarl Rupp 1254b9ad928SBarry Smith yy[5 * i] = diag[0] * x0 + diag[5] * x1 + diag[10] * x2 + diag[15] * x3 + diag[20] * x4; 1264b9ad928SBarry Smith yy[5 * i + 1] = diag[1] * x0 + diag[6] * x1 + diag[11] * x2 + diag[16] * x3 + diag[21] * x4; 1274b9ad928SBarry Smith yy[5 * i + 2] = diag[2] * x0 + diag[7] * x1 + diag[12] * x2 + diag[17] * x3 + diag[22] * x4; 1284b9ad928SBarry Smith yy[5 * i + 3] = diag[3] * x0 + diag[8] * x1 + diag[13] * x2 + diag[18] * x3 + diag[23] * x4; 1294b9ad928SBarry Smith yy[5 * i + 4] = diag[4] * x0 + diag[9] * x1 + diag[14] * x2 + diag[19] * x3 + diag[24] * x4; 1304b9ad928SBarry Smith diag += 25; 1314b9ad928SBarry Smith } 1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 1349566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(45.0 * m)); 1354b9ad928SBarry Smith PetscFunctionReturn(0); 1364b9ad928SBarry Smith } 1379371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_6(PC pc, Vec x, Vec y) { 1380e1b4bd6SMark F. Adams PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 1390e1b4bd6SMark F. Adams PetscInt i, m = jac->mbs; 1400e1b4bd6SMark F. Adams const MatScalar *diag = jac->diag; 1410e666e78SToby Isaac PetscScalar x0, x1, x2, x3, x4, x5, *yy; 1420e666e78SToby Isaac const PetscScalar *xx; 1430e1b4bd6SMark F. Adams 1440e1b4bd6SMark F. Adams PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1469566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 1470e1b4bd6SMark F. Adams for (i = 0; i < m; i++) { 1489371c9d4SSatish Balay x0 = xx[6 * i]; 1499371c9d4SSatish Balay x1 = xx[6 * i + 1]; 1509371c9d4SSatish Balay x2 = xx[6 * i + 2]; 1519371c9d4SSatish Balay x3 = xx[6 * i + 3]; 1529371c9d4SSatish Balay x4 = xx[6 * i + 4]; 1539371c9d4SSatish Balay x5 = xx[6 * i + 5]; 1542fa5cd67SKarl Rupp 1550e1b4bd6SMark F. Adams yy[6 * i] = diag[0] * x0 + diag[6] * x1 + diag[12] * x2 + diag[18] * x3 + diag[24] * x4 + diag[30] * x5; 1560e1b4bd6SMark F. Adams yy[6 * i + 1] = diag[1] * x0 + diag[7] * x1 + diag[13] * x2 + diag[19] * x3 + diag[25] * x4 + diag[31] * x5; 1570e1b4bd6SMark F. Adams yy[6 * i + 2] = diag[2] * x0 + diag[8] * x1 + diag[14] * x2 + diag[20] * x3 + diag[26] * x4 + diag[32] * x5; 1580e1b4bd6SMark F. Adams yy[6 * i + 3] = diag[3] * x0 + diag[9] * x1 + diag[15] * x2 + diag[21] * x3 + diag[27] * x4 + diag[33] * x5; 1590e1b4bd6SMark F. Adams yy[6 * i + 4] = diag[4] * x0 + diag[10] * x1 + diag[16] * x2 + diag[22] * x3 + diag[28] * x4 + diag[34] * x5; 1600e1b4bd6SMark F. Adams yy[6 * i + 5] = diag[5] * x0 + diag[11] * x1 + diag[17] * x2 + diag[23] * x3 + diag[29] * x4 + diag[35] * x5; 1610e1b4bd6SMark F. Adams diag += 36; 1620e1b4bd6SMark F. Adams } 1639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 1659566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(66.0 * m)); 1660e1b4bd6SMark F. Adams PetscFunctionReturn(0); 1670e1b4bd6SMark F. Adams } 1689371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_7(PC pc, Vec x, Vec y) { 1690a4ca348SMatthew G Knepley PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 1700a4ca348SMatthew G Knepley PetscInt i, m = jac->mbs; 1710a4ca348SMatthew G Knepley const MatScalar *diag = jac->diag; 1720e666e78SToby Isaac PetscScalar x0, x1, x2, x3, x4, x5, x6, *yy; 1730e666e78SToby Isaac const PetscScalar *xx; 1740a4ca348SMatthew G Knepley 1750a4ca348SMatthew G Knepley PetscFunctionBegin; 1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1779566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 1780a4ca348SMatthew G Knepley for (i = 0; i < m; i++) { 1799371c9d4SSatish Balay x0 = xx[7 * i]; 1809371c9d4SSatish Balay x1 = xx[7 * i + 1]; 1819371c9d4SSatish Balay x2 = xx[7 * i + 2]; 1829371c9d4SSatish Balay x3 = xx[7 * i + 3]; 1839371c9d4SSatish Balay x4 = xx[7 * i + 4]; 1849371c9d4SSatish Balay x5 = xx[7 * i + 5]; 1859371c9d4SSatish Balay x6 = xx[7 * i + 6]; 1862fa5cd67SKarl Rupp 1870a4ca348SMatthew G Knepley yy[7 * i] = diag[0] * x0 + diag[7] * x1 + diag[14] * x2 + diag[21] * x3 + diag[28] * x4 + diag[35] * x5 + diag[42] * x6; 1880a4ca348SMatthew G Knepley yy[7 * i + 1] = diag[1] * x0 + diag[8] * x1 + diag[15] * x2 + diag[22] * x3 + diag[29] * x4 + diag[36] * x5 + diag[43] * x6; 1890a4ca348SMatthew G Knepley yy[7 * i + 2] = diag[2] * x0 + diag[9] * x1 + diag[16] * x2 + diag[23] * x3 + diag[30] * x4 + diag[37] * x5 + diag[44] * x6; 1900a4ca348SMatthew G Knepley yy[7 * i + 3] = diag[3] * x0 + diag[10] * x1 + diag[17] * x2 + diag[24] * x3 + diag[31] * x4 + diag[38] * x5 + diag[45] * x6; 1910a4ca348SMatthew G Knepley yy[7 * i + 4] = diag[4] * x0 + diag[11] * x1 + diag[18] * x2 + diag[25] * x3 + diag[32] * x4 + diag[39] * x5 + diag[46] * x6; 1920a4ca348SMatthew G Knepley yy[7 * i + 5] = diag[5] * x0 + diag[12] * x1 + diag[19] * x2 + diag[26] * x3 + diag[33] * x4 + diag[40] * x5 + diag[47] * x6; 1930a4ca348SMatthew G Knepley yy[7 * i + 6] = diag[6] * x0 + diag[13] * x1 + diag[20] * x2 + diag[27] * x3 + diag[34] * x4 + diag[41] * x5 + diag[48] * x6; 1940a4ca348SMatthew G Knepley diag += 49; 1950a4ca348SMatthew G Knepley } 1969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 1989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(91.0 * m)); /* 2*bs2 - bs */ 1990a4ca348SMatthew G Knepley PetscFunctionReturn(0); 2000a4ca348SMatthew G Knepley } 2019371c9d4SSatish Balay static PetscErrorCode PCApply_PBJacobi_N(PC pc, Vec x, Vec y) { 202f79daf70SMark Lohry PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 203f79daf70SMark Lohry PetscInt i, ib, jb; 204f79daf70SMark Lohry const PetscInt m = jac->mbs; 205f79daf70SMark Lohry const PetscInt bs = jac->bs; 206f79daf70SMark Lohry const MatScalar *diag = jac->diag; 207f79daf70SMark Lohry PetscScalar *yy; 208f79daf70SMark Lohry const PetscScalar *xx; 209f79daf70SMark Lohry 210f79daf70SMark Lohry PetscFunctionBegin; 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 213f79daf70SMark Lohry for (i = 0; i < m; i++) { 214f79daf70SMark Lohry for (ib = 0; ib < bs; ib++) { 215f79daf70SMark Lohry PetscScalar rowsum = 0; 216ad540459SPierre Jolivet for (jb = 0; jb < bs; jb++) rowsum += diag[ib + jb * bs] * xx[bs * i + jb]; 217f79daf70SMark Lohry yy[bs * i + ib] = rowsum; 218f79daf70SMark Lohry } 219f79daf70SMark Lohry diag += bs * bs; 220f79daf70SMark Lohry } 2219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 2239566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((2.0 * bs * bs - bs) * m)); /* 2*bs2 - bs */ 224f79daf70SMark Lohry PetscFunctionReturn(0); 225f79daf70SMark Lohry } 226faff7d23SJed Brown 2279371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_PBJacobi_N(PC pc, Vec x, Vec y) { 2288cd12025SJed Brown PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 2298cd12025SJed Brown PetscInt i, j, k, m = jac->mbs, bs = jac->bs; 2308cd12025SJed Brown const MatScalar *diag = jac->diag; 2318cd12025SJed Brown const PetscScalar *xx; 2328cd12025SJed Brown PetscScalar *yy; 2338cd12025SJed Brown 2348cd12025SJed Brown PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 2378cd12025SJed Brown for (i = 0; i < m; i++) { 2388cd12025SJed Brown for (j = 0; j < bs; j++) yy[i * bs + j] = 0.; 2398cd12025SJed Brown for (j = 0; j < bs; j++) { 240ad540459SPierre Jolivet for (k = 0; k < bs; k++) yy[i * bs + k] += diag[k * bs + j] * xx[i * bs + j]; 2418cd12025SJed Brown } 2428cd12025SJed Brown diag += bs * bs; 2438cd12025SJed Brown } 2449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 2469566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(m * bs * (2 * bs - 1))); 2478cd12025SJed Brown PetscFunctionReturn(0); 2488cd12025SJed Brown } 2498cd12025SJed Brown 2509371c9d4SSatish Balay static PetscErrorCode PCSetUp_PBJacobi(PC pc) { 2514b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 2524b9ad928SBarry Smith Mat A = pc->pmat; 25300e125f8SBarry Smith MatFactorError err; 254539c167fSBarry Smith PetscInt nlocal; 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith PetscFunctionBegin; 2579566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(A, &jac->diag)); 2589566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(A, &err)); 25900e125f8SBarry Smith if (err) pc->failedreason = (PCFailedReason)err; 260d7806518SHong Zhang 2619566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &jac->bs)); 2629566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &nlocal, NULL)); 263539c167fSBarry Smith jac->mbs = nlocal / jac->bs; 264521d7252SBarry Smith switch (jac->bs) { 2659371c9d4SSatish Balay case 1: pc->ops->apply = PCApply_PBJacobi_1; break; 2669371c9d4SSatish Balay case 2: pc->ops->apply = PCApply_PBJacobi_2; break; 2679371c9d4SSatish Balay case 3: pc->ops->apply = PCApply_PBJacobi_3; break; 2689371c9d4SSatish Balay case 4: pc->ops->apply = PCApply_PBJacobi_4; break; 2699371c9d4SSatish Balay case 5: pc->ops->apply = PCApply_PBJacobi_5; break; 2709371c9d4SSatish Balay case 6: pc->ops->apply = PCApply_PBJacobi_6; break; 2719371c9d4SSatish Balay case 7: pc->ops->apply = PCApply_PBJacobi_7; break; 2729371c9d4SSatish Balay default: pc->ops->apply = PCApply_PBJacobi_N; break; 2734b9ad928SBarry Smith } 2748cd12025SJed Brown pc->ops->applytranspose = PCApplyTranspose_PBJacobi_N; 2754b9ad928SBarry Smith PetscFunctionReturn(0); 2764b9ad928SBarry Smith } 277*f1580f4eSBarry Smith 2789371c9d4SSatish Balay static PetscErrorCode PCDestroy_PBJacobi(PC pc) { 2794b9ad928SBarry Smith PetscFunctionBegin; 2804b9ad928SBarry Smith /* 2814b9ad928SBarry Smith Free the private data structure that was hanging off the PC 2824b9ad928SBarry Smith */ 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2844b9ad928SBarry Smith PetscFunctionReturn(0); 2854b9ad928SBarry Smith } 286fa939db7SJed Brown 2879371c9d4SSatish Balay static PetscErrorCode PCView_PBJacobi(PC pc, PetscViewer viewer) { 288fa939db7SJed Brown PC_PBJacobi *jac = (PC_PBJacobi *)pc->data; 289fa939db7SJed Brown PetscBool iascii; 290fa939db7SJed Brown 291fa939db7SJed Brown PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 29348a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " point-block size %" PetscInt_FMT "\n", jac->bs)); 294fa939db7SJed Brown PetscFunctionReturn(0); 295fa939db7SJed Brown } 296fa939db7SJed Brown 29737a17b4dSBarry Smith /*MC 298422a814eSBarry Smith PCPBJACOBI - Point block Jacobi preconditioner 299422a814eSBarry Smith 30095452b02SPatrick Sanan Notes: 301*f1580f4eSBarry Smith See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable-size point block, and `PCBJACOBI` for large size blocks 302422a814eSBarry Smith 303*f1580f4eSBarry Smith This works for `MATAIJ` and `MATBAIJ` matrices and uses the blocksize provided to the matrix 304422a814eSBarry Smith 305422a814eSBarry Smith Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot 306422a814eSBarry Smith is detected a PETSc error is generated. 307422a814eSBarry Smith 30895452b02SPatrick Sanan Developer Notes: 309*f1580f4eSBarry Smith This should support the `PCSetErrorIfFailure()` flag set to `PETSC_TRUE` to allow 310422a814eSBarry Smith the factorization to continue even after a zero pivot is found resulting in a Nan and hence 311*f1580f4eSBarry Smith terminating `KSP` with a `KSP_DIVERGED_NANORIF` allowing 312422a814eSBarry Smith a nonlinear solver/ODE integrator to recover without stopping the program as currently happens. 313422a814eSBarry Smith 31490dfcc32SBarry Smith Perhaps should provide an option that allows generation of a valid preconditioner 315*f1580f4eSBarry Smith even if a block is singular as the `PCJACOBI` does. 31637a17b4dSBarry Smith 31737a17b4dSBarry Smith Level: beginner 31837a17b4dSBarry Smith 319db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCVPBJACOBI`, `PCBJACOBI` 32037a17b4dSBarry Smith M*/ 32137a17b4dSBarry Smith 3229371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc) { 3234b9ad928SBarry Smith PC_PBJacobi *jac; 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith PetscFunctionBegin; 3264b9ad928SBarry Smith /* 3274b9ad928SBarry Smith Creates the private data structure for this preconditioner and 3284b9ad928SBarry Smith attach it to the PC object. 3294b9ad928SBarry Smith */ 3309566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &jac)); 3314b9ad928SBarry Smith pc->data = (void *)jac; 3324b9ad928SBarry Smith 3334b9ad928SBarry Smith /* 3344b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 3354b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 3364b9ad928SBarry Smith */ 3370a545947SLisandro Dalcin jac->diag = NULL; 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith /* 3404b9ad928SBarry Smith Set the pointers for the functions that are provided above. 3414b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 3424b9ad928SBarry Smith are called, they will automatically call these functions. Note we 3434b9ad928SBarry Smith choose not to provide a couple of these functions since they are 3444b9ad928SBarry Smith not needed. 3454b9ad928SBarry Smith */ 3460a545947SLisandro Dalcin pc->ops->apply = NULL; /*set depending on the block size */ 3470a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 3484b9ad928SBarry Smith pc->ops->setup = PCSetUp_PBJacobi; 3494b9ad928SBarry Smith pc->ops->destroy = PCDestroy_PBJacobi; 3500a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 351fa939db7SJed Brown pc->ops->view = PCView_PBJacobi; 3520a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3530a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 3540a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 3554b9ad928SBarry Smith PetscFunctionReturn(0); 3564b9ad928SBarry Smith } 357