xref: /petsc/src/ksp/pc/impls/none/none.c (revision 4b9ad92859ccb93b5e851e53cb8c4c04ac10e155)
1*4b9ad928SBarry Smith /*$Id: none.c,v 1.30 2001/03/23 23:23:06 balay Exp $*/
2*4b9ad928SBarry Smith /*
3*4b9ad928SBarry Smith     Identity preconditioner, simply copies vector x to y.
4*4b9ad928SBarry Smith */
5*4b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h"          /*I "petscpc.h" I*/
6*4b9ad928SBarry Smith 
7*4b9ad928SBarry Smith #undef __FUNCT__
8*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_None"
9*4b9ad928SBarry Smith int PCApply_None(PC pc,Vec x,Vec y)
10*4b9ad928SBarry Smith {
11*4b9ad928SBarry Smith   int ierr;
12*4b9ad928SBarry Smith 
13*4b9ad928SBarry Smith   PetscFunctionBegin;
14*4b9ad928SBarry Smith   ierr = VecCopy(x,y);CHKERRQ(ierr);
15*4b9ad928SBarry Smith   PetscFunctionReturn(0);
16*4b9ad928SBarry Smith }
17*4b9ad928SBarry Smith 
18*4b9ad928SBarry Smith /*MC
19*4b9ad928SBarry Smith      PCNONE - This is used when you wish to employ a nonpreconditioned
20*4b9ad928SBarry Smith              Krylov method.
21*4b9ad928SBarry Smith 
22*4b9ad928SBarry Smith    Level: beginner
23*4b9ad928SBarry Smith 
24*4b9ad928SBarry Smith   Concepts: preconditioners
25*4b9ad928SBarry Smith 
26*4b9ad928SBarry Smith   Notes: This is implemented by a VecCopy()
27*4b9ad928SBarry Smith 
28*4b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
29*4b9ad928SBarry Smith M*/
30*4b9ad928SBarry Smith 
31*4b9ad928SBarry Smith EXTERN_C_BEGIN
32*4b9ad928SBarry Smith #undef __FUNCT__
33*4b9ad928SBarry Smith #define __FUNCT__ "PCCreate_None"
34*4b9ad928SBarry Smith int PCCreate_None(PC pc)
35*4b9ad928SBarry Smith {
36*4b9ad928SBarry Smith   PetscFunctionBegin;
37*4b9ad928SBarry Smith   pc->ops->apply               = PCApply_None;
38*4b9ad928SBarry Smith   pc->ops->applytranspose      = PCApply_None;
39*4b9ad928SBarry Smith   pc->ops->destroy             = 0;
40*4b9ad928SBarry Smith   pc->ops->setup               = 0;
41*4b9ad928SBarry Smith   pc->ops->view                = 0;
42*4b9ad928SBarry Smith   pc->ops->applysymmetricleft  = PCApply_None;
43*4b9ad928SBarry Smith   pc->ops->applysymmetricright = PCApply_None;
44*4b9ad928SBarry Smith 
45*4b9ad928SBarry Smith   pc->data                     = 0;
46*4b9ad928SBarry Smith   PetscFunctionReturn(0);
47*4b9ad928SBarry Smith }
48*4b9ad928SBarry Smith EXTERN_C_END
49