xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
14b9ad928SBarry Smith 
2*af0996ceSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
34b9ad928SBarry Smith 
44b9ad928SBarry Smith #undef __FUNCT__
54b9ad928SBarry Smith #define __FUNCT__ "PCApply_Mat"
66849ba73SBarry Smith static PetscErrorCode PCApply_Mat(PC pc,Vec x,Vec y)
74b9ad928SBarry Smith {
8dfbe8321SBarry Smith   PetscErrorCode ierr;
94b9ad928SBarry Smith 
104b9ad928SBarry Smith   PetscFunctionBegin;
114b9ad928SBarry Smith   ierr = MatMult(pc->pmat,x,y);CHKERRQ(ierr);
124b9ad928SBarry Smith   PetscFunctionReturn(0);
134b9ad928SBarry Smith }
144b9ad928SBarry Smith 
154b9ad928SBarry Smith #undef __FUNCT__
164b9ad928SBarry Smith #define __FUNCT__ "PCApplyTranspose_Mat"
176849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_Mat(PC pc,Vec x,Vec y)
184b9ad928SBarry Smith {
19dfbe8321SBarry Smith   PetscErrorCode ierr;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
224b9ad928SBarry Smith   ierr = MatMultTranspose(pc->pmat,x,y);CHKERRQ(ierr);
234b9ad928SBarry Smith   PetscFunctionReturn(0);
244b9ad928SBarry Smith }
254b9ad928SBarry Smith 
264b9ad928SBarry Smith #undef __FUNCT__
274b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Mat"
286849ba73SBarry Smith static PetscErrorCode PCDestroy_Mat(PC pc)
294b9ad928SBarry Smith {
304b9ad928SBarry Smith   PetscFunctionBegin;
314b9ad928SBarry Smith   PetscFunctionReturn(0);
324b9ad928SBarry Smith }
334b9ad928SBarry Smith 
344b9ad928SBarry Smith /*MC
354b9ad928SBarry Smith      PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied
364b9ad928SBarry Smith              in PCSetOperators() or KSPSetOperators()
374b9ad928SBarry Smith 
384b9ad928SBarry Smith    Notes:  This one is a little strange. One rarely has an explict matrix that approximates the
394b9ad928SBarry Smith          inverse of the matrix they wish to solve for.
404b9ad928SBarry Smith 
414b9ad928SBarry Smith    Level: intermediate
424b9ad928SBarry Smith 
434b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
444b9ad928SBarry Smith            PCSHELL
454b9ad928SBarry Smith 
464b9ad928SBarry Smith M*/
474b9ad928SBarry Smith 
484b9ad928SBarry Smith #undef __FUNCT__
494b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Mat"
508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
514b9ad928SBarry Smith {
524b9ad928SBarry Smith   PetscFunctionBegin;
534b9ad928SBarry Smith   pc->ops->apply               = PCApply_Mat;
544b9ad928SBarry Smith   pc->ops->applytranspose      = PCApplyTranspose_Mat;
554b9ad928SBarry Smith   pc->ops->setup               = 0;
564b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Mat;
574b9ad928SBarry Smith   pc->ops->setfromoptions      = 0;
584b9ad928SBarry Smith   pc->ops->view                = 0;
594b9ad928SBarry Smith   pc->ops->applyrichardson     = 0;
604b9ad928SBarry Smith   pc->ops->applysymmetricleft  = 0;
614b9ad928SBarry Smith   pc->ops->applysymmetricright = 0;
624b9ad928SBarry Smith   PetscFunctionReturn(0);
634b9ad928SBarry Smith }
644b9ad928SBarry Smith 
65