xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
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 
7*b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
8*b45d2f2cSJed Brown #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 {
1485f4f44aSBarry Smith   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);
33bbead8a2SBarry Smith   for (i=0; i<m; i++) {
34bbead8a2SBarry Smith     yy[i] = diag[i]*xx[i];
35bbead8a2SBarry Smith   }
36bbead8a2SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
37bbead8a2SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
38bbead8a2SBarry Smith   ierr = PetscLogFlops(2.0*m);CHKERRQ(ierr);
39bbead8a2SBarry Smith   PetscFunctionReturn(0);
40bbead8a2SBarry Smith }
41bbead8a2SBarry Smith 
42bbead8a2SBarry Smith #undef __FUNCT__
434b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
446849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
454b9ad928SBarry Smith {
464b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
47dfbe8321SBarry Smith   PetscErrorCode  ierr;
48c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
4985f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
5085f4f44aSBarry Smith   PetscScalar     x0,x1,*xx,*yy;
514b9ad928SBarry Smith 
524b9ad928SBarry Smith   PetscFunctionBegin;
534b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
544b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
554b9ad928SBarry Smith   for (i=0; i<m; i++) {
564b9ad928SBarry Smith     x0 = xx[2*i]; x1 = xx[2*i+1];
574b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
584b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
594b9ad928SBarry Smith     diag     += 4;
604b9ad928SBarry Smith   }
614b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
624b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
63dc0b31edSSatish Balay   ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr);
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith #undef __FUNCT__
674b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
686849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
694b9ad928SBarry Smith {
704b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
71dfbe8321SBarry Smith   PetscErrorCode  ierr;
72c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
7385f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
7485f4f44aSBarry Smith   PetscScalar     x0,x1,x2,*xx,*yy;
754b9ad928SBarry Smith 
764b9ad928SBarry Smith   PetscFunctionBegin;
774b9ad928SBarry Smith   ierr = VecGetArray(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];
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   }
864b9ad928SBarry Smith   ierr = VecRestoreArray(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;
9985f4f44aSBarry Smith   PetscScalar      x0,x1,x2,x3,*xx,*yy;
1004b9ad928SBarry Smith 
1014b9ad928SBarry Smith   PetscFunctionBegin;
1024b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1034b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1044b9ad928SBarry Smith   for (i=0; i<m; i++) {
1054b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
1064b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
1074b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
1084b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
1094b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
1104b9ad928SBarry Smith     diag     += 16;
1114b9ad928SBarry Smith   }
1124b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1134b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
114dc0b31edSSatish Balay   ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr);
1154b9ad928SBarry Smith   PetscFunctionReturn(0);
1164b9ad928SBarry Smith }
1174b9ad928SBarry Smith #undef __FUNCT__
1184b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
1196849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1204b9ad928SBarry Smith {
1214b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
122dfbe8321SBarry Smith   PetscErrorCode  ierr;
123c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
12485f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
12585f4f44aSBarry Smith   PetscScalar     x0,x1,x2,x3,x4,*xx,*yy;
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith   PetscFunctionBegin;
1284b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1294b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1304b9ad928SBarry Smith   for (i=0; i<m; i++) {
1314b9ad928SBarry 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];
1324b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1334b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1344b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1354b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1364b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1374b9ad928SBarry Smith     diag     += 25;
1384b9ad928SBarry Smith   }
1394b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1404b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
141dc0b31edSSatish Balay   ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr);
1424b9ad928SBarry Smith   PetscFunctionReturn(0);
1434b9ad928SBarry Smith }
1440e1b4bd6SMark F. Adams #undef __FUNCT__
1450e1b4bd6SMark F. Adams #define __FUNCT__ "PCApply_PBJacobi_6"
1460e1b4bd6SMark F. Adams static PetscErrorCode PCApply_PBJacobi_6(PC pc,Vec x,Vec y)
1470e1b4bd6SMark F. Adams {
1480e1b4bd6SMark F. Adams   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1490e1b4bd6SMark F. Adams   PetscErrorCode  ierr;
1500e1b4bd6SMark F. Adams   PetscInt        i,m = jac->mbs;
1510e1b4bd6SMark F. Adams   const MatScalar *diag = jac->diag;
1520e1b4bd6SMark F. Adams   PetscScalar     x0,x1,x2,x3,x4,x5,*xx,*yy;
1530e1b4bd6SMark F. Adams 
1540e1b4bd6SMark F. Adams   PetscFunctionBegin;
1550e1b4bd6SMark F. Adams   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1560e1b4bd6SMark F. Adams   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1570e1b4bd6SMark F. Adams   for (i=0; i<m; i++) {
1580e1b4bd6SMark 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];
1590e1b4bd6SMark F. Adams     yy[6*i]   = diag[0]*x0 + diag[6]*x1  + diag[12]*x2  + diag[18]*x3 + diag[24]*x4 + diag[30]*x5;
1600e1b4bd6SMark F. Adams     yy[6*i+1] = diag[1]*x0 + diag[7]*x1  + diag[13]*x2  + diag[19]*x3 + diag[25]*x4 + diag[31]*x5;
1610e1b4bd6SMark F. Adams     yy[6*i+2] = diag[2]*x0 + diag[8]*x1  + diag[14]*x2  + diag[20]*x3 + diag[26]*x4 + diag[32]*x5;
1620e1b4bd6SMark F. Adams     yy[6*i+3] = diag[3]*x0 + diag[9]*x1  + diag[15]*x2  + diag[21]*x3 + diag[27]*x4 + diag[33]*x5;
1630e1b4bd6SMark F. Adams     yy[6*i+4] = diag[4]*x0 + diag[10]*x1 + diag[16]*x2  + diag[22]*x3 + diag[28]*x4 + diag[34]*x5;
1640e1b4bd6SMark F. Adams     yy[6*i+5] = diag[5]*x0 + diag[11]*x1 + diag[17]*x2  + diag[23]*x3 + diag[29]*x4 + diag[35]*x5;
1650e1b4bd6SMark F. Adams     diag     += 36;
1660e1b4bd6SMark F. Adams   }
1670e1b4bd6SMark F. Adams   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1680e1b4bd6SMark F. Adams   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1690e1b4bd6SMark F. Adams   ierr = PetscLogFlops(66.0*m);CHKERRQ(ierr);
1700e1b4bd6SMark F. Adams   PetscFunctionReturn(0);
1710e1b4bd6SMark F. Adams }
1720a4ca348SMatthew G Knepley #undef __FUNCT__
1730a4ca348SMatthew G Knepley #define __FUNCT__ "PCApply_PBJacobi_7"
1740a4ca348SMatthew G Knepley static PetscErrorCode PCApply_PBJacobi_7(PC pc,Vec x,Vec y)
1750a4ca348SMatthew G Knepley {
1760a4ca348SMatthew G Knepley   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
1770a4ca348SMatthew G Knepley   PetscErrorCode  ierr;
1780a4ca348SMatthew G Knepley   PetscInt        i,m = jac->mbs;
1790a4ca348SMatthew G Knepley   const MatScalar *diag = jac->diag;
1800a4ca348SMatthew G Knepley   PetscScalar     x0,x1,x2,x3,x4,x5,x6,*xx,*yy;
1810a4ca348SMatthew G Knepley 
1820a4ca348SMatthew G Knepley   PetscFunctionBegin;
1830a4ca348SMatthew G Knepley   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1840a4ca348SMatthew G Knepley   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1850a4ca348SMatthew G Knepley   for (i=0; i<m; i++) {
1860a4ca348SMatthew 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];
1870a4ca348SMatthew 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;
1880a4ca348SMatthew 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;
1890a4ca348SMatthew 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;
1900a4ca348SMatthew 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;
1910a4ca348SMatthew 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;
1920a4ca348SMatthew 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;
1930a4ca348SMatthew 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;
1940a4ca348SMatthew G Knepley     diag     += 49;
1950a4ca348SMatthew G Knepley   }
1960a4ca348SMatthew G Knepley   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1970a4ca348SMatthew G Knepley   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1980a4ca348SMatthew G Knepley   ierr = PetscLogFlops(80.0*m);CHKERRQ(ierr);
1990a4ca348SMatthew G Knepley   PetscFunctionReturn(0);
2000a4ca348SMatthew G Knepley }
2014b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2024b9ad928SBarry Smith #undef __FUNCT__
2034b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
2046849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
2054b9ad928SBarry Smith {
2064b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
207dfbe8321SBarry Smith   PetscErrorCode ierr;
2084b9ad928SBarry Smith   Mat            A = pc->pmat;
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith   PetscFunctionBegin;
211e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage");
212a1d92eedSBarry Smith 
213bbead8a2SBarry Smith   ierr        = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr);
214d0f46423SBarry Smith   jac->bs     = A->rmap->bs;
215bbead8a2SBarry Smith   jac->mbs    = A->rmap->n/A->rmap->bs;
216521d7252SBarry Smith   switch (jac->bs){
217bbead8a2SBarry Smith     case 1:
218bbead8a2SBarry Smith       pc->ops->apply = PCApply_PBJacobi_1;
219bbead8a2SBarry Smith       break;
2204b9ad928SBarry Smith     case 2:
2214b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_2;
2224b9ad928SBarry Smith       break;
2234b9ad928SBarry Smith     case 3:
2244b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_3;
2254b9ad928SBarry Smith       break;
2264b9ad928SBarry Smith     case 4:
2274b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_4;
2284b9ad928SBarry Smith       break;
2294b9ad928SBarry Smith     case 5:
2304b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_5;
2314b9ad928SBarry Smith       break;
2320e1b4bd6SMark F. Adams     case 6:
2330e1b4bd6SMark F. Adams       pc->ops->apply = PCApply_PBJacobi_6;
2340e1b4bd6SMark F. Adams       break;
2350a4ca348SMatthew G Knepley     case 7:
2360a4ca348SMatthew G Knepley       pc->ops->apply = PCApply_PBJacobi_7;
2370a4ca348SMatthew G Knepley       break;
2384b9ad928SBarry Smith     default:
23965e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
2404b9ad928SBarry Smith   }
2414b9ad928SBarry Smith 
2424b9ad928SBarry Smith   PetscFunctionReturn(0);
2434b9ad928SBarry Smith }
2444b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
2454b9ad928SBarry Smith #undef __FUNCT__
2464b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
2476849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
2484b9ad928SBarry Smith {
249dfbe8321SBarry Smith   PetscErrorCode ierr;
2504b9ad928SBarry Smith 
2514b9ad928SBarry Smith   PetscFunctionBegin;
2524b9ad928SBarry Smith   /*
2534b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
2544b9ad928SBarry Smith   */
255c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2564b9ad928SBarry Smith   PetscFunctionReturn(0);
2574b9ad928SBarry Smith }
2584b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
25937a17b4dSBarry Smith /*MC
26037a17b4dSBarry Smith      PCPBJACOBI - Point block Jacobi
26137a17b4dSBarry Smith 
26237a17b4dSBarry Smith    Level: beginner
26337a17b4dSBarry Smith 
26437a17b4dSBarry Smith   Concepts: point block Jacobi
26537a17b4dSBarry Smith 
26637a17b4dSBarry Smith 
26737a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
26837a17b4dSBarry Smith 
26937a17b4dSBarry Smith M*/
27037a17b4dSBarry Smith 
2714b9ad928SBarry Smith EXTERN_C_BEGIN
2724b9ad928SBarry Smith #undef __FUNCT__
2734b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
2747087cfbeSBarry Smith PetscErrorCode  PCCreate_PBJacobi(PC pc)
2754b9ad928SBarry Smith {
2764b9ad928SBarry Smith   PC_PBJacobi    *jac;
277dfbe8321SBarry Smith   PetscErrorCode ierr;
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith   PetscFunctionBegin;
2804b9ad928SBarry Smith 
2814b9ad928SBarry Smith   /*
2824b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
2834b9ad928SBarry Smith      attach it to the PC object.
2844b9ad928SBarry Smith   */
28538f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(pc,PC_PBJacobi,&jac);CHKERRQ(ierr);
2864b9ad928SBarry Smith   pc->data  = (void*)jac;
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith   /*
2894b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
2904b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
2914b9ad928SBarry Smith   */
2924b9ad928SBarry Smith   jac->diag          = 0;
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith   /*
2954b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
2964b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2974b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
2984b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
2994b9ad928SBarry Smith       not needed.
3004b9ad928SBarry Smith   */
3014b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
3024b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
3034b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
3044b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
3054b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
3064b9ad928SBarry Smith   pc->ops->view                = 0;
3074b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
3084b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
3094b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
3104b9ad928SBarry Smith   PetscFunctionReturn(0);
3114b9ad928SBarry Smith }
3124b9ad928SBarry Smith EXTERN_C_END
3134b9ad928SBarry Smith 
3144b9ad928SBarry Smith 
315