xref: /petsc/src/ksp/pc/impls/pbjacobi/pbjacobi.c (revision dba47a550923b04c7c4ebbb735eb62a1b3e4e9ae)
1*dba47a55SKris Buschelman #define PETSCKSP_DLL
2*dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith    Include files needed for the PBJacobi preconditioner:
54b9ad928SBarry Smith      pcimpl.h - private include file intended for use by all preconditioners
64b9ad928SBarry Smith */
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #include "src/ksp/pc/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 {
144b9ad928SBarry Smith   PetscScalar *diag;
15c1ac3661SBarry Smith   PetscInt    bs,mbs;
164b9ad928SBarry Smith } PC_PBJacobi;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith /*
194b9ad928SBarry Smith    Currently only implemented for baij matrices and directly access baij
204b9ad928SBarry Smith   data structures.
214b9ad928SBarry Smith */
224b9ad928SBarry Smith #include "src/mat/impls/baij/mpi/mpibaij.h"
234b9ad928SBarry Smith #include "src/inline/ilu.h"
244b9ad928SBarry Smith 
254b9ad928SBarry Smith #undef __FUNCT__
264b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2"
276849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y)
284b9ad928SBarry Smith {
294b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
30dfbe8321SBarry Smith   PetscErrorCode ierr;
31c1ac3661SBarry Smith   PetscInt       i,m = jac->mbs;
324b9ad928SBarry Smith   PetscScalar    *diag = jac->diag,x0,x1,*xx,*yy;
334b9ad928SBarry Smith 
344b9ad928SBarry Smith   PetscFunctionBegin;
354b9ad928SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
364b9ad928SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
374b9ad928SBarry Smith   for (i=0; i<m; i++) {
384b9ad928SBarry Smith     x0 = xx[2*i]; x1 = xx[2*i+1];
394b9ad928SBarry Smith     yy[2*i]   = diag[0]*x0 + diag[2]*x1;
404b9ad928SBarry Smith     yy[2*i+1] = diag[1]*x0 + diag[3]*x1;
414b9ad928SBarry Smith     diag     += 4;
424b9ad928SBarry Smith   }
434b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
444b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
45efee365bSSatish Balay   ierr = PetscLogFlops(6*m);CHKERRQ(ierr);
464b9ad928SBarry Smith   PetscFunctionReturn(0);
474b9ad928SBarry Smith }
484b9ad928SBarry Smith #undef __FUNCT__
494b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3"
506849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y)
514b9ad928SBarry Smith {
524b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
54c1ac3661SBarry Smith   PetscInt       i,m = jac->mbs;
554b9ad928SBarry Smith   PetscScalar    *diag = jac->diag,x0,x1,x2,*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[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2];
624b9ad928SBarry Smith     yy[3*i]   = diag[0]*x0 + diag[3]*x1 + diag[6]*x2;
634b9ad928SBarry Smith     yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2;
644b9ad928SBarry Smith     yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2;
654b9ad928SBarry Smith     diag     += 9;
664b9ad928SBarry Smith   }
674b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
684b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
69efee365bSSatish Balay   ierr = PetscLogFlops(15*m);CHKERRQ(ierr);
704b9ad928SBarry Smith   PetscFunctionReturn(0);
714b9ad928SBarry Smith }
724b9ad928SBarry Smith #undef __FUNCT__
734b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4"
746849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y)
754b9ad928SBarry Smith {
764b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
77dfbe8321SBarry Smith   PetscErrorCode ierr;
78c1ac3661SBarry Smith   PetscInt       i,m = jac->mbs;
794b9ad928SBarry Smith   PetscScalar    *diag = jac->diag,x0,x1,x2,x3,*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[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3];
864b9ad928SBarry Smith     yy[4*i]   = diag[0]*x0 + diag[4]*x1 + diag[8]*x2  + diag[12]*x3;
874b9ad928SBarry Smith     yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2  + diag[13]*x3;
884b9ad928SBarry Smith     yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3;
894b9ad928SBarry Smith     yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3;
904b9ad928SBarry Smith     diag     += 16;
914b9ad928SBarry Smith   }
924b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
934b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
94efee365bSSatish Balay   ierr = PetscLogFlops(28*m);CHKERRQ(ierr);
954b9ad928SBarry Smith   PetscFunctionReturn(0);
964b9ad928SBarry Smith }
974b9ad928SBarry Smith #undef __FUNCT__
984b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5"
996849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y)
1004b9ad928SBarry Smith {
1014b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
102dfbe8321SBarry Smith   PetscErrorCode ierr;
103c1ac3661SBarry Smith   PetscInt       i,m = jac->mbs;
1044b9ad928SBarry Smith   PetscScalar    *diag = jac->diag,x0,x1,x2,x3,x4,*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[5*i]; x1 = xx[5*i+1]; x2 = xx[5*i+2]; x3 = xx[5*i+3]; x4 = xx[5*i+4];
1114b9ad928SBarry Smith     yy[5*i]   = diag[0]*x0 + diag[5]*x1 + diag[10]*x2  + diag[15]*x3 + diag[20]*x4;
1124b9ad928SBarry Smith     yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2  + diag[16]*x3 + diag[21]*x4;
1134b9ad928SBarry Smith     yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4;
1144b9ad928SBarry Smith     yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4;
1154b9ad928SBarry Smith     yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4;
1164b9ad928SBarry Smith     diag     += 25;
1174b9ad928SBarry Smith   }
1184b9ad928SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1194b9ad928SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
120efee365bSSatish Balay   ierr = PetscLogFlops(45*m);CHKERRQ(ierr);
1214b9ad928SBarry Smith   PetscFunctionReturn(0);
1224b9ad928SBarry Smith }
1234b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1244b9ad928SBarry Smith #undef __FUNCT__
1254b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi"
1266849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc)
1274b9ad928SBarry Smith {
1284b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
129dfbe8321SBarry Smith   PetscErrorCode ierr;
130c1ac3661SBarry Smith   PetscMPIInt    size;
13121360622SBarry Smith   PetscTruth     seqbaij,mpibaij,baij;
1324b9ad928SBarry Smith   Mat            A = pc->pmat;
1334b9ad928SBarry Smith   Mat_SeqBAIJ    *a;
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith   PetscFunctionBegin;
1364b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);CHKERRQ(ierr);
1374b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);CHKERRQ(ierr);
13821360622SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);CHKERRQ(ierr);
13921360622SBarry Smith   if (!seqbaij && !mpibaij && !baij) {
1401302d50aSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Currently only supports BAIJ matrices");
1414b9ad928SBarry Smith   }
14221360622SBarry Smith   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
14321360622SBarry Smith   if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A;
1441302d50aSBarry Smith   if (A->m != A->n) SETERRQ(PETSC_ERR_SUP,"Supported only for square matrices and square storage");
145a1d92eedSBarry Smith 
14643516a2dSKris Buschelman   ierr        =  MatSeqBAIJInvertBlockDiagonal(A);CHKERRQ(ierr);
1474b9ad928SBarry Smith   a           = (Mat_SeqBAIJ*)A->data;
148a1d92eedSBarry Smith   jac->diag   = a->idiag;
149521d7252SBarry Smith   jac->bs     = A->bs;
1504b9ad928SBarry Smith   jac->mbs    = a->mbs;
151521d7252SBarry Smith   switch (jac->bs){
1524b9ad928SBarry Smith     case 2:
1534b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_2;
1544b9ad928SBarry Smith       break;
1554b9ad928SBarry Smith     case 3:
1564b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_3;
1574b9ad928SBarry Smith       break;
1584b9ad928SBarry Smith     case 4:
1594b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_4;
1604b9ad928SBarry Smith       break;
1614b9ad928SBarry Smith     case 5:
1624b9ad928SBarry Smith       pc->ops->apply = PCApply_PBJacobi_5;
1634b9ad928SBarry Smith       break;
1644b9ad928SBarry Smith     default:
165521d7252SBarry Smith       SETERRQ1(PETSC_ERR_SUP,"not supported for block size %D",jac->bs);
1664b9ad928SBarry Smith   }
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith   PetscFunctionReturn(0);
1694b9ad928SBarry Smith }
1704b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
1714b9ad928SBarry Smith #undef __FUNCT__
1724b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi"
1736849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc)
1744b9ad928SBarry Smith {
1754b9ad928SBarry Smith   PC_PBJacobi    *jac = (PC_PBJacobi*)pc->data;
176dfbe8321SBarry Smith   PetscErrorCode ierr;
1774b9ad928SBarry Smith 
1784b9ad928SBarry Smith   PetscFunctionBegin;
1794b9ad928SBarry Smith   /*
1804b9ad928SBarry Smith       Free the private data structure that was hanging off the PC
1814b9ad928SBarry Smith   */
1824b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1834b9ad928SBarry Smith   PetscFunctionReturn(0);
1844b9ad928SBarry Smith }
1854b9ad928SBarry Smith /* -------------------------------------------------------------------------- */
18637a17b4dSBarry Smith /*MC
18737a17b4dSBarry Smith      PCPBJACOBI - Point block Jacobi
18837a17b4dSBarry Smith 
18937a17b4dSBarry Smith    Level: beginner
19037a17b4dSBarry Smith 
19137a17b4dSBarry Smith   Concepts: point block Jacobi
19237a17b4dSBarry Smith 
19337a17b4dSBarry Smith    Notes: Only implemented for the BAIJ matrix formats.
19437a17b4dSBarry Smith 
19537a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
19637a17b4dSBarry Smith 
19737a17b4dSBarry Smith M*/
19837a17b4dSBarry Smith 
1994b9ad928SBarry Smith EXTERN_C_BEGIN
2004b9ad928SBarry Smith #undef __FUNCT__
2014b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi"
202*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PBJacobi(PC pc)
2034b9ad928SBarry Smith {
2044b9ad928SBarry Smith   PC_PBJacobi    *jac;
205dfbe8321SBarry Smith   PetscErrorCode ierr;
2064b9ad928SBarry Smith 
2074b9ad928SBarry Smith   PetscFunctionBegin;
2084b9ad928SBarry Smith 
2094b9ad928SBarry Smith   /*
2104b9ad928SBarry Smith      Creates the private data structure for this preconditioner and
2114b9ad928SBarry Smith      attach it to the PC object.
2124b9ad928SBarry Smith   */
2134b9ad928SBarry Smith   ierr      = PetscNew(PC_PBJacobi,&jac);CHKERRQ(ierr);
2144b9ad928SBarry Smith   pc->data  = (void*)jac;
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith   /*
2174b9ad928SBarry Smith      Logs the memory usage; this is not needed but allows PETSc to
2184b9ad928SBarry Smith      monitor how much memory is being used for various purposes.
2194b9ad928SBarry Smith   */
22052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_PBJacobi));CHKERRQ(ierr);
2214b9ad928SBarry Smith 
2224b9ad928SBarry Smith   /*
2234b9ad928SBarry Smith      Initialize the pointers to vectors to ZERO; these will be used to store
2244b9ad928SBarry Smith      diagonal entries of the matrix for fast preconditioner application.
2254b9ad928SBarry Smith   */
2264b9ad928SBarry Smith   jac->diag          = 0;
2274b9ad928SBarry Smith 
2284b9ad928SBarry Smith   /*
2294b9ad928SBarry Smith       Set the pointers for the functions that are provided above.
2304b9ad928SBarry Smith       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
2314b9ad928SBarry Smith       are called, they will automatically call these functions.  Note we
2324b9ad928SBarry Smith       choose not to provide a couple of these functions since they are
2334b9ad928SBarry Smith       not needed.
2344b9ad928SBarry Smith   */
2354b9ad928SBarry Smith   pc->ops->apply               = 0; /*set depending on the block size */
2364b9ad928SBarry Smith   pc->ops->applytranspose      = 0;
2374b9ad928SBarry Smith   pc->ops->setup               = PCSetUp_PBJacobi;
2384b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_PBJacobi;
2394b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
2404b9ad928SBarry Smith   pc->ops->view                = 0;
2414b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
2424b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
2434b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
2444b9ad928SBarry Smith   PetscFunctionReturn(0);
2454b9ad928SBarry Smith }
2464b9ad928SBarry Smith EXTERN_C_END
2474b9ad928SBarry Smith 
2484b9ad928SBarry Smith 
249