xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
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 
18bbead8a2SBarry Smith static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y)
19bbead8a2SBarry Smith {
20bbead8a2SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi*)pc->data;
21bbead8a2SBarry Smith   PetscErrorCode    ierr;
22bbead8a2SBarry Smith   PetscInt          i,m = jac->mbs;
23bbead8a2SBarry Smith   const MatScalar   *diag = jac->diag;
24bbead8a2SBarry Smith   const PetscScalar *xx;
25bbead8a2SBarry Smith   PetscScalar       *yy;
26bbead8a2SBarry Smith 
27bbead8a2SBarry Smith   PetscFunctionBegin;
28bbead8a2SBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
29bbead8a2SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
302fa5cd67SKarl Rupp   for (i=0; i<m; i++) yy[i] = diag[i]*xx[i];
31bbead8a2SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
32bbead8a2SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
3379cd6d07SHong Zhang   ierr = PetscLogFlops(m);CHKERRQ(ierr);
34bbead8a2SBarry Smith   PetscFunctionReturn(0);
35bbead8a2SBarry Smith }
36bbead8a2SBarry Smith 
376849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
384b9ad928SBarry Smith {
394b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
40dfbe8321SBarry Smith   PetscErrorCode  ierr;
41c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
4285f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
430e666e78SToby Isaac   PetscScalar     x0,x1,*yy;
440e666e78SToby Isaac   const PetscScalar *xx;
454b9ad928SBarry Smith 
464b9ad928SBarry Smith   PetscFunctionBegin;
470e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
484b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
494b9ad928SBarry Smith   for (i=0; i<m; i++) {
504b9ad928SBarry Smith     x0        = xx[2*i]; x1 = xx[2*i+1];
514b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
524b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
534b9ad928SBarry Smith     diag     += 4;
544b9ad928SBarry Smith   }
550e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
564b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
57dc0b31edSSatish Balay   ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr);
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
606849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
614b9ad928SBarry Smith {
624b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
63dfbe8321SBarry Smith   PetscErrorCode  ierr;
64c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
6585f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
660e666e78SToby Isaac   PetscScalar     x0,x1,x2,*yy;
670e666e78SToby Isaac   const PetscScalar *xx;
684b9ad928SBarry Smith 
694b9ad928SBarry Smith   PetscFunctionBegin;
700e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
714b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
724b9ad928SBarry Smith   for (i=0; i<m; i++) {
734b9ad928SBarry Smith     x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
742fa5cd67SKarl Rupp 
754b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
764b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
774b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
784b9ad928SBarry Smith     diag     += 9;
794b9ad928SBarry Smith   }
800e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
814b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
82dc0b31edSSatish Balay   ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr);
834b9ad928SBarry Smith   PetscFunctionReturn(0);
844b9ad928SBarry Smith }
856849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
864b9ad928SBarry Smith {
874b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
88dfbe8321SBarry Smith   PetscErrorCode  ierr;
89c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
9085f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
910e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,*yy;
920e666e78SToby Isaac   const PetscScalar *xx;
934b9ad928SBarry Smith 
944b9ad928SBarry Smith   PetscFunctionBegin;
950e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
964b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
974b9ad928SBarry Smith   for (i=0; i<m; i++) {
984b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
992fa5cd67SKarl Rupp 
1004b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
1014b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
1024b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
1034b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
1044b9ad928SBarry Smith     diag     += 16;
1054b9ad928SBarry Smith   }
1060e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1074b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
108dc0b31edSSatish Balay   ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr);
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1116849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1124b9ad928SBarry Smith {
1134b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
114dfbe8321SBarry Smith   PetscErrorCode  ierr;
115c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
11685f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
1170e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,*yy;
1180e666e78SToby Isaac   const PetscScalar *xx;
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith   PetscFunctionBegin;
1210e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1224b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1234b9ad928SBarry Smith   for (i=0; i<m; i++) {
1244b9ad928SBarry 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];
1252fa5cd67SKarl Rupp 
1264b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1274b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1284b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1294b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1304b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1314b9ad928SBarry Smith     diag     += 25;
1324b9ad928SBarry Smith   }
1330e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1344b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
135dc0b31edSSatish Balay   ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr);
1364b9ad928SBarry Smith   PetscFunctionReturn(0);
1374b9ad928SBarry Smith }
1380e1b4bd6SMark F. Adams static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y)
1390e1b4bd6SMark F. Adams {
1400e1b4bd6SMark F. Adams   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1410e1b4bd6SMark F. Adams   PetscErrorCode  ierr;
1420e1b4bd6SMark F. Adams   PetscInt        i,m = jac->mbs;
1430e1b4bd6SMark F. Adams   const MatScalar *diag = jac->diag;
1440e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,x5,*yy;
1450e666e78SToby Isaac   const PetscScalar *xx;
1460e1b4bd6SMark F. Adams 
1470e1b4bd6SMark F. Adams   PetscFunctionBegin;
1480e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1490e1b4bd6SMark F. Adams   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1500e1b4bd6SMark F. Adams   for (i=0; i<m; i++) {
1510e1b4bd6SMark 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];
1522fa5cd67SKarl Rupp 
1530e1b4bd6SMark F. Adams     yy[6*i]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
1540e1b4bd6SMark F. Adams     yy[6*i+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
1550e1b4bd6SMark F. Adams     yy[6*i+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
1560e1b4bd6SMark F. Adams     yy[6*i+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
1570e1b4bd6SMark F. Adams     yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
1580e1b4bd6SMark F. Adams     yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
1590e1b4bd6SMark F. Adams     diag     += 36;
1600e1b4bd6SMark F. Adams   }
1610e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1620e1b4bd6SMark F. Adams   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1630e1b4bd6SMark F. Adams   ierr = PetscLogFlops(66.0*m);CHKERRQ(ierr);
1640e1b4bd6SMark F. Adams   PetscFunctionReturn(0);
1650e1b4bd6SMark F. Adams }
1660a4ca348SMatthew G Knepley static PetscErrorCode PCApply_PBJacobi_7(PC pc,Vec x,Vec y)
1670a4ca348SMatthew G Knepley {
1680a4ca348SMatthew G Knepley   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1690a4ca348SMatthew G Knepley   PetscErrorCode  ierr;
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;
1760e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1770a4ca348SMatthew G Knepley   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1780a4ca348SMatthew G Knepley   for (i=0; i<m; i++) {
1790a4ca348SMatthew 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];
1802fa5cd67SKarl Rupp 
1810a4ca348SMatthew 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;
1820a4ca348SMatthew 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;
1830a4ca348SMatthew 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;
1840a4ca348SMatthew 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;
1850a4ca348SMatthew 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;
1860a4ca348SMatthew 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;
1870a4ca348SMatthew 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;
1880a4ca348SMatthew G Knepley     diag     += 49;
1890a4ca348SMatthew G Knepley   }
1900e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1910a4ca348SMatthew G Knepley   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
19279cd6d07SHong Zhang   ierr = PetscLogFlops(91*m);CHKERRQ(ierr); /* 2*bs2 - bs */
1930a4ca348SMatthew G Knepley   PetscFunctionReturn(0);
1940a4ca348SMatthew G Knepley }
1954b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1966849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
1974b9ad928SBarry Smith {
1984b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
199dfbe8321SBarry Smith   PetscErrorCode ierr;
2004b9ad928SBarry Smith   Mat            A = pc->pmat;
20100e125f8SBarry Smith   MatFactorError err;
202539c167fSBarry Smith   PetscInt       nlocal;
2034b9ad928SBarry Smith 
2044b9ad928SBarry Smith   PetscFunctionBegin;
205bbead8a2SBarry Smith   ierr = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr);
20600e125f8SBarry Smith   ierr = MatFactorGetError(A,&err);CHKERRQ(ierr);
20700e125f8SBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
208d7806518SHong Zhang 
20933d57670SJed Brown   ierr = MatGetBlockSize(A,&jac->bs);CHKERRQ(ierr);
210539c167fSBarry Smith   ierr = MatGetLocalSize(A,&nlocal,NULL);CHKERRQ(ierr);
211539c167fSBarry Smith   jac->mbs = nlocal/jac->bs;
212521d7252SBarry Smith   switch (jac->bs) {
213bbead8a2SBarry Smith   case 1:
214bbead8a2SBarry Smith     pc->ops->apply = PCApply_PBJacobi_1;
215bbead8a2SBarry Smith     break;
2164b9ad928SBarry Smith   case 2:
2174b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_2;
2184b9ad928SBarry Smith     break;
2194b9ad928SBarry Smith   case 3:
2204b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_3;
2214b9ad928SBarry Smith     break;
2224b9ad928SBarry Smith   case 4:
2234b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_4;
2244b9ad928SBarry Smith     break;
2254b9ad928SBarry Smith   case 5:
2264b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_5;
2274b9ad928SBarry Smith     break;
2280e1b4bd6SMark F. Adams   case 6:
2290e1b4bd6SMark F. Adams     pc->ops->apply = PCApply_PBJacobi_6;
2300e1b4bd6SMark F. Adams     break;
2310a4ca348SMatthew G Knepley   case 7:
2320a4ca348SMatthew G Knepley     pc->ops->apply = PCApply_PBJacobi_7;
2330a4ca348SMatthew G Knepley     break;
2344b9ad928SBarry Smith   default:
235ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
2364b9ad928SBarry Smith   }
2374b9ad928SBarry Smith   PetscFunctionReturn(0);
2384b9ad928SBarry Smith }
2394b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2406849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
2414b9ad928SBarry Smith {
242dfbe8321SBarry Smith   PetscErrorCode ierr;
2434b9ad928SBarry Smith 
2444b9ad928SBarry Smith   PetscFunctionBegin;
2454b9ad928SBarry Smith   /*
2464b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
2474b9ad928SBarry Smith   */
248c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2494b9ad928SBarry Smith   PetscFunctionReturn(0);
2504b9ad928SBarry Smith }
251fa939db7SJed Brown 
252fa939db7SJed Brown static PetscErrorCode PCView_PBJacobi(PC pc,PetscViewer viewer)
253fa939db7SJed Brown {
254fa939db7SJed Brown   PetscErrorCode ierr;
255fa939db7SJed Brown   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
256fa939db7SJed Brown   PetscBool      iascii;
257fa939db7SJed Brown 
258fa939db7SJed Brown   PetscFunctionBegin;
259251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
260fa939db7SJed Brown   if (iascii) {
261*efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  point-block size %D\n",jac->bs);CHKERRQ(ierr);
26211aeaf0aSBarry Smith   }
263fa939db7SJed Brown   PetscFunctionReturn(0);
264fa939db7SJed Brown }
265fa939db7SJed Brown 
2664b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
26737a17b4dSBarry Smith /*MC
268422a814eSBarry Smith      PCPBJACOBI - Point block Jacobi preconditioner
269422a814eSBarry Smith 
270422a814eSBarry Smith 
271422a814eSBarry Smith    Notes: See PCJACOBI for point Jacobi preconditioning
272422a814eSBarry Smith 
273422a814eSBarry Smith    This works for AIJ and BAIJ matrices and uses the blocksize provided to the matrix
274422a814eSBarry Smith 
275422a814eSBarry Smith    Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
276422a814eSBarry Smith    is detected a PETSc error is generated.
277422a814eSBarry Smith 
278422a814eSBarry Smith    Developer Notes: This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow
279422a814eSBarry Smith    the factorization to continue even after a zero pivot is found resulting in a Nan and hence
280422a814eSBarry Smith    terminating KSP with a KSP_DIVERGED_NANORIF allowing
281422a814eSBarry Smith    a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
282422a814eSBarry Smith 
283422a814eSBarry Smith    Developer Note: Perhaps should provide an option that allows generation of a valid preconditioner
284422a814eSBarry Smith    even if a block is singular as the PCJACOBI does.
28537a17b4dSBarry Smith 
28637a17b4dSBarry Smith    Level: beginner
28737a17b4dSBarry Smith 
28837a17b4dSBarry Smith   Concepts: point block Jacobi
28937a17b4dSBarry Smith 
29037a17b4dSBarry Smith 
291422a814eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI
29237a17b4dSBarry Smith 
29337a17b4dSBarry Smith M*/
29437a17b4dSBarry Smith 
2958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc)
2964b9ad928SBarry Smith {
2974b9ad928SBarry Smith   PC_PBJacobi    *jac;
298dfbe8321SBarry Smith   PetscErrorCode ierr;
2994b9ad928SBarry Smith 
3004b9ad928SBarry Smith   PetscFunctionBegin;
3014b9ad928SBarry Smith   /*
3024b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
3034b9ad928SBarry Smith      attach it to the PC object.
3044b9ad928SBarry Smith   */
305b00a9115SJed Brown   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
3064b9ad928SBarry Smith   pc->data = (void*)jac;
3074b9ad928SBarry Smith 
3084b9ad928SBarry Smith   /*
3094b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
3104b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
3114b9ad928SBarry Smith   */
3124b9ad928SBarry Smith   jac->diag = 0;
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith   /*
3154b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
3164b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
3174b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
3184b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
3194b9ad928SBarry Smith       not needed.
3204b9ad928SBarry Smith   */
3214b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
3224b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
3234b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
3244b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
3254b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
326fa939db7SJed Brown   pc->ops->view                = PCView_PBJacobi;
3274b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
3284b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
3294b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
3304b9ad928SBarry Smith   PetscFunctionReturn(0);
3314b9ad928SBarry Smith }
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith 
334