xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision db7814771ca77b190574494e87b584e981451db0)
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 {
64b9ad928SBarry Smith   PetscFunctionBegin;
79566063dSJacob Faibussowitsch   PetscCall(MatMult(pc->pmat,x,y));
84b9ad928SBarry Smith   PetscFunctionReturn(0);
94b9ad928SBarry Smith }
104b9ad928SBarry Smith 
117b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_Mat(PC pc,Mat X,Mat Y)
127b6e2003SPierre Jolivet {
137b6e2003SPierre Jolivet   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(MatMatMult(pc->pmat,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
157b6e2003SPierre Jolivet   PetscFunctionReturn(0);
167b6e2003SPierre Jolivet }
177b6e2003SPierre Jolivet 
186849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_Mat(PC pc,Vec x,Vec y)
194b9ad928SBarry Smith {
204b9ad928SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(pc->pmat,x,y));
224b9ad928SBarry Smith   PetscFunctionReturn(0);
234b9ad928SBarry Smith }
244b9ad928SBarry Smith 
256849ba73SBarry Smith static PetscErrorCode PCDestroy_Mat(PC pc)
264b9ad928SBarry Smith {
274b9ad928SBarry Smith   PetscFunctionBegin;
284b9ad928SBarry Smith   PetscFunctionReturn(0);
294b9ad928SBarry Smith }
304b9ad928SBarry Smith 
314b9ad928SBarry Smith /*MC
324b9ad928SBarry Smith      PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied
334b9ad928SBarry Smith              in PCSetOperators() or KSPSetOperators()
344b9ad928SBarry Smith 
3595452b02SPatrick Sanan    Notes:
36a5b23f4aSJose E. Roman     This one is a little strange. One rarely has an explicit matrix that approximates the
374b9ad928SBarry Smith          inverse of the matrix they wish to solve for.
384b9ad928SBarry Smith 
394b9ad928SBarry Smith    Level: intermediate
404b9ad928SBarry Smith 
41*db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
42*db781477SPatrick Sanan           `PCSHELL`
434b9ad928SBarry Smith 
444b9ad928SBarry Smith M*/
454b9ad928SBarry Smith 
468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
474b9ad928SBarry Smith {
484b9ad928SBarry Smith   PetscFunctionBegin;
494b9ad928SBarry Smith   pc->ops->apply               = PCApply_Mat;
507b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Mat;
514b9ad928SBarry Smith   pc->ops->applytranspose      = PCApplyTranspose_Mat;
520a545947SLisandro Dalcin   pc->ops->setup               = NULL;
534b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Mat;
540a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
550a545947SLisandro Dalcin   pc->ops->view                = NULL;
560a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
570a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
580a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
594b9ad928SBarry Smith   PetscFunctionReturn(0);
604b9ad928SBarry Smith }
61