xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision a5b23f4acc7afc99d3844ebd5fb65a81c16e8b8c)
14b9ad928SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
34b9ad928SBarry Smith 
46849ba73SBarry Smith static PetscErrorCode PCApply_Mat(PC pc,Vec x,Vec y)
54b9ad928SBarry Smith {
6dfbe8321SBarry Smith   PetscErrorCode ierr;
74b9ad928SBarry Smith 
84b9ad928SBarry Smith   PetscFunctionBegin;
94b9ad928SBarry Smith   ierr = MatMult(pc->pmat,x,y);CHKERRQ(ierr);
104b9ad928SBarry Smith   PetscFunctionReturn(0);
114b9ad928SBarry Smith }
124b9ad928SBarry Smith 
137b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Mat(PC pc,Mat X,Mat Y)
147b6e2003SPierre Jolivet {
157b6e2003SPierre Jolivet   PetscErrorCode ierr;
167b6e2003SPierre Jolivet 
177b6e2003SPierre Jolivet   PetscFunctionBegin;
187b6e2003SPierre Jolivet   ierr = MatMatMult(pc->pmat,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr);
197b6e2003SPierre Jolivet   PetscFunctionReturn(0);
207b6e2003SPierre Jolivet }
217b6e2003SPierre Jolivet 
226849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_Mat(PC pc,Vec x,Vec y)
234b9ad928SBarry Smith {
24dfbe8321SBarry Smith   PetscErrorCode ierr;
254b9ad928SBarry Smith 
264b9ad928SBarry Smith   PetscFunctionBegin;
274b9ad928SBarry Smith   ierr = MatMultTranspose(pc->pmat,x,y);CHKERRQ(ierr);
284b9ad928SBarry Smith   PetscFunctionReturn(0);
294b9ad928SBarry Smith }
304b9ad928SBarry Smith 
316849ba73SBarry Smith static PetscErrorCode PCDestroy_Mat(PC pc)
324b9ad928SBarry Smith {
334b9ad928SBarry Smith   PetscFunctionBegin;
344b9ad928SBarry Smith   PetscFunctionReturn(0);
354b9ad928SBarry Smith }
364b9ad928SBarry Smith 
374b9ad928SBarry Smith /*MC
384b9ad928SBarry Smith      PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied
394b9ad928SBarry Smith              in PCSetOperators() or KSPSetOperators()
404b9ad928SBarry Smith 
4195452b02SPatrick Sanan    Notes:
42*a5b23f4aSJose E. Roman     This one is a little strange. One rarely has an explicit matrix that approximates the
434b9ad928SBarry Smith          inverse of the matrix they wish to solve for.
444b9ad928SBarry Smith 
454b9ad928SBarry Smith    Level: intermediate
464b9ad928SBarry Smith 
474b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
484b9ad928SBarry Smith            PCSHELL
494b9ad928SBarry Smith 
504b9ad928SBarry Smith M*/
514b9ad928SBarry Smith 
528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
534b9ad928SBarry Smith {
544b9ad928SBarry Smith   PetscFunctionBegin;
554b9ad928SBarry Smith   pc->ops->apply               = PCApply_Mat;
567b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Mat;
574b9ad928SBarry Smith   pc->ops->applytranspose      = PCApplyTranspose_Mat;
580a545947SLisandro Dalcin   pc->ops->setup               = NULL;
594b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Mat;
600a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
610a545947SLisandro Dalcin   pc->ops->view                = NULL;
620a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
630a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
640a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
654b9ad928SBarry Smith   PetscFunctionReturn(0);
664b9ad928SBarry Smith }
674b9ad928SBarry Smith 
68