xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision 85e3dda7fddc385e91bc763ba6bbff3d650f88ca)
1 #define PETSCKSP_DLL
2 
3 #include "private/pcimpl.h"   /*I "petscpc.h" I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCApply_Mat"
7 static PetscErrorCode PCApply_Mat(PC pc,Vec x,Vec y)
8 {
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   ierr = MatMult(pc->pmat,x,y);CHKERRQ(ierr);
13   PetscFunctionReturn(0);
14 }
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PCApplyTranspose_Mat"
18 static PetscErrorCode PCApplyTranspose_Mat(PC pc,Vec x,Vec y)
19 {
20   PetscErrorCode ierr;
21 
22   PetscFunctionBegin;
23   ierr = MatMultTranspose(pc->pmat,x,y);CHKERRQ(ierr);
24   PetscFunctionReturn(0);
25 }
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "PCDestroy_Mat"
29 static PetscErrorCode PCDestroy_Mat(PC pc)
30 {
31   PetscFunctionBegin;
32   PetscFunctionReturn(0);
33 }
34 
35 /*MC
36      PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied
37              in PCSetOperators() or KSPSetOperators()
38 
39    Notes:  This one is a little strange. One rarely has an explict matrix that approximates the
40          inverse of the matrix they wish to solve for.
41 
42    Level: intermediate
43 
44 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
45            PCSHELL
46 
47 M*/
48 
49 EXTERN_C_BEGIN
50 #undef __FUNCT__
51 #define __FUNCT__ "PCCreate_Mat"
52 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Mat(PC pc)
53 {
54   PetscFunctionBegin;
55   pc->ops->apply               = PCApply_Mat;
56   pc->ops->applytranspose      = PCApplyTranspose_Mat;
57   pc->ops->setup               = 0;
58   pc->ops->destroy             = PCDestroy_Mat;
59   pc->ops->setfromoptions      = 0;
60   pc->ops->view                = 0;
61   pc->ops->applyrichardson     = 0;
62   pc->ops->applysymmetricleft  = 0;
63   pc->ops->applysymmetricright = 0;
64   PetscFunctionReturn(0);
65 }
66 EXTERN_C_END
67 
68