xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision 539c167f0a31d1b968f995eba09d310cf6f1900b)
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 
174b9ad928SBarry Smith 
184b9ad928SBarry Smith #undef __FUNCT__
19bbead8a2SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_1"
20bbead8a2SBarry Smith static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y)
21bbead8a2SBarry Smith {
22bbead8a2SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi*)pc->data;
23bbead8a2SBarry Smith   PetscErrorCode    ierr;
24bbead8a2SBarry Smith   PetscInt          i,m = jac->mbs;
25bbead8a2SBarry Smith   const MatScalar   *diag = jac->diag;
26bbead8a2SBarry Smith   const PetscScalar *xx;
27bbead8a2SBarry Smith   PetscScalar       *yy;
28bbead8a2SBarry Smith 
29bbead8a2SBarry Smith   PetscFunctionBegin;
30bbead8a2SBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
31bbead8a2SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
322fa5cd67SKarl Rupp   for (i=0; i<m; i++) yy[i] = diag[i]*xx[i];
33bbead8a2SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
34bbead8a2SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
35bbead8a2SBarry Smith   ierr = PetscLogFlops(2.0*m);CHKERRQ(ierr);
36bbead8a2SBarry Smith   PetscFunctionReturn(0);
37bbead8a2SBarry Smith }
38bbead8a2SBarry Smith 
39bbead8a2SBarry Smith #undef __FUNCT__
404b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
416849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
424b9ad928SBarry Smith {
434b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
44dfbe8321SBarry Smith   PetscErrorCode  ierr;
45c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
4685f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
470e666e78SToby Isaac   PetscScalar     x0,x1,*yy;
480e666e78SToby Isaac   const PetscScalar *xx;
494b9ad928SBarry Smith 
504b9ad928SBarry Smith   PetscFunctionBegin;
510e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
524b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
534b9ad928SBarry Smith   for (i=0; i<m; i++) {
544b9ad928SBarry Smith     x0        = xx[2*i]; x1 = xx[2*i+1];
554b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
564b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
574b9ad928SBarry Smith     diag     += 4;
584b9ad928SBarry Smith   }
590e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
604b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
61dc0b31edSSatish Balay   ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr);
624b9ad928SBarry Smith   PetscFunctionReturn(0);
634b9ad928SBarry Smith }
644b9ad928SBarry Smith #undef __FUNCT__
654b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
666849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
674b9ad928SBarry Smith {
684b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
69dfbe8321SBarry Smith   PetscErrorCode  ierr;
70c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
7185f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
720e666e78SToby Isaac   PetscScalar     x0,x1,x2,*yy;
730e666e78SToby Isaac   const PetscScalar *xx;
744b9ad928SBarry Smith 
754b9ad928SBarry Smith   PetscFunctionBegin;
760e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
774b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
784b9ad928SBarry Smith   for (i=0; i<m; i++) {
794b9ad928SBarry Smith     x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
802fa5cd67SKarl Rupp 
814b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
824b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
834b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
844b9ad928SBarry Smith     diag     += 9;
854b9ad928SBarry Smith   }
860e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
874b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
88dc0b31edSSatish Balay   ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr);
894b9ad928SBarry Smith   PetscFunctionReturn(0);
904b9ad928SBarry Smith }
914b9ad928SBarry Smith #undef __FUNCT__
924b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4"
936849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
944b9ad928SBarry Smith {
954b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
96dfbe8321SBarry Smith   PetscErrorCode  ierr;
97c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
9885f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
990e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,*yy;
1000e666e78SToby Isaac   const PetscScalar *xx;
1014b9ad928SBarry Smith 
1024b9ad928SBarry Smith   PetscFunctionBegin;
1030e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1044b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1054b9ad928SBarry Smith   for (i=0; i<m; i++) {
1064b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
1072fa5cd67SKarl Rupp 
1084b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
1094b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
1104b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
1114b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
1124b9ad928SBarry Smith     diag     += 16;
1134b9ad928SBarry Smith   }
1140e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1154b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
116dc0b31edSSatish Balay   ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr);
1174b9ad928SBarry Smith   PetscFunctionReturn(0);
1184b9ad928SBarry Smith }
1194b9ad928SBarry Smith #undef __FUNCT__
1204b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
1216849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1224b9ad928SBarry Smith {
1234b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
124dfbe8321SBarry Smith   PetscErrorCode  ierr;
125c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
12685f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
1270e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,*yy;
1280e666e78SToby Isaac   const PetscScalar *xx;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
1310e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1324b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1334b9ad928SBarry Smith   for (i=0; i<m; i++) {
1344b9ad928SBarry 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];
1352fa5cd67SKarl Rupp 
1364b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1374b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1384b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1394b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1404b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1414b9ad928SBarry Smith     diag     += 25;
1424b9ad928SBarry Smith   }
1430e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1444b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
145dc0b31edSSatish Balay   ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr);
1464b9ad928SBarry Smith   PetscFunctionReturn(0);
1474b9ad928SBarry Smith }
1480e1b4bd6SMark F. Adams #undef __FUNCT__
1490e1b4bd6SMark F. Adams #define __FUNCT__ "PCApply_PBJacobi_6"
1500e1b4bd6SMark F. Adams static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y)
1510e1b4bd6SMark F. Adams {
1520e1b4bd6SMark F. Adams   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1530e1b4bd6SMark F. Adams   PetscErrorCode  ierr;
1540e1b4bd6SMark F. Adams   PetscInt        i,m = jac->mbs;
1550e1b4bd6SMark F. Adams   const MatScalar *diag = jac->diag;
1560e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,x5,*yy;
1570e666e78SToby Isaac   const PetscScalar *xx;
1580e1b4bd6SMark F. Adams 
1590e1b4bd6SMark F. Adams   PetscFunctionBegin;
1600e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1610e1b4bd6SMark F. Adams   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1620e1b4bd6SMark F. Adams   for (i=0; i<m; i++) {
1630e1b4bd6SMark 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];
1642fa5cd67SKarl Rupp 
1650e1b4bd6SMark F. Adams     yy[6*i]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
1660e1b4bd6SMark F. Adams     yy[6*i+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
1670e1b4bd6SMark F. Adams     yy[6*i+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
1680e1b4bd6SMark F. Adams     yy[6*i+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
1690e1b4bd6SMark F. Adams     yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
1700e1b4bd6SMark F. Adams     yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
1710e1b4bd6SMark F. Adams     diag     += 36;
1720e1b4bd6SMark F. Adams   }
1730e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1740e1b4bd6SMark F. Adams   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1750e1b4bd6SMark F. Adams   ierr = PetscLogFlops(66.0*m);CHKERRQ(ierr);
1760e1b4bd6SMark F. Adams   PetscFunctionReturn(0);
1770e1b4bd6SMark F. Adams }
1780a4ca348SMatthew G Knepley #undef __FUNCT__
1790a4ca348SMatthew G Knepley #define __FUNCT__ "PCApply_PBJacobi_7"
1800a4ca348SMatthew G Knepley static PetscErrorCode PCApply_PBJacobi_7(PC pc,Vec x,Vec y)
1810a4ca348SMatthew G Knepley {
1820a4ca348SMatthew G Knepley   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1830a4ca348SMatthew G Knepley   PetscErrorCode  ierr;
1840a4ca348SMatthew G Knepley   PetscInt        i,m = jac->mbs;
1850a4ca348SMatthew G Knepley   const MatScalar *diag = jac->diag;
1860e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,x5,x6,*yy;
1870e666e78SToby Isaac   const PetscScalar *xx;
1880a4ca348SMatthew G Knepley 
1890a4ca348SMatthew G Knepley   PetscFunctionBegin;
1900e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1910a4ca348SMatthew G Knepley   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1920a4ca348SMatthew G Knepley   for (i=0; i<m; i++) {
1930a4ca348SMatthew 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];
1942fa5cd67SKarl Rupp 
1950a4ca348SMatthew 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;
1960a4ca348SMatthew 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;
1970a4ca348SMatthew 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;
1980a4ca348SMatthew 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;
1990a4ca348SMatthew 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;
2000a4ca348SMatthew 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;
2010a4ca348SMatthew 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;
2020a4ca348SMatthew G Knepley     diag     += 49;
2030a4ca348SMatthew G Knepley   }
2040e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
2050a4ca348SMatthew G Knepley   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
2060a4ca348SMatthew G Knepley   ierr = PetscLogFlops(80.0*m);CHKERRQ(ierr);
2070a4ca348SMatthew G Knepley   PetscFunctionReturn(0);
2080a4ca348SMatthew G Knepley }
2094b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2104b9ad928SBarry Smith #undef __FUNCT__
2114b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
2126849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
2134b9ad928SBarry Smith {
2144b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
215dfbe8321SBarry Smith   PetscErrorCode ierr;
2164b9ad928SBarry Smith   Mat            A = pc->pmat;
21700e125f8SBarry Smith   MatFactorError err;
218*539c167fSBarry Smith   PetscInt       nlocal;
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith   PetscFunctionBegin;
221bbead8a2SBarry Smith   ierr = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr);
22200e125f8SBarry Smith   ierr = MatFactorGetError(A,&err);CHKERRQ(ierr);
22300e125f8SBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
224d7806518SHong Zhang 
22533d57670SJed Brown   ierr = MatGetBlockSize(A,&jac->bs);CHKERRQ(ierr);
226*539c167fSBarry Smith   ierr = MatGetLocalSize(A,&nlocal,NULL);CHKERRQ(ierr);
227*539c167fSBarry Smith   jac->mbs = nlocal/jac->bs;
228521d7252SBarry Smith   switch (jac->bs) {
229bbead8a2SBarry Smith   case 1:
230bbead8a2SBarry Smith     pc->ops->apply = PCApply_PBJacobi_1;
231bbead8a2SBarry Smith     break;
2324b9ad928SBarry Smith   case 2:
2334b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_2;
2344b9ad928SBarry Smith     break;
2354b9ad928SBarry Smith   case 3:
2364b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_3;
2374b9ad928SBarry Smith     break;
2384b9ad928SBarry Smith   case 4:
2394b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_4;
2404b9ad928SBarry Smith     break;
2414b9ad928SBarry Smith   case 5:
2424b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_5;
2434b9ad928SBarry Smith     break;
2440e1b4bd6SMark F. Adams   case 6:
2450e1b4bd6SMark F. Adams     pc->ops->apply = PCApply_PBJacobi_6;
2460e1b4bd6SMark F. Adams     break;
2470a4ca348SMatthew G Knepley   case 7:
2480a4ca348SMatthew G Knepley     pc->ops->apply = PCApply_PBJacobi_7;
2490a4ca348SMatthew G Knepley     break;
2504b9ad928SBarry Smith   default:
251ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
2524b9ad928SBarry Smith   }
2534b9ad928SBarry Smith   PetscFunctionReturn(0);
2544b9ad928SBarry Smith }
2554b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2564b9ad928SBarry Smith #undef __FUNCT__
2574b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
2586849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
2594b9ad928SBarry Smith {
260dfbe8321SBarry Smith   PetscErrorCode ierr;
2614b9ad928SBarry Smith 
2624b9ad928SBarry Smith   PetscFunctionBegin;
2634b9ad928SBarry Smith   /*
2644b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
2654b9ad928SBarry Smith   */
266c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2674b9ad928SBarry Smith   PetscFunctionReturn(0);
2684b9ad928SBarry Smith }
269fa939db7SJed Brown 
270fa939db7SJed Brown #undef __FUNCT__
271fa939db7SJed Brown #define __FUNCT__ "PCView_PBJacobi"
272fa939db7SJed Brown static PetscErrorCode PCView_PBJacobi(PC pc,PetscViewer viewer)
273fa939db7SJed Brown {
274fa939db7SJed Brown   PetscErrorCode ierr;
275fa939db7SJed Brown   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
276fa939db7SJed Brown   PetscBool      iascii;
277fa939db7SJed Brown 
278fa939db7SJed Brown   PetscFunctionBegin;
279251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
280fa939db7SJed Brown   if (iascii) {
281fa939db7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  point-block Jacobi: block size %D\n",jac->bs);CHKERRQ(ierr);
28211aeaf0aSBarry Smith   }
283fa939db7SJed Brown   PetscFunctionReturn(0);
284fa939db7SJed Brown }
285fa939db7SJed Brown 
2864b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
28737a17b4dSBarry Smith /*MC
288422a814eSBarry Smith      PCPBJACOBI - Point block Jacobi preconditioner
289422a814eSBarry Smith 
290422a814eSBarry Smith 
291422a814eSBarry Smith    Notes: See PCJACOBI for point Jacobi preconditioning
292422a814eSBarry Smith 
293422a814eSBarry Smith    This works for AIJ and BAIJ matrices and uses the blocksize provided to the matrix
294422a814eSBarry Smith 
295422a814eSBarry Smith    Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
296422a814eSBarry Smith    is detected a PETSc error is generated.
297422a814eSBarry Smith 
298422a814eSBarry Smith    Developer Notes: This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow
299422a814eSBarry Smith    the factorization to continue even after a zero pivot is found resulting in a Nan and hence
300422a814eSBarry Smith    terminating KSP with a KSP_DIVERGED_NANORIF allowing
301422a814eSBarry Smith    a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
302422a814eSBarry Smith 
303422a814eSBarry Smith    Developer Note: Perhaps should provide an option that allows generation of a valid preconditioner
304422a814eSBarry Smith    even if a block is singular as the PCJACOBI does.
30537a17b4dSBarry Smith 
30637a17b4dSBarry Smith    Level: beginner
30737a17b4dSBarry Smith 
30837a17b4dSBarry Smith   Concepts: point block Jacobi
30937a17b4dSBarry Smith 
31037a17b4dSBarry Smith 
311422a814eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI
31237a17b4dSBarry Smith 
31337a17b4dSBarry Smith M*/
31437a17b4dSBarry Smith 
3154b9ad928SBarry Smith #undef __FUNCT__
3164b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
3178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc)
3184b9ad928SBarry Smith {
3194b9ad928SBarry Smith   PC_PBJacobi    *jac;
320dfbe8321SBarry Smith   PetscErrorCode ierr;
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith   PetscFunctionBegin;
3234b9ad928SBarry Smith   /*
3244b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
3254b9ad928SBarry Smith      attach it to the PC object.
3264b9ad928SBarry Smith   */
327b00a9115SJed Brown   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
3284b9ad928SBarry Smith   pc->data = (void*)jac;
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith   /*
3314b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
3324b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
3334b9ad928SBarry Smith   */
3344b9ad928SBarry Smith   jac->diag = 0;
3354b9ad928SBarry Smith 
3364b9ad928SBarry Smith   /*
3374b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
3384b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
3394b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
3404b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
3414b9ad928SBarry Smith       not needed.
3424b9ad928SBarry Smith   */
3434b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
3444b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
3454b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
3464b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
3474b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
348fa939db7SJed Brown   pc->ops->view                = PCView_PBJacobi;
3494b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
3504b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
3514b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
3524b9ad928SBarry Smith   PetscFunctionReturn(0);
3534b9ad928SBarry Smith }
3544b9ad928SBarry Smith 
3554b9ad928SBarry Smith 
356