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 17bbead8a2SBarry Smith static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y) 18bbead8a2SBarry Smith { 19bbead8a2SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 20bbead8a2SBarry Smith PetscInt i,m = jac->mbs; 21bbead8a2SBarry Smith const MatScalar *diag = jac->diag; 22bbead8a2SBarry Smith const PetscScalar *xx; 23bbead8a2SBarry Smith PetscScalar *yy; 24bbead8a2SBarry Smith 25bbead8a2SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 279566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 282fa5cd67SKarl Rupp for (i=0; i<m; i++) yy[i] = diag[i]*xx[i]; 299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 319566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(m)); 32bbead8a2SBarry Smith PetscFunctionReturn(0); 33bbead8a2SBarry Smith } 34bbead8a2SBarry Smith 356849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y) 364b9ad928SBarry Smith { 374b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 38c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 3985f4f44aSBarry Smith const MatScalar *diag = jac->diag; 400e666e78SToby Isaac PetscScalar x0,x1,*yy; 410e666e78SToby Isaac const PetscScalar *xx; 424b9ad928SBarry Smith 434b9ad928SBarry Smith PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 459566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 464b9ad928SBarry Smith for (i=0; i<m; i++) { 474b9ad928SBarry Smith x0 = xx[2*i]; x1 = xx[2*i+1]; 484b9ad928SBarry Smith yy[2*i] = diag[0]*x0 + diag[2]*x1; 494b9ad928SBarry Smith yy[2*i+1] = diag[1]*x0 + diag[3]*x1; 504b9ad928SBarry Smith diag += 4; 514b9ad928SBarry Smith } 529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 549566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(6.0*m)); 554b9ad928SBarry Smith PetscFunctionReturn(0); 564b9ad928SBarry Smith } 576849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y) 584b9ad928SBarry Smith { 594b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 60c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 6185f4f44aSBarry Smith const MatScalar *diag = jac->diag; 620e666e78SToby Isaac PetscScalar x0,x1,x2,*yy; 630e666e78SToby Isaac const PetscScalar *xx; 644b9ad928SBarry Smith 654b9ad928SBarry Smith PetscFunctionBegin; 669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 679566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 684b9ad928SBarry Smith for (i=0; i<m; i++) { 694b9ad928SBarry Smith x0 = xx[3*i]; x1 = xx[3*i+1]; 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 } 816849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y) 824b9ad928SBarry Smith { 834b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 84c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 8585f4f44aSBarry Smith const MatScalar *diag = jac->diag; 860e666e78SToby Isaac PetscScalar x0,x1,x2,x3,*yy; 870e666e78SToby Isaac const PetscScalar *xx; 884b9ad928SBarry Smith 894b9ad928SBarry Smith PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 919566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 924b9ad928SBarry Smith for (i=0; i<m; i++) { 934b9ad928SBarry Smith x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3]; 942fa5cd67SKarl Rupp 954b9ad928SBarry Smith yy[4*i] = diag[0]*x0 + diag[4]*x1 + diag[8]*x2 + diag[12]*x3; 964b9ad928SBarry Smith yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2 + diag[13]*x3; 974b9ad928SBarry Smith yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3; 984b9ad928SBarry Smith yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3; 994b9ad928SBarry Smith diag += 16; 1004b9ad928SBarry Smith } 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 1039566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(28.0*m)); 1044b9ad928SBarry Smith PetscFunctionReturn(0); 1054b9ad928SBarry Smith } 1066849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y) 1074b9ad928SBarry Smith { 1084b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 109c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 11085f4f44aSBarry Smith const MatScalar *diag = jac->diag; 1110e666e78SToby Isaac PetscScalar x0,x1,x2,x3,x4,*yy; 1120e666e78SToby Isaac const PetscScalar *xx; 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1169566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 1174b9ad928SBarry Smith for (i=0; i<m; i++) { 1184b9ad928SBarry Smith x0 = xx[5*i]; x1 = xx[5*i+1]; x2 = xx[5*i+2]; x3 = xx[5*i+3]; x4 = xx[5*i+4]; 1192fa5cd67SKarl Rupp 1204b9ad928SBarry Smith yy[5*i] = diag[0]*x0 + diag[5]*x1 + diag[10]*x2 + diag[15]*x3 + diag[20]*x4; 1214b9ad928SBarry Smith yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2 + diag[16]*x3 + diag[21]*x4; 1224b9ad928SBarry Smith yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4; 1234b9ad928SBarry Smith yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4; 1244b9ad928SBarry Smith yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4; 1254b9ad928SBarry Smith diag += 25; 1264b9ad928SBarry Smith } 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 1299566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(45.0*m)); 1304b9ad928SBarry Smith PetscFunctionReturn(0); 1314b9ad928SBarry Smith } 1320e1b4bd6SMark F. Adams static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y) 1330e1b4bd6SMark F. Adams { 1340e1b4bd6SMark F. Adams PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 1350e1b4bd6SMark F. Adams PetscInt i,m = jac->mbs; 1360e1b4bd6SMark F. Adams const MatScalar *diag = jac->diag; 1370e666e78SToby Isaac PetscScalar x0,x1,x2,x3,x4,x5,*yy; 1380e666e78SToby Isaac const PetscScalar *xx; 1390e1b4bd6SMark F. Adams 1400e1b4bd6SMark F. Adams PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1429566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 1430e1b4bd6SMark F. Adams for (i=0; i<m; i++) { 1440e1b4bd6SMark F. Adams x0 = xx[6*i]; x1 = xx[6*i+1]; x2 = xx[6*i+2]; x3 = xx[6*i+3]; x4 = xx[6*i+4]; x5 = xx[6*i+5]; 1452fa5cd67SKarl Rupp 1460e1b4bd6SMark F. Adams yy[6*i] = diag[0]*x0 + diag[6]*x1 + diag[12]*x2 + diag[18]*x3 + diag[24]*x4 + diag[30]*x5; 1470e1b4bd6SMark F. Adams yy[6*i+1] = diag[1]*x0 + diag[7]*x1 + diag[13]*x2 + diag[19]*x3 + diag[25]*x4 + diag[31]*x5; 1480e1b4bd6SMark F. Adams yy[6*i+2] = diag[2]*x0 + diag[8]*x1 + diag[14]*x2 + diag[20]*x3 + diag[26]*x4 + diag[32]*x5; 1490e1b4bd6SMark F. Adams yy[6*i+3] = diag[3]*x0 + diag[9]*x1 + diag[15]*x2 + diag[21]*x3 + diag[27]*x4 + diag[33]*x5; 1500e1b4bd6SMark F. Adams yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2 + diag[22]*x3 + diag[28]*x4 + diag[34]*x5; 1510e1b4bd6SMark F. Adams yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2 + diag[23]*x3 + diag[29]*x4 + diag[35]*x5; 1520e1b4bd6SMark F. Adams diag += 36; 1530e1b4bd6SMark F. Adams } 1549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 1569566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(66.0*m)); 1570e1b4bd6SMark F. Adams PetscFunctionReturn(0); 1580e1b4bd6SMark F. Adams } 1590a4ca348SMatthew G Knepley static PetscErrorCode PCApply_PBJacobi_7(PC pc,Vec x,Vec y) 1600a4ca348SMatthew G Knepley { 1610a4ca348SMatthew G Knepley PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 1620a4ca348SMatthew G Knepley PetscInt i,m = jac->mbs; 1630a4ca348SMatthew G Knepley const MatScalar *diag = jac->diag; 1640e666e78SToby Isaac PetscScalar x0,x1,x2,x3,x4,x5,x6,*yy; 1650e666e78SToby Isaac const PetscScalar *xx; 1660a4ca348SMatthew G Knepley 1670a4ca348SMatthew G Knepley PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1699566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 1700a4ca348SMatthew G Knepley for (i=0; i<m; i++) { 1710a4ca348SMatthew G Knepley x0 = xx[7*i]; x1 = xx[7*i+1]; x2 = xx[7*i+2]; x3 = xx[7*i+3]; x4 = xx[7*i+4]; x5 = xx[7*i+5]; x6 = xx[7*i+6]; 1722fa5cd67SKarl Rupp 1730a4ca348SMatthew 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; 1740a4ca348SMatthew 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; 1750a4ca348SMatthew 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; 1760a4ca348SMatthew 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; 1770a4ca348SMatthew 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; 1780a4ca348SMatthew 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; 1790a4ca348SMatthew 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; 1800a4ca348SMatthew G Knepley diag += 49; 1810a4ca348SMatthew G Knepley } 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 1849566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(91.0*m)); /* 2*bs2 - bs */ 1850a4ca348SMatthew G Knepley PetscFunctionReturn(0); 1860a4ca348SMatthew G Knepley } 187f79daf70SMark Lohry static PetscErrorCode PCApply_PBJacobi_N(PC pc,Vec x,Vec y) 188f79daf70SMark Lohry { 189f79daf70SMark Lohry PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 190f79daf70SMark Lohry PetscInt i,ib,jb; 191f79daf70SMark Lohry const PetscInt m = jac->mbs; 192f79daf70SMark Lohry const PetscInt bs = jac->bs; 193f79daf70SMark Lohry const MatScalar *diag = jac->diag; 194f79daf70SMark Lohry PetscScalar *yy; 195f79daf70SMark Lohry const PetscScalar *xx; 196f79daf70SMark Lohry 197f79daf70SMark Lohry PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1999566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 200f79daf70SMark Lohry for (i=0; i<m; i++) { 201f79daf70SMark Lohry for (ib=0; ib<bs; ib++) { 202f79daf70SMark Lohry PetscScalar rowsum = 0; 203f79daf70SMark Lohry for (jb=0; jb<bs; jb++) { 204f79daf70SMark Lohry rowsum += diag[ib+jb*bs] * xx[bs*i+jb]; 205f79daf70SMark Lohry } 206f79daf70SMark Lohry yy[bs*i+ib] = rowsum; 207f79daf70SMark Lohry } 208f79daf70SMark Lohry diag += bs*bs; 209f79daf70SMark Lohry } 2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 2129566063dSJacob Faibussowitsch PetscCall(PetscLogFlops((2.0*bs*bs-bs)*m)); /* 2*bs2 - bs */ 213f79daf70SMark Lohry PetscFunctionReturn(0); 214f79daf70SMark Lohry } 215faff7d23SJed Brown 2168cd12025SJed Brown static PetscErrorCode PCApplyTranspose_PBJacobi_N(PC pc,Vec x,Vec y) 2178cd12025SJed Brown { 2188cd12025SJed Brown PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 2198cd12025SJed Brown PetscInt i,j,k,m = jac->mbs,bs=jac->bs; 2208cd12025SJed Brown const MatScalar *diag = jac->diag; 2218cd12025SJed Brown const PetscScalar *xx; 2228cd12025SJed Brown PetscScalar *yy; 2238cd12025SJed Brown 2248cd12025SJed Brown PetscFunctionBegin; 2259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&yy)); 2278cd12025SJed Brown for (i=0; i<m; i++) { 2288cd12025SJed Brown for (j=0; j<bs; j++) yy[i*bs+j] = 0.; 2298cd12025SJed Brown for (j=0; j<bs; j++) { 2308cd12025SJed Brown for (k=0; k<bs; k++) { 2318cd12025SJed Brown yy[i*bs+k] += diag[k*bs+j]*xx[i*bs+j]; 2328cd12025SJed Brown } 2338cd12025SJed Brown } 2348cd12025SJed Brown diag += bs*bs; 2358cd12025SJed Brown } 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&yy)); 2389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(m*bs*(2*bs-1))); 2398cd12025SJed Brown PetscFunctionReturn(0); 2408cd12025SJed Brown } 2418cd12025SJed Brown 2424b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2436849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc) 2444b9ad928SBarry Smith { 2454b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 2464b9ad928SBarry Smith Mat A = pc->pmat; 24700e125f8SBarry Smith MatFactorError err; 248539c167fSBarry Smith PetscInt nlocal; 2494b9ad928SBarry Smith 2504b9ad928SBarry Smith PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(A,&jac->diag)); 2529566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(A,&err)); 25300e125f8SBarry Smith if (err) pc->failedreason = (PCFailedReason)err; 254d7806518SHong Zhang 2559566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&jac->bs)); 2569566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&nlocal,NULL)); 257539c167fSBarry Smith jac->mbs = nlocal/jac->bs; 258521d7252SBarry Smith switch (jac->bs) { 259bbead8a2SBarry Smith case 1: 260bbead8a2SBarry Smith pc->ops->apply = PCApply_PBJacobi_1; 261bbead8a2SBarry Smith break; 2624b9ad928SBarry Smith case 2: 2634b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_2; 2644b9ad928SBarry Smith break; 2654b9ad928SBarry Smith case 3: 2664b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_3; 2674b9ad928SBarry Smith break; 2684b9ad928SBarry Smith case 4: 2694b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_4; 2704b9ad928SBarry Smith break; 2714b9ad928SBarry Smith case 5: 2724b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_5; 2734b9ad928SBarry Smith break; 2740e1b4bd6SMark F. Adams case 6: 2750e1b4bd6SMark F. Adams pc->ops->apply = PCApply_PBJacobi_6; 2760e1b4bd6SMark F. Adams break; 2770a4ca348SMatthew G Knepley case 7: 2780a4ca348SMatthew G Knepley pc->ops->apply = PCApply_PBJacobi_7; 2790a4ca348SMatthew G Knepley break; 2804b9ad928SBarry Smith default: 281f79daf70SMark Lohry pc->ops->apply = PCApply_PBJacobi_N; 282f79daf70SMark Lohry break; 2834b9ad928SBarry Smith } 2848cd12025SJed Brown pc->ops->applytranspose = PCApplyTranspose_PBJacobi_N; 2854b9ad928SBarry Smith PetscFunctionReturn(0); 2864b9ad928SBarry Smith } 2874b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 2886849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc) 2894b9ad928SBarry Smith { 2904b9ad928SBarry Smith PetscFunctionBegin; 2914b9ad928SBarry Smith /* 2924b9ad928SBarry Smith Free the private data structure that was hanging off the PC 2934b9ad928SBarry Smith */ 2949566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 297fa939db7SJed Brown 298fa939db7SJed Brown static PetscErrorCode PCView_PBJacobi(PC pc,PetscViewer viewer) 299fa939db7SJed Brown { 300fa939db7SJed Brown PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 301fa939db7SJed Brown PetscBool iascii; 302fa939db7SJed Brown 303fa939db7SJed Brown PetscFunctionBegin; 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 305fa939db7SJed Brown if (iascii) { 30663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," point-block size %" PetscInt_FMT "\n",jac->bs)); 30711aeaf0aSBarry Smith } 308fa939db7SJed Brown PetscFunctionReturn(0); 309fa939db7SJed Brown } 310fa939db7SJed Brown 3114b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 31237a17b4dSBarry Smith /*MC 313422a814eSBarry Smith PCPBJACOBI - Point block Jacobi preconditioner 314422a814eSBarry Smith 31595452b02SPatrick Sanan Notes: 316*468ae2e8SBarry Smith See PCJACOBI for diagonal Jacobi, PCVPBJACOBI for variable-size point block, and PCBJACOBI for large size blocks 317422a814eSBarry Smith 318422a814eSBarry Smith This works for AIJ and BAIJ matrices and uses the blocksize provided to the matrix 319422a814eSBarry Smith 320422a814eSBarry Smith Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot 321422a814eSBarry Smith is detected a PETSc error is generated. 322422a814eSBarry Smith 32395452b02SPatrick Sanan Developer Notes: 32495452b02SPatrick Sanan This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow 325422a814eSBarry Smith the factorization to continue even after a zero pivot is found resulting in a Nan and hence 326422a814eSBarry Smith terminating KSP with a KSP_DIVERGED_NANORIF allowing 327422a814eSBarry Smith a nonlinear solver/ODE integrator to recover without stopping the program as currently happens. 328422a814eSBarry Smith 32990dfcc32SBarry Smith Perhaps should provide an option that allows generation of a valid preconditioner 330422a814eSBarry Smith even if a block is singular as the PCJACOBI does. 33137a17b4dSBarry Smith 33237a17b4dSBarry Smith Level: beginner 33337a17b4dSBarry Smith 33490dfcc32SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI, PCVPBJACOBI, PCBJACOBI 33537a17b4dSBarry Smith 33637a17b4dSBarry Smith M*/ 33737a17b4dSBarry Smith 3388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc) 3394b9ad928SBarry Smith { 3404b9ad928SBarry Smith PC_PBJacobi *jac; 3414b9ad928SBarry Smith 3424b9ad928SBarry Smith PetscFunctionBegin; 3434b9ad928SBarry Smith /* 3444b9ad928SBarry Smith Creates the private data structure for this preconditioner and 3454b9ad928SBarry Smith attach it to the PC object. 3464b9ad928SBarry Smith */ 3479566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&jac)); 3484b9ad928SBarry Smith pc->data = (void*)jac; 3494b9ad928SBarry Smith 3504b9ad928SBarry Smith /* 3514b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 3524b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 3534b9ad928SBarry Smith */ 3540a545947SLisandro Dalcin jac->diag = NULL; 3554b9ad928SBarry Smith 3564b9ad928SBarry Smith /* 3574b9ad928SBarry Smith Set the pointers for the functions that are provided above. 3584b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 3594b9ad928SBarry Smith are called, they will automatically call these functions. Note we 3604b9ad928SBarry Smith choose not to provide a couple of these functions since they are 3614b9ad928SBarry Smith not needed. 3624b9ad928SBarry Smith */ 3630a545947SLisandro Dalcin pc->ops->apply = NULL; /*set depending on the block size */ 3640a545947SLisandro Dalcin pc->ops->applytranspose = NULL; 3654b9ad928SBarry Smith pc->ops->setup = PCSetUp_PBJacobi; 3664b9ad928SBarry Smith pc->ops->destroy = PCDestroy_PBJacobi; 3670a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 368fa939db7SJed Brown pc->ops->view = PCView_PBJacobi; 3690a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3700a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 3710a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 3724b9ad928SBarry Smith PetscFunctionReturn(0); 3734b9ad928SBarry Smith } 374