xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith    Include files needed for the PBJacobi preconditioner:
34b9ad928SBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
44b9ad928SBarry Smith */
54b9ad928SBarry Smith 
64b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h"   /*I "petscpc.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith /*
94b9ad928SBarry Smith    Private context (data structure) for the PBJacobi preconditioner.
104b9ad928SBarry Smith */
114b9ad928SBarry Smith typedef struct {
124b9ad928SBarry Smith   PetscScalar *diag;
134b9ad928SBarry Smith   int         bs,mbs;
144b9ad928SBarry Smith } PC_PBJacobi;
154b9ad928SBarry Smith 
164b9ad928SBarry Smith /*
174b9ad928SBarry Smith    Currently only implemented for baij matrices and directly access baij
184b9ad928SBarry Smith   data structures.
194b9ad928SBarry Smith */
204b9ad928SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
214b9ad928SBarry Smith #include "src/inline/ilu.h"
224b9ad928SBarry Smith 
234b9ad928SBarry Smith #undef __FUNCT__
244b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
25*6849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
264b9ad928SBarry Smith {
274b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
28dfbe8321SBarry Smith   PetscErrorCode ierr;
29dfbe8321SBarry Smith   int i,m = jac->mbs;
304b9ad928SBarry Smith   PetscScalar *diag = jac->diag,x0,x1,*xx,*yy;
314b9ad928SBarry Smith 
324b9ad928SBarry Smith   PetscFunctionBegin;
334b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
344b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
354b9ad928SBarry Smith   for (i=0; i<m; i++) {
364b9ad928SBarry Smith     x0 = xx[2*i]; x1 = xx[2*i+1];
374b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
384b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
394b9ad928SBarry Smith     diag     += 4;
404b9ad928SBarry Smith   }
414b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
424b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
434b9ad928SBarry Smith   PetscLogFlops(6*m);
444b9ad928SBarry Smith   PetscFunctionReturn(0);
454b9ad928SBarry Smith }
464b9ad928SBarry Smith #undef __FUNCT__
474b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
48*6849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
494b9ad928SBarry Smith {
504b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
51dfbe8321SBarry Smith   PetscErrorCode ierr;
52dfbe8321SBarry Smith   int i,m = jac->mbs;
534b9ad928SBarry Smith   PetscScalar *diag = jac->diag,x0,x1,x2,*xx,*yy;
544b9ad928SBarry Smith 
554b9ad928SBarry Smith   PetscFunctionBegin;
564b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
574b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
584b9ad928SBarry Smith   for (i=0; i<m; i++) {
594b9ad928SBarry Smith     x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
604b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
614b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
624b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
634b9ad928SBarry Smith     diag     += 9;
644b9ad928SBarry Smith   }
654b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
664b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
674b9ad928SBarry Smith   PetscLogFlops(15*m);
684b9ad928SBarry Smith   PetscFunctionReturn(0);
694b9ad928SBarry Smith }
704b9ad928SBarry Smith #undef __FUNCT__
714b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4"
72*6849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
734b9ad928SBarry Smith {
744b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
75dfbe8321SBarry Smith   PetscErrorCode ierr;
76dfbe8321SBarry Smith   int i,m = jac->mbs;
774b9ad928SBarry Smith   PetscScalar *diag = jac->diag,x0,x1,x2,x3,*xx,*yy;
784b9ad928SBarry Smith 
794b9ad928SBarry Smith   PetscFunctionBegin;
804b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
814b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
824b9ad928SBarry Smith   for (i=0; i<m; i++) {
834b9ad928SBarry Smith     x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
844b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
854b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
864b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
874b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
884b9ad928SBarry Smith     diag     += 16;
894b9ad928SBarry Smith   }
904b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
914b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
924b9ad928SBarry Smith   PetscLogFlops(28*m);
934b9ad928SBarry Smith   PetscFunctionReturn(0);
944b9ad928SBarry Smith }
954b9ad928SBarry Smith #undef __FUNCT__
964b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
97*6849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
984b9ad928SBarry Smith {
994b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
100dfbe8321SBarry Smith   PetscErrorCode ierr;
101dfbe8321SBarry Smith   int i,m = jac->mbs;
1024b9ad928SBarry Smith   PetscScalar *diag = jac->diag,x0,x1,x2,x3,x4,*xx,*yy;
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith   PetscFunctionBegin;
1054b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1064b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1074b9ad928SBarry Smith   for (i=0; i<m; i++) {
1084b9ad928SBarry 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];
1094b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1104b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1114b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1124b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1134b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1144b9ad928SBarry Smith     diag     += 25;
1154b9ad928SBarry Smith   }
1164b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1174b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1184b9ad928SBarry Smith   PetscLogFlops(45*m);
1194b9ad928SBarry Smith   PetscFunctionReturn(0);
1204b9ad928SBarry Smith }
1214b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
122dfbe8321SBarry Smith EXTERN PetscErrorCode MatInvertBlockDiagonal_SeqBAIJ(Mat);
1234b9ad928SBarry Smith #undef __FUNCT__
1244b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
125*6849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
1264b9ad928SBarry Smith {
1274b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
128dfbe8321SBarry Smith   PetscErrorCode ierr;
129dfbe8321SBarry Smith   int size;
13021360622SBarry Smith   PetscTruth  seqbaij,mpibaij,baij;
1314b9ad928SBarry Smith   Mat         A = pc->pmat;
1324b9ad928SBarry Smith   Mat_SeqBAIJ *a;
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   PetscFunctionBegin;
1354b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);CHKERRQ(ierr);
1364b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);CHKERRQ(ierr);
13721360622SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);CHKERRQ(ierr);
13821360622SBarry Smith   if (!seqbaij && !mpibaij && !baij) {
1394b9ad928SBarry Smith     SETERRQ(1,"Currently only supports BAIJ matrices");
1404b9ad928SBarry Smith   }
14121360622SBarry Smith   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
14221360622SBarry Smith   if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A;
1434b9ad928SBarry Smith   if (A->m != A->n) SETERRQ(1,"Supported only for square matrices and square storage");
144a1d92eedSBarry Smith 
145a1d92eedSBarry Smith   ierr        =  MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
1464b9ad928SBarry Smith   a           = (Mat_SeqBAIJ*)A->data;
147a1d92eedSBarry Smith   jac->diag   = a->idiag;
1484b9ad928SBarry Smith   jac->bs     = a->bs;
1494b9ad928SBarry Smith   jac->mbs    = a->mbs;
1504b9ad928SBarry Smith   switch (a->bs){
1514b9ad928SBarry Smith     case 2:
1524b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_2;
1534b9ad928SBarry Smith       break;
1544b9ad928SBarry Smith     case 3:
1554b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_3;
1564b9ad928SBarry Smith       break;
1574b9ad928SBarry Smith     case 4:
1584b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_4;
1594b9ad928SBarry Smith       break;
1604b9ad928SBarry Smith     case 5:
1614b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_5;
1624b9ad928SBarry Smith       break;
1634b9ad928SBarry Smith     default:
1644b9ad928SBarry Smith       SETERRQ1(1,"not supported for block size %d",a->bs);
1654b9ad928SBarry Smith   }
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith   PetscFunctionReturn(0);
1684b9ad928SBarry Smith }
1694b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1704b9ad928SBarry Smith #undef __FUNCT__
1714b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
172*6849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
1734b9ad928SBarry Smith {
1744b9ad928SBarry Smith   PC_PBJacobi *jac = (PC_PBJacobi*)pc->data;
175dfbe8321SBarry Smith   PetscErrorCode ierr;
1764b9ad928SBarry Smith 
1774b9ad928SBarry Smith   PetscFunctionBegin;
1784b9ad928SBarry Smith   /*
1794b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
1804b9ad928SBarry Smith   */
1814b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1824b9ad928SBarry Smith   PetscFunctionReturn(0);
1834b9ad928SBarry Smith }
1844b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
18537a17b4dSBarry Smith /*MC
18637a17b4dSBarry Smith      PCPBJACOBI - Point block Jacobi
18737a17b4dSBarry Smith 
18837a17b4dSBarry Smith    Level: beginner
18937a17b4dSBarry Smith 
19037a17b4dSBarry Smith   Concepts: point block Jacobi
19137a17b4dSBarry Smith 
19237a17b4dSBarry Smith    Notes: Only implemented for the BAIJ matrix formats.
19337a17b4dSBarry Smith 
19437a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
19537a17b4dSBarry Smith 
19637a17b4dSBarry Smith M*/
19737a17b4dSBarry Smith 
1984b9ad928SBarry Smith EXTERN_C_BEGIN
1994b9ad928SBarry Smith #undef __FUNCT__
2004b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
201dfbe8321SBarry Smith PetscErrorCode PCCreate_PBJacobi(PC pc)
2024b9ad928SBarry Smith {
2034b9ad928SBarry Smith   PC_PBJacobi *jac;
204dfbe8321SBarry Smith   PetscErrorCode ierr;
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith   PetscFunctionBegin;
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith   /*
2094b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
2104b9ad928SBarry Smith      attach it to the PC object.
2114b9ad928SBarry Smith   */
2124b9ad928SBarry Smith   ierr      = PetscNew(PC_PBJacobi,&jac);CHKERRQ(ierr);
2134b9ad928SBarry Smith   pc->data  = (void*)jac;
2144b9ad928SBarry Smith 
2154b9ad928SBarry Smith   /*
2164b9ad928SBarry Smith      Logs the memory usage; this is not needed but allows PETSc to
2174b9ad928SBarry Smith      monitor how much memory is being used for various purposes.
2184b9ad928SBarry Smith   */
2194b9ad928SBarry Smith   PetscLogObjectMemory(pc,sizeof(PC_PBJacobi));
2204b9ad928SBarry Smith 
2214b9ad928SBarry Smith   /*
2224b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
2234b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
2244b9ad928SBarry Smith   */
2254b9ad928SBarry Smith   jac->diag          = 0;
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith   /*
2284b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
2294b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2304b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
2314b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
2324b9ad928SBarry Smith       not needed.
2334b9ad928SBarry Smith   */
2344b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
2354b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
2364b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
2374b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
2384b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
2394b9ad928SBarry Smith   pc->ops->view                = 0;
2404b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
2414b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
2424b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
2434b9ad928SBarry Smith   PetscFunctionReturn(0);
2444b9ad928SBarry Smith }
2454b9ad928SBarry Smith EXTERN_C_END
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith 
248