1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris 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 86356e834SBarry Smith #include "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 { 14*85f4f44aSBarry Smith MatScalar *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 */ 227c4f633dSBarry Smith #include "../src/mat/impls/baij/mpi/mpibaij.h" 23c60f0209SBarry Smith #include "../src/mat/blockinvert.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; 32*85f4f44aSBarry Smith const MatScalar *diag = jac->diag; 33*85f4f44aSBarry Smith PetscScalar x0,x1,*xx,*yy; 344b9ad928SBarry Smith 354b9ad928SBarry Smith PetscFunctionBegin; 364b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 374b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 384b9ad928SBarry Smith for (i=0; i<m; i++) { 394b9ad928SBarry Smith x0 = xx[2*i]; x1 = xx[2*i+1]; 404b9ad928SBarry Smith yy[2*i] = diag[0]*x0 + diag[2]*x1; 414b9ad928SBarry Smith yy[2*i+1] = diag[1]*x0 + diag[3]*x1; 424b9ad928SBarry Smith diag += 4; 434b9ad928SBarry Smith } 444b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 454b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 46dc0b31edSSatish Balay ierr = PetscLogFlops(6.0*m);CHKERRQ(ierr); 474b9ad928SBarry Smith PetscFunctionReturn(0); 484b9ad928SBarry Smith } 494b9ad928SBarry Smith #undef __FUNCT__ 504b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_3" 516849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_3(PC pc,Vec x,Vec y) 524b9ad928SBarry Smith { 534b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 54dfbe8321SBarry Smith PetscErrorCode ierr; 55c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 56*85f4f44aSBarry Smith const MatScalar *diag = jac->diag; 57*85f4f44aSBarry Smith PetscScalar x0,x1,x2,*xx,*yy; 584b9ad928SBarry Smith 594b9ad928SBarry Smith PetscFunctionBegin; 604b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 614b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 624b9ad928SBarry Smith for (i=0; i<m; i++) { 634b9ad928SBarry Smith x0 = xx[3*i]; x1 = xx[3*i+1]; x2 = xx[3*i+2]; 644b9ad928SBarry Smith yy[3*i] = diag[0]*x0 + diag[3]*x1 + diag[6]*x2; 654b9ad928SBarry Smith yy[3*i+1] = diag[1]*x0 + diag[4]*x1 + diag[7]*x2; 664b9ad928SBarry Smith yy[3*i+2] = diag[2]*x0 + diag[5]*x1 + diag[8]*x2; 674b9ad928SBarry Smith diag += 9; 684b9ad928SBarry Smith } 694b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 704b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 71dc0b31edSSatish Balay ierr = PetscLogFlops(15.0*m);CHKERRQ(ierr); 724b9ad928SBarry Smith PetscFunctionReturn(0); 734b9ad928SBarry Smith } 744b9ad928SBarry Smith #undef __FUNCT__ 754b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_4" 766849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_4(PC pc,Vec x,Vec y) 774b9ad928SBarry Smith { 784b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 79dfbe8321SBarry Smith PetscErrorCode ierr; 80c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 81*85f4f44aSBarry Smith const MatScalar *diag = jac->diag; 82*85f4f44aSBarry Smith PetscScalar x0,x1,x2,x3,*xx,*yy; 834b9ad928SBarry Smith 844b9ad928SBarry Smith PetscFunctionBegin; 854b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 864b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 874b9ad928SBarry Smith for (i=0; i<m; i++) { 884b9ad928SBarry Smith x0 = xx[4*i]; x1 = xx[4*i+1]; x2 = xx[4*i+2]; x3 = xx[4*i+3]; 894b9ad928SBarry Smith yy[4*i] = diag[0]*x0 + diag[4]*x1 + diag[8]*x2 + diag[12]*x3; 904b9ad928SBarry Smith yy[4*i+1] = diag[1]*x0 + diag[5]*x1 + diag[9]*x2 + diag[13]*x3; 914b9ad928SBarry Smith yy[4*i+2] = diag[2]*x0 + diag[6]*x1 + diag[10]*x2 + diag[14]*x3; 924b9ad928SBarry Smith yy[4*i+3] = diag[3]*x0 + diag[7]*x1 + diag[11]*x2 + diag[15]*x3; 934b9ad928SBarry Smith diag += 16; 944b9ad928SBarry Smith } 954b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 964b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 97dc0b31edSSatish Balay ierr = PetscLogFlops(28.0*m);CHKERRQ(ierr); 984b9ad928SBarry Smith PetscFunctionReturn(0); 994b9ad928SBarry Smith } 1004b9ad928SBarry Smith #undef __FUNCT__ 1014b9ad928SBarry Smith #define __FUNCT__ "PCApply_PBJacobi_5" 1026849ba73SBarry Smith static PetscErrorCode PCApply_PBJacobi_5(PC pc,Vec x,Vec y) 1034b9ad928SBarry Smith { 1044b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 105dfbe8321SBarry Smith PetscErrorCode ierr; 106c1ac3661SBarry Smith PetscInt i,m = jac->mbs; 107*85f4f44aSBarry Smith const MatScalar *diag = jac->diag; 108*85f4f44aSBarry Smith PetscScalar x0,x1,x2,x3,x4,*xx,*yy; 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith PetscFunctionBegin; 1114b9ad928SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1124b9ad928SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 1134b9ad928SBarry Smith for (i=0; i<m; i++) { 1144b9ad928SBarry 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]; 1154b9ad928SBarry Smith yy[5*i] = diag[0]*x0 + diag[5]*x1 + diag[10]*x2 + diag[15]*x3 + diag[20]*x4; 1164b9ad928SBarry Smith yy[5*i+1] = diag[1]*x0 + diag[6]*x1 + diag[11]*x2 + diag[16]*x3 + diag[21]*x4; 1174b9ad928SBarry Smith yy[5*i+2] = diag[2]*x0 + diag[7]*x1 + diag[12]*x2 + diag[17]*x3 + diag[22]*x4; 1184b9ad928SBarry Smith yy[5*i+3] = diag[3]*x0 + diag[8]*x1 + diag[13]*x2 + diag[18]*x3 + diag[23]*x4; 1194b9ad928SBarry Smith yy[5*i+4] = diag[4]*x0 + diag[9]*x1 + diag[14]*x2 + diag[19]*x3 + diag[24]*x4; 1204b9ad928SBarry Smith diag += 25; 1214b9ad928SBarry Smith } 1224b9ad928SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1234b9ad928SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 124dc0b31edSSatish Balay ierr = PetscLogFlops(45.0*m);CHKERRQ(ierr); 1254b9ad928SBarry Smith PetscFunctionReturn(0); 1264b9ad928SBarry Smith } 1274b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1284b9ad928SBarry Smith #undef __FUNCT__ 1294b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_PBJacobi" 1306849ba73SBarry Smith static PetscErrorCode PCSetUp_PBJacobi(PC pc) 1314b9ad928SBarry Smith { 1324b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 133dfbe8321SBarry Smith PetscErrorCode ierr; 134c1ac3661SBarry Smith PetscMPIInt size; 13521360622SBarry Smith PetscTruth seqbaij,mpibaij,baij; 1364b9ad928SBarry Smith Mat A = pc->pmat; 1374b9ad928SBarry Smith Mat_SeqBAIJ *a; 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith PetscFunctionBegin; 1404b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);CHKERRQ(ierr); 1414b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);CHKERRQ(ierr); 14221360622SBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);CHKERRQ(ierr); 143e7e72b3dSBarry Smith if (!seqbaij && !mpibaij && !baij) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only supports BAIJ matrices"); 1447adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 14521360622SBarry Smith if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A; 146e32f2f54SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage"); 147a1d92eedSBarry Smith 14843516a2dSKris Buschelman ierr = MatSeqBAIJInvertBlockDiagonal(A);CHKERRQ(ierr); 1494b9ad928SBarry Smith a = (Mat_SeqBAIJ*)A->data; 150a1d92eedSBarry Smith jac->diag = a->idiag; 151d0f46423SBarry Smith jac->bs = A->rmap->bs; 1524b9ad928SBarry Smith jac->mbs = a->mbs; 153521d7252SBarry Smith switch (jac->bs){ 1544b9ad928SBarry Smith case 2: 1554b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_2; 1564b9ad928SBarry Smith break; 1574b9ad928SBarry Smith case 3: 1584b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_3; 1594b9ad928SBarry Smith break; 1604b9ad928SBarry Smith case 4: 1614b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_4; 1624b9ad928SBarry Smith break; 1634b9ad928SBarry Smith case 5: 1644b9ad928SBarry Smith pc->ops->apply = PCApply_PBJacobi_5; 1654b9ad928SBarry Smith break; 1664b9ad928SBarry Smith default: 16765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"not supported for block size %D",jac->bs); 1684b9ad928SBarry Smith } 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith PetscFunctionReturn(0); 1714b9ad928SBarry Smith } 1724b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 1734b9ad928SBarry Smith #undef __FUNCT__ 1744b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_PBJacobi" 1756849ba73SBarry Smith static PetscErrorCode PCDestroy_PBJacobi(PC pc) 1764b9ad928SBarry Smith { 1774b9ad928SBarry Smith PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; 178dfbe8321SBarry Smith PetscErrorCode ierr; 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith PetscFunctionBegin; 1814b9ad928SBarry Smith /* 1824b9ad928SBarry Smith Free the private data structure that was hanging off the PC 1834b9ad928SBarry Smith */ 1844b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 1854b9ad928SBarry Smith PetscFunctionReturn(0); 1864b9ad928SBarry Smith } 1874b9ad928SBarry Smith /* -------------------------------------------------------------------------- */ 18837a17b4dSBarry Smith /*MC 18937a17b4dSBarry Smith PCPBJACOBI - Point block Jacobi 19037a17b4dSBarry Smith 19137a17b4dSBarry Smith Level: beginner 19237a17b4dSBarry Smith 19337a17b4dSBarry Smith Concepts: point block Jacobi 19437a17b4dSBarry Smith 19537a17b4dSBarry Smith Notes: Only implemented for the BAIJ matrix formats. 19637a17b4dSBarry Smith 19737a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 19837a17b4dSBarry Smith 19937a17b4dSBarry Smith M*/ 20037a17b4dSBarry Smith 2014b9ad928SBarry Smith EXTERN_C_BEGIN 2024b9ad928SBarry Smith #undef __FUNCT__ 2034b9ad928SBarry Smith #define __FUNCT__ "PCCreate_PBJacobi" 204dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PBJacobi(PC pc) 2054b9ad928SBarry Smith { 2064b9ad928SBarry Smith PC_PBJacobi *jac; 207dfbe8321SBarry Smith PetscErrorCode ierr; 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith PetscFunctionBegin; 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith /* 2124b9ad928SBarry Smith Creates the private data structure for this preconditioner and 2134b9ad928SBarry Smith attach it to the PC object. 2144b9ad928SBarry Smith */ 21538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_PBJacobi,&jac);CHKERRQ(ierr); 2164b9ad928SBarry Smith pc->data = (void*)jac; 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith /* 2194b9ad928SBarry Smith Initialize the pointers to vectors to ZERO; these will be used to store 2204b9ad928SBarry Smith diagonal entries of the matrix for fast preconditioner application. 2214b9ad928SBarry Smith */ 2224b9ad928SBarry Smith jac->diag = 0; 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith /* 2254b9ad928SBarry Smith Set the pointers for the functions that are provided above. 2264b9ad928SBarry Smith Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 2274b9ad928SBarry Smith are called, they will automatically call these functions. Note we 2284b9ad928SBarry Smith choose not to provide a couple of these functions since they are 2294b9ad928SBarry Smith not needed. 2304b9ad928SBarry Smith */ 2314b9ad928SBarry Smith pc->ops->apply = 0; /*set depending on the block size */ 2324b9ad928SBarry Smith pc->ops->applytranspose = 0; 2334b9ad928SBarry Smith pc->ops->setup = PCSetUp_PBJacobi; 2344b9ad928SBarry Smith pc->ops->destroy = PCDestroy_PBJacobi; 2354b9ad928SBarry Smith pc->ops->setfromoptions = 0; 2364b9ad928SBarry Smith pc->ops->view = 0; 2374b9ad928SBarry Smith pc->ops->applyrichardson = 0; 2384b9ad928SBarry Smith pc->ops->applysymmetricleft = 0; 2394b9ad928SBarry Smith pc->ops->applysymmetricright = 0; 2404b9ad928SBarry Smith PetscFunctionReturn(0); 2414b9ad928SBarry Smith } 2424b9ad928SBarry Smith EXTERN_C_END 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith 245