xref: /petsc/src/ksp/pc/impls/none/none.c (revision 347806e75ae86c9025d175dd7fbb3f1396681ded)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith     Identity preconditioner, simply copies vector x to y.
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
54b9ad928SBarry Smith 
PCApply_None(PC pc,Vec x,Vec y)666976f2fSJacob Faibussowitsch static PetscErrorCode PCApply_None(PC pc, Vec x, Vec y)
7d71ae5a4SJacob Faibussowitsch {
84b9ad928SBarry Smith   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, y));
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114b9ad928SBarry Smith }
124b9ad928SBarry Smith 
PCMatApply_None(PC pc,Mat X,Mat Y)1366976f2fSJacob Faibussowitsch static PetscErrorCode PCMatApply_None(PC pc, Mat X, Mat Y)
14d71ae5a4SJacob Faibussowitsch {
157b6e2003SPierre Jolivet   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(MatCopy(X, Y, SAME_NONZERO_PATTERN));
173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187b6e2003SPierre Jolivet }
197b6e2003SPierre Jolivet 
204b9ad928SBarry Smith /*MC
214b9ad928SBarry Smith      PCNONE - This is used when you wish to employ a nonpreconditioned
224b9ad928SBarry Smith              Krylov method.
234b9ad928SBarry Smith 
244b9ad928SBarry Smith    Level: beginner
254b9ad928SBarry Smith 
26f1580f4eSBarry Smith   Developer Note:
27f1580f4eSBarry Smith   This is implemented by a `VecCopy()`. It would be nice if the `KSP` implementations could be organized to avoid this copy without making them
28f1580f4eSBarry Smith   more complex.
294b9ad928SBarry Smith 
30562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`
314b9ad928SBarry Smith M*/
324b9ad928SBarry Smith 
PCCreate_None(PC pc)33d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_None(PC pc)
34d71ae5a4SJacob Faibussowitsch {
354b9ad928SBarry Smith   PetscFunctionBegin;
364b9ad928SBarry Smith   pc->ops->apply               = PCApply_None;
377b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_None;
384b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_None;
39*4dbf25a8SPierre Jolivet   pc->ops->matapplytranspose   = PCMatApply_None;
400a545947SLisandro Dalcin   pc->ops->destroy             = NULL;
410a545947SLisandro Dalcin   pc->ops->setup               = NULL;
420a545947SLisandro Dalcin   pc->ops->view                = NULL;
434b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApply_None;
444b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApply_None;
454b9ad928SBarry Smith 
460a545947SLisandro Dalcin   pc->data = NULL;
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
484b9ad928SBarry Smith }
49