xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision 00e125f844f1e8da6dd97b92f2bb419dde0a009f)
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/matimpl.h>
8af0996ceSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
94b9ad928SBarry Smith 
104b9ad928SBarry Smith /*
114b9ad928SBarry Smith    Private context (data structure) for the PBJacobi preconditioner.
124b9ad928SBarry Smith */
134b9ad928SBarry Smith typedef struct {
14713ccfa9SJed Brown   const MatScalar *diag;
15c1ac3661SBarry Smith   PetscInt        bs,mbs;
164b9ad928SBarry Smith } PC_PBJacobi;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith 
194b9ad928SBarry Smith #undef __FUNCT__
20bbead8a2SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_1"
21bbead8a2SBarry Smith static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y)
22bbead8a2SBarry Smith {
23bbead8a2SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi*)pc->data;
24bbead8a2SBarry Smith   PetscErrorCode    ierr;
25bbead8a2SBarry Smith   PetscInt          i,m = jac->mbs;
26bbead8a2SBarry Smith   const MatScalar   *diag = jac->diag;
27bbead8a2SBarry Smith   const PetscScalar *xx;
28bbead8a2SBarry Smith   PetscScalar       *yy;
29bbead8a2SBarry Smith 
30bbead8a2SBarry Smith   PetscFunctionBegin;
31bbead8a2SBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
32bbead8a2SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
332fa5cd67SKarl Rupp   for (i=0; i<m; i++) yy[i] = diag[i]*xx[i];
34bbead8a2SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
35bbead8a2SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
36bbead8a2SBarry Smith   ierr = PetscLogFlops(2.0*m);CHKERRQ(ierr);
37bbead8a2SBarry Smith   PetscFunctionReturn(0);
38bbead8a2SBarry Smith }
39bbead8a2SBarry Smith 
40bbead8a2SBarry Smith #undef __FUNCT__
414b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
426849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
434b9ad928SBarry Smith {
444b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
45dfbe8321SBarry Smith   PetscErrorCode  ierr;
46c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
4785f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
480e666e78SToby Isaac   PetscScalar     x0,x1,*yy;
490e666e78SToby Isaac   const PetscScalar *xx;
504b9ad928SBarry Smith 
514b9ad928SBarry Smith   PetscFunctionBegin;
520e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
534b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
544b9ad928SBarry Smith   for (i=0; i<m; i++) {
554b9ad928SBarry Smith     x0        = xx[2*i]; x1 = xx[2*i+1];
564b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
574b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
584b9ad928SBarry Smith     diag     += 4;
594b9ad928SBarry Smith   }
600e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
614b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
62dc0b31edSSatish Balay   ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr);
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith #undef __FUNCT__
664b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
676849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
684b9ad928SBarry Smith {
694b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
70dfbe8321SBarry Smith   PetscErrorCode  ierr;
71c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
7285f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
730e666e78SToby Isaac   PetscScalar     x0,x1,x2,*yy;
740e666e78SToby Isaac   const PetscScalar *xx;
754b9ad928SBarry Smith 
764b9ad928SBarry Smith   PetscFunctionBegin;
770e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
784b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
794b9ad928SBarry Smith   for (i=0; i<m; i++) {
804b9ad928SBarry Smith     x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
812fa5cd67SKarl Rupp 
824b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
834b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
844b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
854b9ad928SBarry Smith     diag     += 9;
864b9ad928SBarry Smith   }
870e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
884b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
89dc0b31edSSatish Balay   ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr);
904b9ad928SBarry Smith   PetscFunctionReturn(0);
914b9ad928SBarry Smith }
924b9ad928SBarry Smith #undef __FUNCT__
934b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4"
946849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
954b9ad928SBarry Smith {
964b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
97dfbe8321SBarry Smith   PetscErrorCode  ierr;
98c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
9985f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
1000e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,*yy;
1010e666e78SToby Isaac   const PetscScalar *xx;
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith   PetscFunctionBegin;
1040e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1054b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1064b9ad928SBarry Smith   for (i=0; i<m; i++) {
1074b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
1082fa5cd67SKarl Rupp 
1094b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
1104b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
1114b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
1124b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
1134b9ad928SBarry Smith     diag     += 16;
1144b9ad928SBarry Smith   }
1150e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1164b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
117dc0b31edSSatish Balay   ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr);
1184b9ad928SBarry Smith   PetscFunctionReturn(0);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith #undef __FUNCT__
1214b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
1226849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1234b9ad928SBarry Smith {
1244b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
125dfbe8321SBarry Smith   PetscErrorCode  ierr;
126c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
12785f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
1280e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,*yy;
1290e666e78SToby Isaac   const PetscScalar *xx;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
1320e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1334b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1344b9ad928SBarry Smith   for (i=0; i<m; i++) {
1354b9ad928SBarry 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];
1362fa5cd67SKarl Rupp 
1374b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1384b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1394b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1404b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1414b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1424b9ad928SBarry Smith     diag     += 25;
1434b9ad928SBarry Smith   }
1440e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1454b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
146dc0b31edSSatish Balay   ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr);
1474b9ad928SBarry Smith   PetscFunctionReturn(0);
1484b9ad928SBarry Smith }
1490e1b4bd6SMark F. Adams #undef __FUNCT__
1500e1b4bd6SMark F. Adams #define __FUNCT__ "PCApply_PBJacobi_6"
1510e1b4bd6SMark F. Adams static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y)
1520e1b4bd6SMark F. Adams {
1530e1b4bd6SMark F. Adams   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1540e1b4bd6SMark F. Adams   PetscErrorCode  ierr;
1550e1b4bd6SMark F. Adams   PetscInt        i,m = jac->mbs;
1560e1b4bd6SMark F. Adams   const MatScalar *diag = jac->diag;
1570e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,x5,*yy;
1580e666e78SToby Isaac   const PetscScalar *xx;
1590e1b4bd6SMark F. Adams 
1600e1b4bd6SMark F. Adams   PetscFunctionBegin;
1610e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1620e1b4bd6SMark F. Adams   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1630e1b4bd6SMark F. Adams   for (i=0; i<m; i++) {
1640e1b4bd6SMark 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];
1652fa5cd67SKarl Rupp 
1660e1b4bd6SMark F. Adams     yy[6*i]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
1670e1b4bd6SMark F. Adams     yy[6*i+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
1680e1b4bd6SMark F. Adams     yy[6*i+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
1690e1b4bd6SMark F. Adams     yy[6*i+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
1700e1b4bd6SMark F. Adams     yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
1710e1b4bd6SMark F. Adams     yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
1720e1b4bd6SMark F. Adams     diag     += 36;
1730e1b4bd6SMark F. Adams   }
1740e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1750e1b4bd6SMark F. Adams   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1760e1b4bd6SMark F. Adams   ierr = PetscLogFlops(66.0*m);CHKERRQ(ierr);
1770e1b4bd6SMark F. Adams   PetscFunctionReturn(0);
1780e1b4bd6SMark F. Adams }
1790a4ca348SMatthew G Knepley #undef __FUNCT__
1800a4ca348SMatthew G Knepley #define __FUNCT__ "PCApply_PBJacobi_7"
1810a4ca348SMatthew G Knepley static PetscErrorCode PCApply_PBJacobi_7(PC pc,Vec x,Vec y)
1820a4ca348SMatthew G Knepley {
1830a4ca348SMatthew G Knepley   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1840a4ca348SMatthew G Knepley   PetscErrorCode  ierr;
1850a4ca348SMatthew G Knepley   PetscInt        i,m = jac->mbs;
1860a4ca348SMatthew G Knepley   const MatScalar *diag = jac->diag;
1870e666e78SToby Isaac   PetscScalar     x0,x1,x2,x3,x4,x5,x6,*yy;
1880e666e78SToby Isaac   const PetscScalar *xx;
1890a4ca348SMatthew G Knepley 
1900a4ca348SMatthew G Knepley   PetscFunctionBegin;
1910e666e78SToby Isaac   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1920a4ca348SMatthew G Knepley   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1930a4ca348SMatthew G Knepley   for (i=0; i<m; i++) {
1940a4ca348SMatthew 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];
1952fa5cd67SKarl Rupp 
1960a4ca348SMatthew 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;
1970a4ca348SMatthew 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;
1980a4ca348SMatthew 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;
1990a4ca348SMatthew 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;
2000a4ca348SMatthew 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;
2010a4ca348SMatthew 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;
2020a4ca348SMatthew 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;
2030a4ca348SMatthew G Knepley     diag     += 49;
2040a4ca348SMatthew G Knepley   }
2050e666e78SToby Isaac   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
2060a4ca348SMatthew G Knepley   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
2070a4ca348SMatthew G Knepley   ierr = PetscLogFlops(80.0*m);CHKERRQ(ierr);
2080a4ca348SMatthew G Knepley   PetscFunctionReturn(0);
2090a4ca348SMatthew G Knepley }
2104b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2114b9ad928SBarry Smith #undef __FUNCT__
2124b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
2136849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
2144b9ad928SBarry Smith {
2154b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
216dfbe8321SBarry Smith   PetscErrorCode ierr;
2174b9ad928SBarry Smith   Mat            A = pc->pmat;
218*00e125f8SBarry Smith   MatFactorError err;
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith   PetscFunctionBegin;
221e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage");
222a1d92eedSBarry Smith 
223bbead8a2SBarry Smith   ierr = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr);
224*00e125f8SBarry Smith   ierr = MatFactorGetError(A,&err);CHKERRQ(ierr);
225*00e125f8SBarry Smith   if (err) pc->failedreason = (PCFailedReason)err;
226d7806518SHong Zhang 
22733d57670SJed Brown   ierr = MatGetBlockSize(A,&jac->bs);CHKERRQ(ierr);
22833d57670SJed Brown   jac->mbs = A->rmap->n/jac->bs;
229521d7252SBarry Smith   switch (jac->bs) {
230bbead8a2SBarry Smith   case 1:
231bbead8a2SBarry Smith     pc->ops->apply = PCApply_PBJacobi_1;
232bbead8a2SBarry Smith     break;
2334b9ad928SBarry Smith   case 2:
2344b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_2;
2354b9ad928SBarry Smith     break;
2364b9ad928SBarry Smith   case 3:
2374b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_3;
2384b9ad928SBarry Smith     break;
2394b9ad928SBarry Smith   case 4:
2404b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_4;
2414b9ad928SBarry Smith     break;
2424b9ad928SBarry Smith   case 5:
2434b9ad928SBarry Smith     pc->ops->apply = PCApply_PBJacobi_5;
2444b9ad928SBarry Smith     break;
2450e1b4bd6SMark F. Adams   case 6:
2460e1b4bd6SMark F. Adams     pc->ops->apply = PCApply_PBJacobi_6;
2470e1b4bd6SMark F. Adams     break;
2480a4ca348SMatthew G Knepley   case 7:
2490a4ca348SMatthew G Knepley     pc->ops->apply = PCApply_PBJacobi_7;
2500a4ca348SMatthew G Knepley     break;
2514b9ad928SBarry Smith   default:
252ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
2534b9ad928SBarry Smith   }
2544b9ad928SBarry Smith   PetscFunctionReturn(0);
2554b9ad928SBarry Smith }
2564b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2574b9ad928SBarry Smith #undef __FUNCT__
2584b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
2596849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
2604b9ad928SBarry Smith {
261dfbe8321SBarry Smith   PetscErrorCode ierr;
2624b9ad928SBarry Smith 
2634b9ad928SBarry Smith   PetscFunctionBegin;
2644b9ad928SBarry Smith   /*
2654b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
2664b9ad928SBarry Smith   */
267c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2684b9ad928SBarry Smith   PetscFunctionReturn(0);
2694b9ad928SBarry Smith }
270fa939db7SJed Brown 
271fa939db7SJed Brown #undef __FUNCT__
272fa939db7SJed Brown #define __FUNCT__ "PCView_PBJacobi"
273fa939db7SJed Brown static PetscErrorCode PCView_PBJacobi(PC pc,PetscViewer viewer)
274fa939db7SJed Brown {
275fa939db7SJed Brown   PetscErrorCode ierr;
276fa939db7SJed Brown   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
277fa939db7SJed Brown   PetscBool      iascii;
278fa939db7SJed Brown 
279fa939db7SJed Brown   PetscFunctionBegin;
280251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
281fa939db7SJed Brown   if (iascii) {
282fa939db7SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  point-block Jacobi: block size %D\n",jac->bs);CHKERRQ(ierr);
28311aeaf0aSBarry Smith   }
284fa939db7SJed Brown   PetscFunctionReturn(0);
285fa939db7SJed Brown }
286fa939db7SJed Brown 
2874b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
28837a17b4dSBarry Smith /*MC
289422a814eSBarry Smith      PCPBJACOBI - Point block Jacobi preconditioner
290422a814eSBarry Smith 
291422a814eSBarry Smith 
292422a814eSBarry Smith    Notes: See PCJACOBI for point Jacobi preconditioning
293422a814eSBarry Smith 
294422a814eSBarry Smith    This works for AIJ and BAIJ matrices and uses the blocksize provided to the matrix
295422a814eSBarry Smith 
296422a814eSBarry Smith    Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot
297422a814eSBarry Smith    is detected a PETSc error is generated.
298422a814eSBarry Smith 
299422a814eSBarry Smith    Developer Notes: This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow
300422a814eSBarry Smith    the factorization to continue even after a zero pivot is found resulting in a Nan and hence
301422a814eSBarry Smith    terminating KSP with a KSP_DIVERGED_NANORIF allowing
302422a814eSBarry Smith    a nonlinear solver/ODE integrator to recover without stopping the program as currently happens.
303422a814eSBarry Smith 
304422a814eSBarry Smith    Developer Note: Perhaps should provide an option that allows generation of a valid preconditioner
305422a814eSBarry Smith    even if a block is singular as the PCJACOBI does.
30637a17b4dSBarry Smith 
30737a17b4dSBarry Smith    Level: beginner
30837a17b4dSBarry Smith 
30937a17b4dSBarry Smith   Concepts: point block Jacobi
31037a17b4dSBarry Smith 
31137a17b4dSBarry Smith 
312422a814eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI
31337a17b4dSBarry Smith 
31437a17b4dSBarry Smith M*/
31537a17b4dSBarry Smith 
3164b9ad928SBarry Smith #undef __FUNCT__
3174b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
3188cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PBJacobi(PC pc)
3194b9ad928SBarry Smith {
3204b9ad928SBarry Smith   PC_PBJacobi    *jac;
321dfbe8321SBarry Smith   PetscErrorCode ierr;
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith   PetscFunctionBegin;
3244b9ad928SBarry Smith   /*
3254b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
3264b9ad928SBarry Smith      attach it to the PC object.
3274b9ad928SBarry Smith   */
328b00a9115SJed Brown   ierr     = PetscNewLog(pc,&jac);CHKERRQ(ierr);
3294b9ad928SBarry Smith   pc->data = (void*)jac;
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith   /*
3324b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
3334b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
3344b9ad928SBarry Smith   */
3354b9ad928SBarry Smith   jac->diag = 0;
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith   /*
3384b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
3394b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
3404b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
3414b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
3424b9ad928SBarry Smith       not needed.
3434b9ad928SBarry Smith   */
3444b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
3454b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
3464b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
3474b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
3484b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
349fa939db7SJed Brown   pc->ops->view                = PCView_PBJacobi;
3504b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
3514b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
3524b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
3534b9ad928SBarry Smith   PetscFunctionReturn(0);
3544b9ad928SBarry Smith }
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith 
357