14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 34b9ad928SBarry Smith 4*9371c9d4SSatish Balay static PetscErrorCode PCApply_Mat(PC pc, Vec x, Vec y) { 54b9ad928SBarry Smith PetscFunctionBegin; 69566063dSJacob Faibussowitsch PetscCall(MatMult(pc->pmat, x, y)); 74b9ad928SBarry Smith PetscFunctionReturn(0); 84b9ad928SBarry Smith } 94b9ad928SBarry Smith 10*9371c9d4SSatish Balay static PetscErrorCode PCMatApply_Mat(PC pc, Mat X, Mat Y) { 117b6e2003SPierre Jolivet PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(MatMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 137b6e2003SPierre Jolivet PetscFunctionReturn(0); 147b6e2003SPierre Jolivet } 157b6e2003SPierre Jolivet 16*9371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_Mat(PC pc, Vec x, Vec y) { 174b9ad928SBarry Smith PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->pmat, x, y)); 194b9ad928SBarry Smith PetscFunctionReturn(0); 204b9ad928SBarry Smith } 214b9ad928SBarry Smith 22*9371c9d4SSatish Balay static PetscErrorCode PCDestroy_Mat(PC pc) { 234b9ad928SBarry Smith PetscFunctionBegin; 244b9ad928SBarry Smith PetscFunctionReturn(0); 254b9ad928SBarry Smith } 264b9ad928SBarry Smith 274b9ad928SBarry Smith /*MC 284b9ad928SBarry Smith PCMAT - A preconditioner obtained by multiplying by the preconditioner matrix supplied 294b9ad928SBarry Smith in PCSetOperators() or KSPSetOperators() 304b9ad928SBarry Smith 3195452b02SPatrick Sanan Notes: 32a5b23f4aSJose E. Roman This one is a little strange. One rarely has an explicit matrix that approximates the 334b9ad928SBarry Smith inverse of the matrix they wish to solve for. 344b9ad928SBarry Smith 354b9ad928SBarry Smith Level: intermediate 364b9ad928SBarry Smith 37db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 38db781477SPatrick Sanan `PCSHELL` 394b9ad928SBarry Smith 404b9ad928SBarry Smith M*/ 414b9ad928SBarry Smith 42*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc) { 434b9ad928SBarry Smith PetscFunctionBegin; 444b9ad928SBarry Smith pc->ops->apply = PCApply_Mat; 457b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Mat; 464b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Mat; 470a545947SLisandro Dalcin pc->ops->setup = NULL; 484b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Mat; 490a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 500a545947SLisandro Dalcin pc->ops->view = NULL; 510a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 520a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 530a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 544b9ad928SBarry Smith PetscFunctionReturn(0); 554b9ad928SBarry Smith } 56