Actual source code: pcmat.c
2: #include <petsc/private/pcimpl.h>
4: static PetscErrorCode PCApply_Mat(PC pc,Vec x,Vec y)
5: {
9: MatMult(pc->pmat,x,y);
10: return(0);
11: }
13: static PetscErrorCode PCMatApply_Mat(PC pc,Mat X,Mat Y)
14: {
18: MatMatMult(pc->pmat,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);
19: return(0);
20: }
22: static PetscErrorCode PCApplyTranspose_Mat(PC pc,Vec x,Vec y)
23: {
27: MatMultTranspose(pc->pmat,x,y);
28: return(0);
29: }
31: static PetscErrorCode PCDestroy_Mat(PC pc)
32: {
34: return(0);
35: }
37: /*MC
38: PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied
39: in PCSetOperators() or KSPSetOperators()
41: Notes:
42: This one is a little strange. One rarely has an explict matrix that approximates the
43: inverse of the matrix they wish to solve for.
45: Level: intermediate
47: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
48: PCSHELL
50: M*/
52: PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
53: {
55: pc->ops->apply = PCApply_Mat;
56: pc->ops->matapply = PCMatApply_Mat;
57: pc->ops->applytranspose = PCApplyTranspose_Mat;
58: pc->ops->setup = NULL;
59: pc->ops->destroy = PCDestroy_Mat;
60: pc->ops->setfromoptions = NULL;
61: pc->ops->view = NULL;
62: pc->ops->applyrichardson = NULL;
63: pc->ops->applysymmetricleft = NULL;
64: pc->ops->applysymmetricright = NULL;
65: return(0);
66: }