xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
14b9ad928SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
34b9ad928SBarry Smith 
4*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Mat(PC pc, Vec x, Vec y)
5*d71ae5a4SJacob Faibussowitsch {
64b9ad928SBarry Smith   PetscFunctionBegin;
79566063dSJacob Faibussowitsch   PetscCall(MatMult(pc->pmat, x, y));
84b9ad928SBarry Smith   PetscFunctionReturn(0);
94b9ad928SBarry Smith }
104b9ad928SBarry Smith 
11*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Mat(PC pc, Mat X, Mat Y)
12*d71ae5a4SJacob Faibussowitsch {
137b6e2003SPierre Jolivet   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(MatMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
157b6e2003SPierre Jolivet   PetscFunctionReturn(0);
167b6e2003SPierre Jolivet }
177b6e2003SPierre Jolivet 
18*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Mat(PC pc, Vec x, Vec y)
19*d71ae5a4SJacob Faibussowitsch {
204b9ad928SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(pc->pmat, x, y));
224b9ad928SBarry Smith   PetscFunctionReturn(0);
234b9ad928SBarry Smith }
244b9ad928SBarry Smith 
25*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Mat(PC pc)
26*d71ae5a4SJacob Faibussowitsch {
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
33f1580f4eSBarry Smith              in `PCSetOperators()` or `KSPSetOperators()`
344b9ad928SBarry Smith 
35f1580f4eSBarry Smith    Note:
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 
41db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
42db781477SPatrick Sanan           `PCSHELL`
434b9ad928SBarry Smith M*/
444b9ad928SBarry Smith 
45*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
46*d71ae5a4SJacob Faibussowitsch {
474b9ad928SBarry Smith   PetscFunctionBegin;
484b9ad928SBarry Smith   pc->ops->apply               = PCApply_Mat;
497b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Mat;
504b9ad928SBarry Smith   pc->ops->applytranspose      = PCApplyTranspose_Mat;
510a545947SLisandro Dalcin   pc->ops->setup               = NULL;
524b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Mat;
530a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
540a545947SLisandro Dalcin   pc->ops->view                = NULL;
550a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
560a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
570a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
60