/* Include files needed for the variable size block PBJacobi preconditioner: pcimpl.h - private include file intended for use by all preconditioners */ #include /*I "petscpc.h" I*/ /* Private context (data structure) for the VPBJacobi preconditioner. */ typedef struct { MatScalar *diag; } PC_VPBJacobi; static PetscErrorCode PCApply_VPBJacobi(PC pc,Vec x,Vec y) { PC_VPBJacobi *jac = (PC_VPBJacobi*)pc->data; PetscErrorCode ierr; PetscInt i,ncnt = 0; const MatScalar *diag = jac->diag; PetscInt ib,jb,bs; const PetscScalar *xx; PetscScalar *yy,x0,x1,x2,x3,x4,x5,x6; PetscInt nblocks; const PetscInt *bsizes; PetscFunctionBegin; ierr = MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes);CHKERRQ(ierr); ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); ierr = VecGetArray(y,&yy);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; Mat A = pc->pmat; MatFactorError err; PetscInt i,nsize = 0,nlocal; PetscInt nblocks; const PetscInt *bsizes; PetscFunctionBegin; ierr = MatGetVariableBlockSizes(pc->pmat,&nblocks,&bsizes);CHKERRQ(ierr); ierr = MatGetLocalSize(pc->pmat,&nlocal,NULL);CHKERRQ(ierr); if (nlocal && !nblocks) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatSetVariableBlockSizes() before using PCVPBJACOBI"); if (!jac->diag) { for (i=0; idiag);CHKERRQ(ierr); } ierr = MatInvertVariableBlockDiagonal(A,nblocks,bsizes,jac->diag);CHKERRQ(ierr); ierr = MatFactorGetError(A,&err);CHKERRQ(ierr); if (err) pc->failedreason = (PCFailedReason)err; pc->ops->apply = PCApply_VPBJacobi; PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ static PetscErrorCode PCDestroy_VPBJacobi(PC pc) { PC_VPBJacobi *jac = (PC_VPBJacobi*)pc->data; PetscErrorCode ierr; PetscFunctionBegin; /* Free the private data structure that was hanging off the PC */ ierr = PetscFree(jac->diag);CHKERRQ(ierr); ierr = PetscFree(pc->data);CHKERRQ(ierr); PetscFunctionReturn(0); } /* -------------------------------------------------------------------------- */ /*MC PCVPBJACOBI - Variable size point block Jacobi preconditioner Notes: See PCJACOBI for point Jacobi preconditioning, PCPBJACOBI for fixed point block size, and PCBJACOBI for large size blocks This works for AIJ matrices Uses dense LU factorization with partial pivoting to invert the blocks; if a zero pivot is detected a PETSc error is generated. One must call MatSetVariableBlockSizes() to use this preconditioner 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. Perhaps should provide an option that allows generation of a valid preconditioner even if a block is singular as the PCJACOBI does. Level: beginner .seealso: MatSetVariableBlockSizes(), PCCreate(), PCSetType(), PCType (for list of available types), PC, PCJACOBI M*/ PETSC_EXTERN PetscErrorCode PCCreate_VPBJacobi(PC pc) { PC_VPBJacobi *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 = NULL; /* 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 = PCApply_VPBJacobi; pc->ops->applytranspose = NULL; pc->ops->setup = PCSetUp_VPBJacobi; pc->ops->destroy = PCDestroy_VPBJacobi; pc->ops->setfromoptions = NULL; pc->ops->applyrichardson = NULL; pc->ops->applysymmetricleft = NULL; pc->ops->applysymmetricright = NULL; PetscFunctionReturn(0); }