1f1be3500SJunchao Zhang #include <../src/ksp/pc/impls/vpbjacobi/vpbjacobi.h> 20da83c2eSBarry Smith 3d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_VPBJacobi(PC pc, Vec x, Vec y) 4d71ae5a4SJacob Faibussowitsch { 50da83c2eSBarry Smith PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data; 60da83c2eSBarry Smith PetscInt i, ncnt = 0; 70da83c2eSBarry Smith const MatScalar *diag = jac->diag; 80da83c2eSBarry Smith PetscInt ib, jb, bs; 90da83c2eSBarry Smith const PetscScalar *xx; 100da83c2eSBarry Smith PetscScalar *yy, x0, x1, x2, x3, x4, x5, x6; 110da83c2eSBarry Smith PetscInt nblocks; 120da83c2eSBarry Smith const PetscInt *bsizes; 130da83c2eSBarry Smith 140da83c2eSBarry Smith PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes)); 169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 179566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &yy)); 180da83c2eSBarry Smith for (i = 0; i < nblocks; i++) { 190da83c2eSBarry Smith bs = bsizes[i]; 200da83c2eSBarry Smith switch (bs) { 21d71ae5a4SJacob Faibussowitsch case 1: 22d71ae5a4SJacob Faibussowitsch yy[ncnt] = *diag * xx[ncnt]; 23d71ae5a4SJacob Faibussowitsch break; 240da83c2eSBarry Smith case 2: 259371c9d4SSatish Balay x0 = xx[ncnt]; 269371c9d4SSatish Balay x1 = xx[ncnt + 1]; 270da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[2] * x1; 280da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[3] * x1; 290da83c2eSBarry Smith break; 300da83c2eSBarry Smith case 3: 319371c9d4SSatish Balay x0 = xx[ncnt]; 329371c9d4SSatish Balay x1 = xx[ncnt + 1]; 339371c9d4SSatish Balay x2 = xx[ncnt + 2]; 340da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[3] * x1 + diag[6] * x2; 350da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[4] * x1 + diag[7] * x2; 360da83c2eSBarry Smith yy[ncnt + 2] = diag[2] * x0 + diag[5] * x1 + diag[8] * x2; 370da83c2eSBarry Smith break; 380da83c2eSBarry Smith case 4: 399371c9d4SSatish Balay x0 = xx[ncnt]; 409371c9d4SSatish Balay x1 = xx[ncnt + 1]; 419371c9d4SSatish Balay x2 = xx[ncnt + 2]; 429371c9d4SSatish Balay x3 = xx[ncnt + 3]; 430da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[4] * x1 + diag[8] * x2 + diag[12] * x3; 440da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[5] * x1 + diag[9] * x2 + diag[13] * x3; 450da83c2eSBarry Smith yy[ncnt + 2] = diag[2] * x0 + diag[6] * x1 + diag[10] * x2 + diag[14] * x3; 460da83c2eSBarry Smith yy[ncnt + 3] = diag[3] * x0 + diag[7] * x1 + diag[11] * x2 + diag[15] * x3; 470da83c2eSBarry Smith break; 480da83c2eSBarry Smith case 5: 499371c9d4SSatish Balay x0 = xx[ncnt]; 509371c9d4SSatish Balay x1 = xx[ncnt + 1]; 519371c9d4SSatish Balay x2 = xx[ncnt + 2]; 529371c9d4SSatish Balay x3 = xx[ncnt + 3]; 539371c9d4SSatish Balay x4 = xx[ncnt + 4]; 540da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[5] * x1 + diag[10] * x2 + diag[15] * x3 + diag[20] * x4; 550da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[6] * x1 + diag[11] * x2 + diag[16] * x3 + diag[21] * x4; 560da83c2eSBarry Smith yy[ncnt + 2] = diag[2] * x0 + diag[7] * x1 + diag[12] * x2 + diag[17] * x3 + diag[22] * x4; 570da83c2eSBarry Smith yy[ncnt + 3] = diag[3] * x0 + diag[8] * x1 + diag[13] * x2 + diag[18] * x3 + diag[23] * x4; 580da83c2eSBarry Smith yy[ncnt + 4] = diag[4] * x0 + diag[9] * x1 + diag[14] * x2 + diag[19] * x3 + diag[24] * x4; 590da83c2eSBarry Smith break; 600da83c2eSBarry Smith case 6: 619371c9d4SSatish Balay x0 = xx[ncnt]; 629371c9d4SSatish Balay x1 = xx[ncnt + 1]; 639371c9d4SSatish Balay x2 = xx[ncnt + 2]; 649371c9d4SSatish Balay x3 = xx[ncnt + 3]; 659371c9d4SSatish Balay x4 = xx[ncnt + 4]; 669371c9d4SSatish Balay x5 = xx[ncnt + 5]; 670da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[6] * x1 + diag[12] * x2 + diag[18] * x3 + diag[24] * x4 + diag[30] * x5; 680da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[7] * x1 + diag[13] * x2 + diag[19] * x3 + diag[25] * x4 + diag[31] * x5; 690da83c2eSBarry Smith yy[ncnt + 2] = diag[2] * x0 + diag[8] * x1 + diag[14] * x2 + diag[20] * x3 + diag[26] * x4 + diag[32] * x5; 700da83c2eSBarry Smith yy[ncnt + 3] = diag[3] * x0 + diag[9] * x1 + diag[15] * x2 + diag[21] * x3 + diag[27] * x4 + diag[33] * x5; 710da83c2eSBarry Smith yy[ncnt + 4] = diag[4] * x0 + diag[10] * x1 + diag[16] * x2 + diag[22] * x3 + diag[28] * x4 + diag[34] * x5; 720da83c2eSBarry Smith yy[ncnt + 5] = diag[5] * x0 + diag[11] * x1 + diag[17] * x2 + diag[23] * x3 + diag[29] * x4 + diag[35] * x5; 730da83c2eSBarry Smith break; 740da83c2eSBarry Smith case 7: 759371c9d4SSatish Balay x0 = xx[ncnt]; 769371c9d4SSatish Balay x1 = xx[ncnt + 1]; 779371c9d4SSatish Balay x2 = xx[ncnt + 2]; 789371c9d4SSatish Balay x3 = xx[ncnt + 3]; 799371c9d4SSatish Balay x4 = xx[ncnt + 4]; 809371c9d4SSatish Balay x5 = xx[ncnt + 5]; 819371c9d4SSatish Balay x6 = xx[ncnt + 6]; 820da83c2eSBarry Smith yy[ncnt] = diag[0] * x0 + diag[7] * x1 + diag[14] * x2 + diag[21] * x3 + diag[28] * x4 + diag[35] * x5 + diag[42] * x6; 830da83c2eSBarry Smith yy[ncnt + 1] = diag[1] * x0 + diag[8] * x1 + diag[15] * x2 + diag[22] * x3 + diag[29] * x4 + diag[36] * x5 + diag[43] * x6; 840da83c2eSBarry Smith yy[ncnt + 2] = diag[2] * x0 + diag[9] * x1 + diag[16] * x2 + diag[23] * x3 + diag[30] * x4 + diag[37] * x5 + diag[44] * x6; 850da83c2eSBarry Smith yy[ncnt + 3] = diag[3] * x0 + diag[10] * x1 + diag[17] * x2 + diag[24] * x3 + diag[31] * x4 + diag[38] * x5 + diag[45] * x6; 860da83c2eSBarry Smith yy[ncnt + 4] = diag[4] * x0 + diag[11] * x1 + diag[18] * x2 + diag[25] * x3 + diag[32] * x4 + diag[39] * x5 + diag[46] * x6; 870da83c2eSBarry Smith yy[ncnt + 5] = diag[5] * x0 + diag[12] * x1 + diag[19] * x2 + diag[26] * x3 + diag[33] * x4 + diag[40] * x5 + diag[47] * x6; 880da83c2eSBarry Smith yy[ncnt + 6] = diag[6] * x0 + diag[13] * x1 + diag[20] * x2 + diag[27] * x3 + diag[34] * x4 + diag[41] * x5 + diag[48] * x6; 890da83c2eSBarry Smith break; 900da83c2eSBarry Smith default: 910da83c2eSBarry Smith for (ib = 0; ib < bs; ib++) { 920da83c2eSBarry Smith PetscScalar rowsum = 0; 93ad540459SPierre Jolivet for (jb = 0; jb < bs; jb++) rowsum += diag[ib + jb * bs] * xx[ncnt + jb]; 940da83c2eSBarry Smith yy[ncnt + ib] = rowsum; 950da83c2eSBarry Smith } 960da83c2eSBarry Smith } 970da83c2eSBarry Smith ncnt += bsizes[i]; 980da83c2eSBarry Smith diag += bsizes[i] * bsizes[i]; 990da83c2eSBarry Smith } 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &yy)); 1020da83c2eSBarry Smith PetscFunctionReturn(0); 1030da83c2eSBarry Smith } 1040da83c2eSBarry Smith 105*69eda9daSJed Brown static PetscErrorCode PCApplyTranspose_VPBJacobi(PC pc, Vec x, Vec y) 106*69eda9daSJed Brown { 107*69eda9daSJed Brown PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data; 108*69eda9daSJed Brown PetscInt i, ncnt = 0; 109*69eda9daSJed Brown const MatScalar *diag = jac->diag; 110*69eda9daSJed Brown PetscInt ib, jb, bs; 111*69eda9daSJed Brown const PetscScalar *xx; 112*69eda9daSJed Brown PetscScalar *yy, x0, x1, x2, x3, x4, x5, x6; 113*69eda9daSJed Brown PetscInt nblocks; 114*69eda9daSJed Brown const PetscInt *bsizes; 115*69eda9daSJed Brown 116*69eda9daSJed Brown PetscFunctionBegin; 117*69eda9daSJed Brown PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes)); 118*69eda9daSJed Brown PetscCall(VecGetArrayRead(x, &xx)); 119*69eda9daSJed Brown PetscCall(VecGetArray(y, &yy)); 120*69eda9daSJed Brown for (i = 0; i < nblocks; i++) { 121*69eda9daSJed Brown bs = bsizes[i]; 122*69eda9daSJed Brown switch (bs) { 123*69eda9daSJed Brown case 1: 124*69eda9daSJed Brown yy[ncnt] = *diag * xx[ncnt]; 125*69eda9daSJed Brown break; 126*69eda9daSJed Brown case 2: 127*69eda9daSJed Brown x0 = xx[ncnt]; 128*69eda9daSJed Brown x1 = xx[ncnt + 1]; 129*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1; 130*69eda9daSJed Brown yy[ncnt + 1] = diag[2] * x0 + diag[3] * x1; 131*69eda9daSJed Brown break; 132*69eda9daSJed Brown case 3: 133*69eda9daSJed Brown x0 = xx[ncnt]; 134*69eda9daSJed Brown x1 = xx[ncnt + 1]; 135*69eda9daSJed Brown x2 = xx[ncnt + 2]; 136*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1 + diag[2] * x2; 137*69eda9daSJed Brown yy[ncnt + 1] = diag[3] * x0 + diag[4] * x1 + diag[5] * x2; 138*69eda9daSJed Brown yy[ncnt + 2] = diag[6] * x0 + diag[7] * x1 + diag[8] * x2; 139*69eda9daSJed Brown break; 140*69eda9daSJed Brown case 4: 141*69eda9daSJed Brown x0 = xx[ncnt]; 142*69eda9daSJed Brown x1 = xx[ncnt + 1]; 143*69eda9daSJed Brown x2 = xx[ncnt + 2]; 144*69eda9daSJed Brown x3 = xx[ncnt + 3]; 145*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1 + diag[2] * x2 + diag[3] * x3; 146*69eda9daSJed Brown yy[ncnt + 1] = diag[4] * x0 + diag[5] * x1 + diag[6] * x2 + diag[7] * x3; 147*69eda9daSJed Brown yy[ncnt + 2] = diag[8] * x0 + diag[9] * x1 + diag[10] * x2 + diag[11] * x3; 148*69eda9daSJed Brown yy[ncnt + 3] = diag[12] * x0 + diag[13] * x1 + diag[14] * x2 + diag[15] * x3; 149*69eda9daSJed Brown break; 150*69eda9daSJed Brown case 5: 151*69eda9daSJed Brown x0 = xx[ncnt]; 152*69eda9daSJed Brown x1 = xx[ncnt + 1]; 153*69eda9daSJed Brown x2 = xx[ncnt + 2]; 154*69eda9daSJed Brown x3 = xx[ncnt + 3]; 155*69eda9daSJed Brown x4 = xx[ncnt + 4]; 156*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1 + diag[2] * x2 + diag[3] * x3 + diag[4] * x4; 157*69eda9daSJed Brown yy[ncnt + 1] = diag[5] * x0 + diag[6] * x1 + diag[7] * x2 + diag[8] * x3 + diag[9] * x4; 158*69eda9daSJed Brown yy[ncnt + 2] = diag[10] * x0 + diag[11] * x1 + diag[12] * x2 + diag[13] * x3 + diag[14] * x4; 159*69eda9daSJed Brown yy[ncnt + 3] = diag[15] * x0 + diag[16] * x1 + diag[17] * x2 + diag[18] * x3 + diag[19] * x4; 160*69eda9daSJed Brown yy[ncnt + 4] = diag[20] * x0 + diag[21] * x1 + diag[22] * x2 + diag[23] * x3 + diag[24] * x4; 161*69eda9daSJed Brown break; 162*69eda9daSJed Brown case 6: 163*69eda9daSJed Brown x0 = xx[ncnt]; 164*69eda9daSJed Brown x1 = xx[ncnt + 1]; 165*69eda9daSJed Brown x2 = xx[ncnt + 2]; 166*69eda9daSJed Brown x3 = xx[ncnt + 3]; 167*69eda9daSJed Brown x4 = xx[ncnt + 4]; 168*69eda9daSJed Brown x5 = xx[ncnt + 5]; 169*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1 + diag[2] * x2 + diag[3] * x3 + diag[4] * x4 + diag[5] * x5; 170*69eda9daSJed Brown yy[ncnt + 1] = diag[6] * x0 + diag[7] * x1 + diag[8] * x2 + diag[9] * x3 + diag[10] * x4 + diag[11] * x5; 171*69eda9daSJed Brown yy[ncnt + 2] = diag[12] * x0 + diag[13] * x1 + diag[14] * x2 + diag[15] * x3 + diag[16] * x4 + diag[17] * x5; 172*69eda9daSJed Brown yy[ncnt + 3] = diag[18] * x0 + diag[19] * x1 + diag[20] * x2 + diag[21] * x3 + diag[22] * x4 + diag[23] * x5; 173*69eda9daSJed Brown yy[ncnt + 4] = diag[24] * x0 + diag[25] * x1 + diag[26] * x2 + diag[27] * x3 + diag[28] * x4 + diag[29] * x5; 174*69eda9daSJed Brown yy[ncnt + 5] = diag[30] * x0 + diag[31] * x1 + diag[32] * x2 + diag[33] * x3 + diag[34] * x4 + diag[35] * x5; 175*69eda9daSJed Brown break; 176*69eda9daSJed Brown case 7: 177*69eda9daSJed Brown x0 = xx[ncnt]; 178*69eda9daSJed Brown x1 = xx[ncnt + 1]; 179*69eda9daSJed Brown x2 = xx[ncnt + 2]; 180*69eda9daSJed Brown x3 = xx[ncnt + 3]; 181*69eda9daSJed Brown x4 = xx[ncnt + 4]; 182*69eda9daSJed Brown x5 = xx[ncnt + 5]; 183*69eda9daSJed Brown x6 = xx[ncnt + 6]; 184*69eda9daSJed Brown yy[ncnt] = diag[0] * x0 + diag[1] * x1 + diag[2] * x2 + diag[3] * x3 + diag[4] * x4 + diag[5] * x5 + diag[6] * x6; 185*69eda9daSJed Brown yy[ncnt + 1] = diag[7] * x0 + diag[8] * x1 + diag[9] * x2 + diag[10] * x3 + diag[11] * x4 + diag[12] * x5 + diag[13] * x6; 186*69eda9daSJed Brown yy[ncnt + 2] = diag[14] * x0 + diag[15] * x1 + diag[16] * x2 + diag[17] * x3 + diag[18] * x4 + diag[19] * x5 + diag[20] * x6; 187*69eda9daSJed Brown yy[ncnt + 3] = diag[21] * x0 + diag[22] * x1 + diag[23] * x2 + diag[24] * x3 + diag[25] * x4 + diag[26] * x5 + diag[27] * x6; 188*69eda9daSJed Brown yy[ncnt + 4] = diag[28] * x0 + diag[29] * x1 + diag[30] * x2 + diag[31] * x3 + diag[32] * x4 + diag[33] * x5 + diag[34] * x6; 189*69eda9daSJed Brown yy[ncnt + 5] = diag[35] * x0 + diag[36] * x1 + diag[37] * x2 + diag[38] * x3 + diag[39] * x4 + diag[40] * x5 + diag[41] * x6; 190*69eda9daSJed Brown yy[ncnt + 6] = diag[42] * x0 + diag[43] * x1 + diag[44] * x2 + diag[45] * x3 + diag[46] * x4 + diag[47] * x5 + diag[48] * x6; 191*69eda9daSJed Brown break; 192*69eda9daSJed Brown default: 193*69eda9daSJed Brown for (ib = 0; ib < bs; ib++) { 194*69eda9daSJed Brown PetscScalar rowsum = 0; 195*69eda9daSJed Brown for (jb = 0; jb < bs; jb++) rowsum += diag[ib * bs + jb] * xx[ncnt + jb]; 196*69eda9daSJed Brown yy[ncnt + ib] = rowsum; 197*69eda9daSJed Brown } 198*69eda9daSJed Brown } 199*69eda9daSJed Brown ncnt += bsizes[i]; 200*69eda9daSJed Brown diag += bsizes[i] * bsizes[i]; 201*69eda9daSJed Brown } 202*69eda9daSJed Brown PetscCall(VecRestoreArrayRead(x, &xx)); 203*69eda9daSJed Brown PetscCall(VecRestoreArray(y, &yy)); 204*69eda9daSJed Brown PetscFunctionReturn(0); 205*69eda9daSJed Brown } 206*69eda9daSJed Brown 207d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCSetUp_VPBJacobi_Host(PC pc) 208d71ae5a4SJacob Faibussowitsch { 2090da83c2eSBarry Smith PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data; 2100da83c2eSBarry Smith Mat A = pc->pmat; 2110da83c2eSBarry Smith MatFactorError err; 2120da83c2eSBarry Smith PetscInt i, nsize = 0, nlocal; 2130da83c2eSBarry Smith PetscInt nblocks; 2140da83c2eSBarry Smith const PetscInt *bsizes; 2150da83c2eSBarry Smith 2160da83c2eSBarry Smith PetscFunctionBegin; 2179566063dSJacob Faibussowitsch PetscCall(MatGetVariableBlockSizes(pc->pmat, &nblocks, &bsizes)); 2189566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &nlocal, NULL)); 21908401ef6SPierre Jolivet PetscCheck(!nlocal || nblocks, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call MatSetVariableBlockSizes() before using PCVPBJACOBI"); 2200da83c2eSBarry Smith if (!jac->diag) { 2210da83c2eSBarry Smith for (i = 0; i < nblocks; i++) nsize += bsizes[i] * bsizes[i]; 2229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nsize, &jac->diag)); 2230da83c2eSBarry Smith } 2249566063dSJacob Faibussowitsch PetscCall(MatInvertVariableBlockDiagonal(A, nblocks, bsizes, jac->diag)); 2259566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(A, &err)); 2260da83c2eSBarry Smith if (err) pc->failedreason = (PCFailedReason)err; 2270da83c2eSBarry Smith pc->ops->apply = PCApply_VPBJacobi; 228*69eda9daSJed Brown pc->ops->applytranspose = PCApplyTranspose_VPBJacobi; 2290da83c2eSBarry Smith PetscFunctionReturn(0); 2300da83c2eSBarry Smith } 2318a9c020eSBarry Smith 232d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_VPBJacobi(PC pc) 233d71ae5a4SJacob Faibussowitsch { 234f1be3500SJunchao Zhang PetscFunctionBegin; 235f1be3500SJunchao Zhang /* In PCCreate_VPBJacobi() pmat might have not been set, so we wait to the last minute to do the dispatch */ 236f1be3500SJunchao Zhang #if defined(PETSC_HAVE_CUDA) 237f1be3500SJunchao Zhang PetscBool isCuda; 238f1be3500SJunchao Zhang PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isCuda, MATSEQAIJCUSPARSE, MATMPIAIJCUSPARSE, "")); 239f1be3500SJunchao Zhang #endif 240f1be3500SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS_KERNELS) 241f1be3500SJunchao Zhang PetscBool isKok; 242f1be3500SJunchao Zhang PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isKok, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, "")); 243f1be3500SJunchao Zhang #endif 244f1be3500SJunchao Zhang 245f1be3500SJunchao Zhang #if defined(PETSC_HAVE_CUDA) 2469371c9d4SSatish Balay if (isCuda) PetscCall(PCSetUp_VPBJacobi_CUDA(pc)); 2479371c9d4SSatish Balay else 248f1be3500SJunchao Zhang #endif 249f1be3500SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS_KERNELS) 2509371c9d4SSatish Balay if (isKok) 2519371c9d4SSatish Balay PetscCall(PCSetUp_VPBJacobi_Kokkos(pc)); 2529371c9d4SSatish Balay else 253f1be3500SJunchao Zhang #endif 2549371c9d4SSatish Balay { 2559371c9d4SSatish Balay PetscCall(PCSetUp_VPBJacobi_Host(pc)); 2569371c9d4SSatish Balay } 257f1be3500SJunchao Zhang PetscFunctionReturn(0); 258f1be3500SJunchao Zhang } 259f1be3500SJunchao Zhang 260d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCDestroy_VPBJacobi(PC pc) 261d71ae5a4SJacob Faibussowitsch { 2620da83c2eSBarry Smith PC_VPBJacobi *jac = (PC_VPBJacobi *)pc->data; 2630da83c2eSBarry Smith 2640da83c2eSBarry Smith PetscFunctionBegin; 2650da83c2eSBarry Smith /* 2660da83c2eSBarry Smith Free the private data structure that was hanging off the PC 2670da83c2eSBarry Smith */ 2689566063dSJacob Faibussowitsch PetscCall(PetscFree(jac->diag)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2700da83c2eSBarry Smith PetscFunctionReturn(0); 2710da83c2eSBarry Smith } 2720da83c2eSBarry Smith 2730da83c2eSBarry Smith /*MC 2740da83c2eSBarry Smith PCVPBJACOBI - Variable size point block Jacobi preconditioner 2750da83c2eSBarry Smith 276f1580f4eSBarry Smith Level: beginner 2770da83c2eSBarry Smith 278f1580f4eSBarry Smith Notes: 279f1580f4eSBarry Smith See `PCJACOBI` for point Jacobi preconditioning, `PCPBJACOBI` for fixed point block size, and `PCBJACOBI` for large size blocks 280f1580f4eSBarry Smith 281f1580f4eSBarry Smith This works for `MATAIJ` matrices 2820da83c2eSBarry Smith 2830da83c2eSBarry Smith Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot 2840da83c2eSBarry Smith is detected a PETSc error is generated. 2850da83c2eSBarry Smith 286f1580f4eSBarry Smith One must call `MatSetVariableBlockSizes()` to use this preconditioner 28790dfcc32SBarry Smith 2880da83c2eSBarry Smith Developer Notes: 289f1580f4eSBarry Smith This should support the `PCSetErrorIfFailure()` flag set to `PETSC_TRUE` to allow 2900da83c2eSBarry Smith the factorization to continue even after a zero pivot is found resulting in a Nan and hence 291f1580f4eSBarry Smith terminating `KSP` with a `KSP_DIVERGED_NANORIF` allowing 2920da83c2eSBarry Smith a nonlinear solver/ODE integrator to recover without stopping the program as currently happens. 2930da83c2eSBarry Smith 2940da83c2eSBarry Smith Perhaps should provide an option that allows generation of a valid preconditioner 295f1580f4eSBarry Smith even if a block is singular as the `PCJACOBI` does. 2960da83c2eSBarry Smith 297db781477SPatrick Sanan .seealso: `MatSetVariableBlockSizes()`, `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, `PCPBJACOBI`, `PCBJACOBI` 2980da83c2eSBarry Smith M*/ 2990da83c2eSBarry Smith 300d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc) 301d71ae5a4SJacob Faibussowitsch { 3020da83c2eSBarry Smith PC_VPBJacobi *jac; 3030da83c2eSBarry Smith 3040da83c2eSBarry Smith PetscFunctionBegin; 3050da83c2eSBarry Smith /* 3060da83c2eSBarry Smith Creates the private data structure for this preconditioner and 3070da83c2eSBarry Smith attach it to the PC object. 3080da83c2eSBarry Smith */ 3094dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 3100da83c2eSBarry Smith pc->data = (void *)jac; 3110da83c2eSBarry Smith 3120da83c2eSBarry Smith /* 3130da83c2eSBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 3140da83c2eSBarry Smith diagonal entries of the matrix for fast preconditioner application. 3150da83c2eSBarry Smith */ 3160da83c2eSBarry Smith jac->diag = NULL; 3170da83c2eSBarry Smith 3180da83c2eSBarry Smith /* 3190da83c2eSBarry Smith Set the pointers for the functions that are provided above. 3200da83c2eSBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 3210da83c2eSBarry Smith are called, they will automatically call these functions. Note we 3220da83c2eSBarry Smith choose not to provide a couple of these functions since they are 3230da83c2eSBarry Smith not needed. 3240da83c2eSBarry Smith */ 3250da83c2eSBarry Smith pc->ops->apply = PCApply_VPBJacobi; 3260a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 3270da83c2eSBarry Smith pc->ops->setup = PCSetUp_VPBJacobi; 3280da83c2eSBarry Smith pc->ops->destroy = PCDestroy_VPBJacobi; 3290a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 3300a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3310a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 3320a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 3330da83c2eSBarry Smith PetscFunctionReturn(0); 3340da83c2eSBarry Smith } 335