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