/* Include files needed for the PBJacobi preconditioner: pcimpl.h - private include file intended for use by all preconditioners */ #include #include /*I "petscpc.h" I*/ /* Private context (data structure) for the PBJacobi preconditioner. */ typedef struct { const MatScalar *diag; PetscInt bs,mbs; } PC_PBJacobi; #undef __FUNCT__ #define __FUNCT__ "PCApply_PBJacobi_1" static PetscErrorCode PCApply_PBJacobi_1(PC pc,Vec x,Vec y) { PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; const PetscScalar *xx; PetscScalar *yy; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,x2,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,x2,x3,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,x2,x3,x4,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,x2,x3,x4,x5,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,m = jac->mbs; const MatScalar *diag = jac->diag; PetscScalar x0,x1,x2,x3,x4,x5,x6,*yy; const PetscScalar *xx; PetscFunctionBegin; ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; Mat A = pc->pmat; PetscFunctionBegin; if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supported only for square matrices and square storage"); ierr = MatInvertBlockDiagonal(A,&jac->diag);CHKERRQ(ierr); ierr = MatGetBlockSize(A,&jac->bs);CHKERRQ(ierr); jac->mbs = A->rmap->n/jac->bs; switch (jac->bs) { case 1: pc->ops->apply = PCApply_PBJacobi_1; break; 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; case 6: pc->ops->apply = PCApply_PBJacobi_6; break; case 7: pc->ops->apply = PCApply_PBJacobi_7; break; default: SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"not supported for block size %D",jac->bs); } PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ #undef __FUNCT__ #define __FUNCT__ "PCDestroy_PBJacobi" static PetscErrorCode PCDestroy_PBJacobi(PC pc) { PetscErrorCode ierr; PetscFunctionBegin; /* Free the private data structure that was hanging off the PC */ ierr = PetscFree(pc->data);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "PCView_PBJacobi" static PetscErrorCode PCView_PBJacobi(PC pc,PetscViewer viewer) { PetscErrorCode ierr; PC_PBJacobi *jac = (PC_PBJacobi*)pc->data; PetscBool iascii; PetscFunctionBegin; ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); if (iascii) { ierr = PetscViewerASCIIPrintf(viewer," point-block Jacobi: block size %D\n",jac->bs);CHKERRQ(ierr); } PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ /*MC PCPBJACOBI - Point block Jacobi preconditioner Notes: See PCJACOBI for point Jacobi preconditioning This works for AIJ and BAIJ matrices and uses the blocksize provided to the matrix Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot is detected a PETSc error is generated. Developer Notes: This should support the PCSetErrorIfFailure() flag set to PETSC_TRUE to allow the factorization to continue even after a zero pivot is found resulting in a Nan and hence terminating KSP with a KSP_DIVERGED_NANORIF allowing a nonlinear solver/ODE integrator to recover without stopping the program as currently happens. Developer Note: Perhaps should provide an option that allows generation of a valid preconditioner even if a block is singular as the PCJACOBI does. Level: beginner Concepts: point block Jacobi .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI M*/ #undef __FUNCT__ #define __FUNCT__ "PCCreate_PBJacobi" PETSC_EXTERN 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 = PetscNewLog(pc,&jac);CHKERRQ(ierr); pc->data = (void*)jac; /* 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 = PCView_PBJacobi; pc->ops->applyrichardson = 0; pc->ops->applysymmetricleft = 0; pc->ops->applysymmetricright = 0; PetscFunctionReturn(0); }