/* Include files needed for the PBJacobi preconditioner: pcimpl.h - private include file intended for use by all preconditioners */ #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ /* Private context (data structure) for the PBJacobi preconditioner. */ typedef struct { PetscScalar *diag; int bs,mbs; } PC_PBJacobi; /* Currently only implemented for baij matrices and directly access baij data structures. */ #include "src/mat/impls/baij/mpi/mpibaij.h" #include "src/inline/ilu.h" #undef __FUNCT__ #define __FUNCT__ "PCApply_PBJacobi_2" static PetscErrorCode PCApply_PBJacobi_2(PC pc,Vec x,Vec y) { PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; PetscErrorCode ierr; int i,m = jac->mbs; PetscScalar *diag = jac->diag,x0,x1,*xx,*yy; PetscFunctionBegin; ierr = VecGetArray(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; int i,m = jac->mbs; PetscScalar *diag = jac->diag,x0,x1,x2,*xx,*yy; PetscFunctionBegin; ierr = VecGetArray(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; int i,m = jac->mbs; PetscScalar *diag = jac->diag,x0,x1,x2,x3,*xx,*yy; PetscFunctionBegin; ierr = VecGetArray(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; int i,m = jac->mbs; PetscScalar *diag = jac->diag,x0,x1,x2,x3,x4,*xx,*yy; PetscFunctionBegin; ierr = VecGetArray(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; int size; PetscTruth seqbaij,mpibaij,baij; Mat A = pc->pmat; Mat_SeqBAIJ *a; PetscFunctionBegin; ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQBAIJ,&seqbaij);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIBAIJ,&mpibaij);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)pc->pmat,MATBAIJ,&baij);CHKERRQ(ierr); if (!seqbaij && !mpibaij && !baij) { SETERRQ(1,"Currently only supports BAIJ matrices"); } ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); if (mpibaij || (baij && (size > 1))) A = ((Mat_MPIBAIJ*)A->data)->A; if (A->m != A->n) SETERRQ(1,"Supported only for square matrices and square storage"); ierr = MatInvertBlockDiagonal_SeqBAIJ(A);CHKERRQ(ierr); a = (Mat_SeqBAIJ*)A->data; jac->diag = a->idiag; jac->bs = a->bs; jac->mbs = a->mbs; switch (a->bs){ case 2: pc->ops->apply = PCApply_PBJacobi_2; break; case 3: pc->ops->apply = PCApply_PBJacobi_3; break; case 4: pc->ops->apply = PCApply_PBJacobi_4; break; case 5: pc->ops->apply = PCApply_PBJacobi_5; break; default: SETERRQ1(1,"not supported for block size %d",a->bs); } PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "PCDestroy_PBJacobi" static PetscErrorCode PCDestroy_PBJacobi(PC pc) { PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; PetscErrorCode ierr; PetscFunctionBegin; /* Free the private data structure that was hanging off the PC */ ierr = PetscFree(jac);CHKERRQ(ierr); PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ /*MC PCPBJACOBI - Point block Jacobi Level: beginner Concepts: point block Jacobi Notes: Only implemented for the BAIJ matrix formats. .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC M*/ EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "PCCreate_PBJacobi" PetscErrorCode PCCreate_PBJacobi(PC pc) { PC_PBJacobi *jac; PetscErrorCode ierr; PetscFunctionBegin; /* Creates the private data structure for this preconditioner and attach it to the PC object. */ ierr = PetscNew(PC_PBJacobi,&jac);CHKERRQ(ierr); pc->data = (void*)jac; /* Logs the memory usage; this is not needed but allows PETSc to monitor how much memory is being used for various purposes. */ PetscLogObjectMemory(pc,sizeof(PC_PBJacobi)); /* Initialize the pointers to vectors to ZERO; these will be used to store diagonal entries of the matrix for fast preconditioner application. */ jac->diag = 0; /* Set the pointers for the functions that are provided above. Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) are called, they will automatically call these functions. Note we choose not to provide a couple of these functions since they are not needed. */ pc->ops->apply = 0; /*set depending on the block size */ pc->ops->applytranspose = 0; pc->ops->setup = PCSetUp_PBJacobi; pc->ops->destroy = PCDestroy_PBJacobi; pc->ops->setfromoptions = 0; pc->ops->view = 0; pc->ops->applyrichardson = 0; pc->ops->applysymmetricleft = 0; pc->ops->applysymmetricright = 0; PetscFunctionReturn(0); } EXTERN_C_END