xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision bbead8a294046a48bbbe8fbcc40ad4dafe54ecae)
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 
7c6db04a5SJed Brown #include <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 {
1385f4f44aSBarry Smith   MatScalar   *diag;
14c1ac3661SBarry Smith   PetscInt    bs,mbs;
154b9ad928SBarry Smith } PC_PBJacobi;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith /*
184b9ad928SBarry Smith    Currently only implemented for baij matrices and directly access baij
194b9ad928SBarry Smith   data structures.
204b9ad928SBarry Smith */
21c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h>
22c6db04a5SJed Brown #include <../src/mat/blockinvert.h>
234b9ad928SBarry Smith 
244b9ad928SBarry Smith #undef __FUNCT__
25*bbead8a2SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_1"
26*bbead8a2SBarry Smith static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y)
27*bbead8a2SBarry Smith {
28*bbead8a2SBarry Smith   PC_PBJacobi       *jac = (PC_PBJacobi*)pc->data;
29*bbead8a2SBarry Smith   PetscErrorCode    ierr;
30*bbead8a2SBarry Smith   PetscInt          i,m = jac->mbs;
31*bbead8a2SBarry Smith   const MatScalar   *diag = jac->diag;
32*bbead8a2SBarry Smith   const PetscScalar *xx;
33*bbead8a2SBarry Smith   PetscScalar       *yy;
34*bbead8a2SBarry Smith 
35*bbead8a2SBarry Smith   PetscFunctionBegin;
36*bbead8a2SBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
37*bbead8a2SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
38*bbead8a2SBarry Smith   for (i=0; i<m; i++) {
39*bbead8a2SBarry Smith     yy[i] = diag[i]*xx[i];
40*bbead8a2SBarry Smith   }
41*bbead8a2SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
42*bbead8a2SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
43*bbead8a2SBarry Smith   ierr = PetscLogFlops(2.0*m);CHKERRQ(ierr);
44*bbead8a2SBarry Smith   PetscFunctionReturn(0);
45*bbead8a2SBarry Smith }
46*bbead8a2SBarry Smith 
47*bbead8a2SBarry Smith #undef __FUNCT__
484b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
496849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
504b9ad928SBarry Smith {
514b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
52dfbe8321SBarry Smith   PetscErrorCode  ierr;
53c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
5485f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
5585f4f44aSBarry Smith   PetscScalar     x0,x1,*xx,*yy;
564b9ad928SBarry Smith 
574b9ad928SBarry Smith   PetscFunctionBegin;
584b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
594b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
604b9ad928SBarry Smith   for (i=0; i<m; i++) {
614b9ad928SBarry Smith     x0 = xx[2*i]; x1 = xx[2*i+1];
624b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
634b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
644b9ad928SBarry Smith     diag     += 4;
654b9ad928SBarry Smith   }
664b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
674b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
68dc0b31edSSatish Balay   ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr);
694b9ad928SBarry Smith   PetscFunctionReturn(0);
704b9ad928SBarry Smith }
714b9ad928SBarry Smith #undef __FUNCT__
724b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
736849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
744b9ad928SBarry Smith {
754b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
76dfbe8321SBarry Smith   PetscErrorCode  ierr;
77c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
7885f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
7985f4f44aSBarry Smith   PetscScalar     x0,x1,x2,*xx,*yy;
804b9ad928SBarry Smith 
814b9ad928SBarry Smith   PetscFunctionBegin;
824b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
834b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
844b9ad928SBarry Smith   for (i=0; i<m; i++) {
854b9ad928SBarry Smith     x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
864b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
874b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
884b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
894b9ad928SBarry Smith     diag     += 9;
904b9ad928SBarry Smith   }
914b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
924b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
93dc0b31edSSatish Balay   ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr);
944b9ad928SBarry Smith   PetscFunctionReturn(0);
954b9ad928SBarry Smith }
964b9ad928SBarry Smith #undef __FUNCT__
974b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4"
986849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
994b9ad928SBarry Smith {
1004b9ad928SBarry Smith   PC_PBJacobi      *jac = (PC_PBJacobi*)pc->data;
101dfbe8321SBarry Smith   PetscErrorCode   ierr;
102c1ac3661SBarry Smith   PetscInt         i,m = jac->mbs;
10385f4f44aSBarry Smith   const MatScalar  *diag = jac->diag;
10485f4f44aSBarry Smith   PetscScalar      x0,x1,x2,x3,*xx,*yy;
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith   PetscFunctionBegin;
1074b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1084b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1094b9ad928SBarry Smith   for (i=0; i<m; i++) {
1104b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
1114b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
1124b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
1134b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
1144b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
1154b9ad928SBarry Smith     diag     += 16;
1164b9ad928SBarry Smith   }
1174b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1184b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
119dc0b31edSSatish Balay   ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr);
1204b9ad928SBarry Smith   PetscFunctionReturn(0);
1214b9ad928SBarry Smith }
1224b9ad928SBarry Smith #undef __FUNCT__
1234b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
1246849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1254b9ad928SBarry Smith {
1264b9ad928SBarry Smith   PC_PBJacobi     *jac = (PC_PBJacobi*)pc->data;
127dfbe8321SBarry Smith   PetscErrorCode  ierr;
128c1ac3661SBarry Smith   PetscInt        i,m = jac->mbs;
12985f4f44aSBarry Smith   const MatScalar *diag = jac->diag;
13085f4f44aSBarry Smith   PetscScalar     x0,x1,x2,x3,x4,*xx,*yy;
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith   PetscFunctionBegin;
1334b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1344b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1354b9ad928SBarry Smith   for (i=0; i<m; i++) {
1364b9ad928SBarry 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];
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   }
1444b9ad928SBarry Smith   ierr = VecRestoreArray(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 }
1494b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1504b9ad928SBarry Smith #undef __FUNCT__
1514b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
1526849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
1534b9ad928SBarry Smith {
1544b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
155dfbe8321SBarry Smith   PetscErrorCode ierr;
1564b9ad928SBarry Smith   Mat            A = pc->pmat;
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith   PetscFunctionBegin;
159e32f2f54SBarry Smith   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage");
160a1d92eedSBarry Smith 
161*bbead8a2SBarry Smith   ierr        = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr);
162d0f46423SBarry Smith   jac->bs     = A->rmap->bs;
163*bbead8a2SBarry Smith   jac->mbs    = A->rmap->n/A->rmap->bs;
164521d7252SBarry Smith   switch (jac->bs){
165*bbead8a2SBarry Smith     case 1:
166*bbead8a2SBarry Smith       pc->ops->apply = PCApply_PBJacobi_1;
167*bbead8a2SBarry Smith       break;
1684b9ad928SBarry Smith     case 2:
1694b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_2;
1704b9ad928SBarry Smith       break;
1714b9ad928SBarry Smith     case 3:
1724b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_3;
1734b9ad928SBarry Smith       break;
1744b9ad928SBarry Smith     case 4:
1754b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_4;
1764b9ad928SBarry Smith       break;
1774b9ad928SBarry Smith     case 5:
1784b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_5;
1794b9ad928SBarry Smith       break;
1804b9ad928SBarry Smith     default:
18165e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
1824b9ad928SBarry Smith   }
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   PetscFunctionReturn(0);
1854b9ad928SBarry Smith }
1864b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1874b9ad928SBarry Smith #undef __FUNCT__
1884b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
1896849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
1904b9ad928SBarry Smith {
191dfbe8321SBarry Smith   PetscErrorCode ierr;
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith   PetscFunctionBegin;
1944b9ad928SBarry Smith   /*
1954b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
1964b9ad928SBarry Smith   */
197c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1984b9ad928SBarry Smith   PetscFunctionReturn(0);
1994b9ad928SBarry Smith }
2004b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
20137a17b4dSBarry Smith /*MC
20237a17b4dSBarry Smith      PCPBJACOBI - Point block Jacobi
20337a17b4dSBarry Smith 
20437a17b4dSBarry Smith    Level: beginner
20537a17b4dSBarry Smith 
20637a17b4dSBarry Smith   Concepts: point block Jacobi
20737a17b4dSBarry Smith 
20837a17b4dSBarry Smith    Notes: Only implemented for the BAIJ matrix formats.
20937a17b4dSBarry Smith 
21037a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
21137a17b4dSBarry Smith 
21237a17b4dSBarry Smith M*/
21337a17b4dSBarry Smith 
2144b9ad928SBarry Smith EXTERN_C_BEGIN
2154b9ad928SBarry Smith #undef __FUNCT__
2164b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
2177087cfbeSBarry Smith PetscErrorCode  PCCreate_PBJacobi(PC pc)
2184b9ad928SBarry Smith {
2194b9ad928SBarry Smith   PC_PBJacobi    *jac;
220dfbe8321SBarry Smith   PetscErrorCode ierr;
2214b9ad928SBarry Smith 
2224b9ad928SBarry Smith   PetscFunctionBegin;
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith   /*
2254b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
2264b9ad928SBarry Smith      attach it to the PC object.
2274b9ad928SBarry Smith   */
22838f2d2fdSLisandro Dalcin   ierr      = PetscNewLog(pc,PC_PBJacobi,&jac);CHKERRQ(ierr);
2294b9ad928SBarry Smith   pc->data  = (void*)jac;
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith   /*
2324b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
2334b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
2344b9ad928SBarry Smith   */
2354b9ad928SBarry Smith   jac->diag          = 0;
2364b9ad928SBarry Smith 
2374b9ad928SBarry Smith   /*
2384b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
2394b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2404b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
2414b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
2424b9ad928SBarry Smith       not needed.
2434b9ad928SBarry Smith   */
2444b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
2454b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
2464b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
2474b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
2484b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
2494b9ad928SBarry Smith   pc->ops->view                = 0;
2504b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
2514b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
2524b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
2534b9ad928SBarry Smith   PetscFunctionReturn(0);
2544b9ad928SBarry Smith }
2554b9ad928SBarry Smith EXTERN_C_END
2564b9ad928SBarry Smith 
2574b9ad928SBarry Smith 
258