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*c6db04a5SJed 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 */ 21*c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 22*c6db04a5SJed Brown #include <../src/mat/blockinvert.h> 234b9ad928SBarry Smith 244b9ad928SBarry Smith #undef __FUNCT__ 254b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_2" 266849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y) 274b9ad928SBarry Smith { 284b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 29dfbe8321SBarry Smith PetscErrorCode ierr; 30c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 3185f4f44aSBarry Smith const MatScalar *diag = jac->diag; 3285f4f44aSBarry Smith PetscScalar 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); 45dc0b31edSSatish Balay ierr = PetscLogFlops(6.0*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; 5585f4f44aSBarry Smith const MatScalar *diag = jac->diag; 5685f4f44aSBarry Smith PetscScalar x0,x1,x2,*xx,*yy; 574b9ad928SBarry Smith 584b9ad928SBarry Smith PetscFunctionBegin; 594b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 604b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 614b9ad928SBarry Smith for (i=0; i<m; i++) { 624b9ad928SBarry Smith x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2]; 634b9ad928SBarry Smith yy[3*i] = diag[0]*x0 + diag[3]*x1 + diag[6]*x2; 644b9ad928SBarry Smith yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2; 654b9ad928SBarry Smith yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2; 664b9ad928SBarry Smith diag += 9; 674b9ad928SBarry Smith } 684b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 694b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 70dc0b31edSSatish Balay ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr); 714b9ad928SBarry Smith PetscFunctionReturn(0); 724b9ad928SBarry Smith } 734b9ad928SBarry Smith #undef __FUNCT__ 744b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4" 756849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y) 764b9ad928SBarry Smith { 774b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 78dfbe8321SBarry Smith PetscErrorCode ierr; 79c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 8085f4f44aSBarry Smith const MatScalar *diag = jac->diag; 8185f4f44aSBarry Smith PetscScalar x0,x1,x2,x3,*xx,*yy; 824b9ad928SBarry Smith 834b9ad928SBarry Smith PetscFunctionBegin; 844b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 854b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 864b9ad928SBarry Smith for (i=0; i<m; i++) { 874b9ad928SBarry Smith x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3]; 884b9ad928SBarry Smith yy[4*i] = diag[0]*x0 + diag[4]*x1 + diag[8]*x2 + diag[12]*x3; 894b9ad928SBarry Smith yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2 + diag[13]*x3; 904b9ad928SBarry Smith yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3; 914b9ad928SBarry Smith yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3; 924b9ad928SBarry Smith diag += 16; 934b9ad928SBarry Smith } 944b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 954b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 96dc0b31edSSatish Balay ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr); 974b9ad928SBarry Smith PetscFunctionReturn(0); 984b9ad928SBarry Smith } 994b9ad928SBarry Smith #undef __FUNCT__ 1004b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5" 1016849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y) 1024b9ad928SBarry Smith { 1034b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 104dfbe8321SBarry Smith PetscErrorCode ierr; 105c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 10685f4f44aSBarry Smith const MatScalar *diag = jac->diag; 10785f4f44aSBarry Smith PetscScalar x0,x1,x2,x3,x4,*xx,*yy; 1084b9ad928SBarry Smith 1094b9ad928SBarry Smith PetscFunctionBegin; 1104b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1114b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1124b9ad928SBarry Smith for (i=0; i<m; i++) { 1134b9ad928SBarry 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]; 1144b9ad928SBarry Smith yy[5*i] = diag[0]*x0 + diag[5]*x1 + diag[10]*x2 + diag[15]*x3 + diag[20]*x4; 1154b9ad928SBarry Smith yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2 + diag[16]*x3 + diag[21]*x4; 1164b9ad928SBarry Smith yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4; 1174b9ad928SBarry Smith yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4; 1184b9ad928SBarry Smith yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4; 1194b9ad928SBarry Smith diag += 25; 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1224b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 123dc0b31edSSatish Balay ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr); 1244b9ad928SBarry Smith PetscFunctionReturn(0); 1254b9ad928SBarry Smith } 1264b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1274b9ad928SBarry Smith #undef __FUNCT__ 1284b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi" 1296849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc) 1304b9ad928SBarry Smith { 1314b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 132dfbe8321SBarry Smith PetscErrorCode ierr; 133c1ac3661SBarry Smith PetscMPIInt size; 134ace3abfcSBarry Smith PetscBool seqbaij,mpibaij,baij; 1354b9ad928SBarry Smith Mat A = pc->pmat; 1364b9ad928SBarry Smith Mat_SeqBAIJ *a; 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith PetscFunctionBegin; 1394b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);CHKERRQ(ierr); 1404b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);CHKERRQ(ierr); 14121360622SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);CHKERRQ(ierr); 142e7e72b3dSBarry Smith if (!seqbaij && !mpibaij && !baij) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only supports BAIJ matrices"); 1437adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 14421360622SBarry Smith if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A; 145e32f2f54SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage"); 146a1d92eedSBarry Smith 14743516a2dSKris Buschelman ierr = MatSeqBAIJInvertBlockDiagonal(A);CHKERRQ(ierr); 1484b9ad928SBarry Smith a = (Mat_SeqBAIJ*)A->data; 149a1d92eedSBarry Smith jac->diag = a->idiag; 150d0f46423SBarry Smith jac->bs = A->rmap->bs; 1514b9ad928SBarry Smith jac->mbs = a->mbs; 152521d7252SBarry Smith switch (jac->bs){ 1534b9ad928SBarry Smith case 2: 1544b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_2; 1554b9ad928SBarry Smith break; 1564b9ad928SBarry Smith case 3: 1574b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_3; 1584b9ad928SBarry Smith break; 1594b9ad928SBarry Smith case 4: 1604b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_4; 1614b9ad928SBarry Smith break; 1624b9ad928SBarry Smith case 5: 1634b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_5; 1644b9ad928SBarry Smith break; 1654b9ad928SBarry Smith default: 16665e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"not supported for block size %D",jac->bs); 1674b9ad928SBarry Smith } 1684b9ad928SBarry Smith 1694b9ad928SBarry Smith PetscFunctionReturn(0); 1704b9ad928SBarry Smith } 1714b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1724b9ad928SBarry Smith #undef __FUNCT__ 1734b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi" 1746849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc) 1754b9ad928SBarry Smith { 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 */ 182c31cb41cSBarry Smith ierr = PetscFree(pc->data);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" 2027087cfbeSBarry Smith PetscErrorCode 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 */ 21338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_PBJacobi,&jac);CHKERRQ(ierr); 2144b9ad928SBarry Smith pc->data = (void*)jac; 2154b9ad928SBarry Smith 2164b9ad928SBarry Smith /* 2174b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 2184b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 2194b9ad928SBarry Smith */ 2204b9ad928SBarry Smith jac->diag = 0; 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith /* 2234b9ad928SBarry Smith Set the pointers for the functions that are provided above. 2244b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 2254b9ad928SBarry Smith are called, they will automatically call these functions. Note we 2264b9ad928SBarry Smith choose not to provide a couple of these functions since they are 2274b9ad928SBarry Smith not needed. 2284b9ad928SBarry Smith */ 2294b9ad928SBarry Smith pc->ops->apply = 0; /*set depending on the block size */ 2304b9ad928SBarry Smith pc->ops->applytranspose = 0; 2314b9ad928SBarry Smith pc->ops->setup = PCSetUp_PBJacobi; 2324b9ad928SBarry Smith pc->ops->destroy = PCDestroy_PBJacobi; 2334b9ad928SBarry Smith pc->ops->setfromoptions = 0; 2344b9ad928SBarry Smith pc->ops->view = 0; 2354b9ad928SBarry Smith pc->ops->applyrichardson = 0; 2364b9ad928SBarry Smith pc->ops->applysymmetricleft = 0; 2374b9ad928SBarry Smith pc->ops->applysymmetricright = 0; 2384b9ad928SBarry Smith PetscFunctionReturn(0); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith EXTERN_C_END 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith 243