xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision e662bc509b40f26ba5a3dffd4204028ba0846ecb)
137eeb815SJakub Kruzik 
237eeb815SJakub Kruzik /*  --------------------------------------------------------------------
337eeb815SJakub Kruzik 
437eeb815SJakub Kruzik      This file implements a Deflation preconditioner in PETSc as part of PC.
537eeb815SJakub Kruzik      You can use this as a starting point for implementing your own
637eeb815SJakub Kruzik      preconditioner that is not provided with PETSc. (You might also consider
737eeb815SJakub Kruzik      just using PCSHELL)
837eeb815SJakub Kruzik 
937eeb815SJakub Kruzik      The following basic routines are required for each preconditioner.
1037eeb815SJakub Kruzik           PCCreate_XXX()          - Creates a preconditioner context
1137eeb815SJakub Kruzik           PCSetFromOptions_XXX()  - Sets runtime options
1237eeb815SJakub Kruzik           PCApply_XXX()           - Applies the preconditioner
1337eeb815SJakub Kruzik           PCDestroy_XXX()         - Destroys the preconditioner context
1437eeb815SJakub Kruzik      where the suffix "_XXX" denotes a particular implementation, in
1537eeb815SJakub Kruzik      this case we use _Deflation (e.g., PCCreate_Deflation, PCApply_Deflation).
1637eeb815SJakub Kruzik      These routines are actually called via the common user interface
1737eeb815SJakub Kruzik      routines PCCreate(), PCSetFromOptions(), PCApply(), and PCDestroy(),
1837eeb815SJakub Kruzik      so the application code interface remains identical for all
1937eeb815SJakub Kruzik      preconditioners.
2037eeb815SJakub Kruzik 
2137eeb815SJakub Kruzik      Another key routine is:
2237eeb815SJakub Kruzik           PCSetUp_XXX()           - Prepares for the use of a preconditioner
2337eeb815SJakub Kruzik      by setting data structures and options.   The interface routine PCSetUp()
2437eeb815SJakub Kruzik      is not usually called directly by the user, but instead is called by
2537eeb815SJakub Kruzik      PCApply() if necessary.
2637eeb815SJakub Kruzik 
2737eeb815SJakub Kruzik      Additional basic routines are:
2837eeb815SJakub Kruzik           PCView_XXX()            - Prints details of runtime options that
2937eeb815SJakub Kruzik                                     have actually been used.
3037eeb815SJakub Kruzik      These are called by application codes via the interface routines
3137eeb815SJakub Kruzik      PCView().
3237eeb815SJakub Kruzik 
3337eeb815SJakub Kruzik      The various types of solvers (preconditioners, Krylov subspace methods,
3437eeb815SJakub Kruzik      nonlinear solvers, timesteppers) are all organized similarly, so the
3537eeb815SJakub Kruzik      above description applies to these categories also.  One exception is
3637eeb815SJakub Kruzik      that the analogues of PCApply() for these components are KSPSolve(),
3737eeb815SJakub Kruzik      SNESSolve(), and TSSolve().
3837eeb815SJakub Kruzik 
3937eeb815SJakub Kruzik      Additional optional functionality unique to preconditioners is left and
4037eeb815SJakub Kruzik      right symmetric preconditioner application via PCApplySymmetricLeft()
4137eeb815SJakub Kruzik      and PCApplySymmetricRight().  The Deflation implementation is
4237eeb815SJakub Kruzik      PCApplySymmetricLeftOrRight_Deflation().
4337eeb815SJakub Kruzik 
4437eeb815SJakub Kruzik     -------------------------------------------------------------------- */
4537eeb815SJakub Kruzik 
4637eeb815SJakub Kruzik /*
4737eeb815SJakub Kruzik    Include files needed for the Deflation preconditioner:
4837eeb815SJakub Kruzik      pcimpl.h - private include file intended for use by all preconditioners
4937eeb815SJakub Kruzik */
5037eeb815SJakub Kruzik 
5137eeb815SJakub Kruzik #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
5237eeb815SJakub Kruzik 
5337eeb815SJakub Kruzik const char *const PCDeflationTypes[]    = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0};
5437eeb815SJakub Kruzik 
5537eeb815SJakub Kruzik /*
5637eeb815SJakub Kruzik    Private context (data structure) for the deflation preconditioner.
5737eeb815SJakub Kruzik */
5837eeb815SJakub Kruzik typedef struct {
5937eeb815SJakub Kruzik   PetscBool init;            /* do only init step - error correction of direction is omitted */
6037eeb815SJakub Kruzik   PetscBool pre;             /* start with x0 being the solution in the deflation space */
6137eeb815SJakub Kruzik   PetscBool correct;         /* add CP (Qr) correction to descent direction */
6237eeb815SJakub Kruzik   PetscBool truenorm;
6337eeb815SJakub Kruzik   PetscBool adaptiveconv;
6437eeb815SJakub Kruzik   PetscReal adaptiveconst;
6537eeb815SJakub Kruzik   PetscInt  reductionfact;
6637eeb815SJakub Kruzik   Mat       W,Wt,AW,WtAW;    /* deflation space, coarse problem mats */
6737eeb815SJakub Kruzik   KSP       WtAWinv;         /* deflation coarse problem */
6837eeb815SJakub Kruzik   Vec       *work;
6937eeb815SJakub Kruzik 
7037eeb815SJakub Kruzik   PCDEFLATIONSpaceType spacetype;
7137eeb815SJakub Kruzik   PetscInt             spacesize;
7237eeb815SJakub Kruzik   PetscInt             nestedlvl;
7337eeb815SJakub Kruzik   PetscInt             maxnestedlvl;
7437eeb815SJakub Kruzik   PetscBool            extendsp;
7537eeb815SJakub Kruzik } PC_Deflation;
7637eeb815SJakub Kruzik 
7737eeb815SJakub Kruzik static PetscErrorCode  PCDeflationSetType_Deflation(PC pc,PCDeflationType type)
7837eeb815SJakub Kruzik {
7937eeb815SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
8037eeb815SJakub Kruzik 
8137eeb815SJakub Kruzik   PetscFunctionBegin;
8237eeb815SJakub Kruzik   def->init = PETSC_FALSE;
8337eeb815SJakub Kruzik   def->pre = PETSC_FALSE;
8437eeb815SJakub Kruzik   if (type == PC_DEFLATION_INIT) {
8537eeb815SJakub Kruzik     def->init = PETSC_TRUE;
8637eeb815SJakub Kruzik     def->pre  = PETSC_TRUE;
8737eeb815SJakub Kruzik   } else if (type == PC_DEFLATION_PRE) {
8837eeb815SJakub Kruzik     def->pre  = PETSC_TRUE;
8937eeb815SJakub Kruzik   }
9037eeb815SJakub Kruzik   PetscFunctionReturn(0);
9137eeb815SJakub Kruzik }
9237eeb815SJakub Kruzik 
9337eeb815SJakub Kruzik static PetscErrorCode  PCDeflationGetType_Deflation(PC pc,PCDeflationType *type)
9437eeb815SJakub Kruzik {
9537eeb815SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
9637eeb815SJakub Kruzik 
9737eeb815SJakub Kruzik   PetscFunctionBegin;
9837eeb815SJakub Kruzik   if (def->init) {
9937eeb815SJakub Kruzik     *type = PC_DEFLATION_INIT;
10037eeb815SJakub Kruzik   } else if (def->pre) {
10137eeb815SJakub Kruzik     *type = PC_DEFLATION_PRE;
10237eeb815SJakub Kruzik   } else {
10337eeb815SJakub Kruzik     *type = PC_DEFLATION_POST;
10437eeb815SJakub Kruzik   }
10537eeb815SJakub Kruzik   PetscFunctionReturn(0);
10637eeb815SJakub Kruzik }
10737eeb815SJakub Kruzik 
108*e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
109*e662bc50SJakub Kruzik {
110*e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
111*e662bc50SJakub Kruzik   PetscErrorCode ierr;
112*e662bc50SJakub Kruzik 
113*e662bc50SJakub Kruzik   PetscFunctionBegin;
114*e662bc50SJakub Kruzik   if (transpose) {
115*e662bc50SJakub Kruzik     def->Wt = W;
116*e662bc50SJakub Kruzik     def->W = NULL;
117*e662bc50SJakub Kruzik   } else {
118*e662bc50SJakub Kruzik     def->W = W;
119*e662bc50SJakub Kruzik   }
120*e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
121*e662bc50SJakub Kruzik   PetscFunctionReturn(0);
122*e662bc50SJakub Kruzik }
123*e662bc50SJakub Kruzik 
124*e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */
125*e662bc50SJakub Kruzik /*@
126*e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
127*e662bc50SJakub Kruzik 
128*e662bc50SJakub Kruzik    Logically Collective on PC
129*e662bc50SJakub Kruzik 
130*e662bc50SJakub Kruzik    Input Parameters:
131*e662bc50SJakub Kruzik +  pc - the preconditioner context
132*e662bc50SJakub Kruzik .  W  - deflation matrix
133*e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
134*e662bc50SJakub Kruzik 
135*e662bc50SJakub Kruzik    Level: intermediate
136*e662bc50SJakub Kruzik 
137*e662bc50SJakub Kruzik .seealso: PCDEFLATION
138*e662bc50SJakub Kruzik @*/
139*e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
140*e662bc50SJakub Kruzik {
141*e662bc50SJakub Kruzik   PetscErrorCode ierr;
142*e662bc50SJakub Kruzik 
143*e662bc50SJakub Kruzik   PetscFunctionBegin;
144*e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
145*e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
146*e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
147*e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
148*e662bc50SJakub Kruzik   PetscFunctionReturn(0);
149*e662bc50SJakub Kruzik }
150*e662bc50SJakub Kruzik 
151*e662bc50SJakub Kruzik 
15237eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
15337eeb815SJakub Kruzik /*
15437eeb815SJakub Kruzik    PCSetUp_Deflation - Prepares for the use of the Deflation preconditioner
15537eeb815SJakub Kruzik                     by setting data structures and options.
15637eeb815SJakub Kruzik 
15737eeb815SJakub Kruzik    Input Parameter:
15837eeb815SJakub Kruzik .  pc - the preconditioner context
15937eeb815SJakub Kruzik 
16037eeb815SJakub Kruzik    Application Interface Routine: PCSetUp()
16137eeb815SJakub Kruzik 
16237eeb815SJakub Kruzik    Notes:
16337eeb815SJakub Kruzik    The interface routine PCSetUp() is not usually called directly by
16437eeb815SJakub Kruzik    the user, but instead is called by PCApply() if necessary.
16537eeb815SJakub Kruzik */
16637eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
16737eeb815SJakub Kruzik {
16837eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
16937eeb815SJakub Kruzik   Vec            diag,diagsqrt;
17037eeb815SJakub Kruzik   PetscErrorCode ierr;
17137eeb815SJakub Kruzik   PetscInt       n,i;
17237eeb815SJakub Kruzik   PetscScalar    *x;
17337eeb815SJakub Kruzik   PetscBool      zeroflag = PETSC_FALSE;
17437eeb815SJakub Kruzik 
17537eeb815SJakub Kruzik   PetscFunctionBegin;
17637eeb815SJakub Kruzik   /*
17737eeb815SJakub Kruzik        For most preconditioners the code would begin here something like
17837eeb815SJakub Kruzik 
17937eeb815SJakub Kruzik   if (pc->setupcalled == 0) { allocate space the first time this is ever called
18037eeb815SJakub Kruzik     ierr = MatCreateVecs(pc->mat,&jac->diag);CHKERRQ(ierr);
18137eeb815SJakub Kruzik     PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->diag);
18237eeb815SJakub Kruzik   }
18337eeb815SJakub Kruzik 
18437eeb815SJakub Kruzik     But for this preconditioner we want to support use of both the matrix' diagonal
18537eeb815SJakub Kruzik     elements (for left or right preconditioning) and square root of diagonal elements
18637eeb815SJakub Kruzik     (for symmetric preconditioning).  Hence we do not allocate space here, since we
18737eeb815SJakub Kruzik     don't know at this point which will be needed (diag and/or diagsqrt) until the user
18837eeb815SJakub Kruzik     applies the preconditioner, and we don't want to allocate BOTH unless we need
18937eeb815SJakub Kruzik     them both.  Thus, the diag and diagsqrt are allocated in PCSetUp_Deflation_NonSymmetric()
19037eeb815SJakub Kruzik     and PCSetUp_Deflation_Symmetric(), respectively.
19137eeb815SJakub Kruzik   */
19237eeb815SJakub Kruzik 
19337eeb815SJakub Kruzik   /*
19437eeb815SJakub Kruzik     Here we set up the preconditioner; that is, we copy the diagonal values from
19537eeb815SJakub Kruzik     the matrix and put them into a format to make them quick to apply as a preconditioner.
19637eeb815SJakub Kruzik   */
19737eeb815SJakub Kruzik   if (zeroflag) {
19837eeb815SJakub Kruzik     ierr = PetscInfo(pc,"Zero detected in diagonal of matrix, using 1 at those locations\n");CHKERRQ(ierr);
19937eeb815SJakub Kruzik   }
20037eeb815SJakub Kruzik   PetscFunctionReturn(0);
20137eeb815SJakub Kruzik }
20237eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
20337eeb815SJakub Kruzik /*
20437eeb815SJakub Kruzik    PCApply_Deflation - Applies the Deflation preconditioner to a vector.
20537eeb815SJakub Kruzik 
20637eeb815SJakub Kruzik    Input Parameters:
20737eeb815SJakub Kruzik .  pc - the preconditioner context
20837eeb815SJakub Kruzik .  x - input vector
20937eeb815SJakub Kruzik 
21037eeb815SJakub Kruzik    Output Parameter:
21137eeb815SJakub Kruzik .  y - output vector
21237eeb815SJakub Kruzik 
21337eeb815SJakub Kruzik    Application Interface Routine: PCApply()
21437eeb815SJakub Kruzik  */
21537eeb815SJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec x,Vec y)
21637eeb815SJakub Kruzik {
21737eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
21837eeb815SJakub Kruzik   PetscErrorCode ierr;
21937eeb815SJakub Kruzik 
22037eeb815SJakub Kruzik   PetscFunctionBegin;
22137eeb815SJakub Kruzik   PetscFunctionReturn(0);
22237eeb815SJakub Kruzik }
22337eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
22437eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
22537eeb815SJakub Kruzik {
22637eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
22737eeb815SJakub Kruzik   PetscErrorCode ierr;
22837eeb815SJakub Kruzik 
22937eeb815SJakub Kruzik   PetscFunctionBegin;
23037eeb815SJakub Kruzik   PetscFunctionReturn(0);
23137eeb815SJakub Kruzik }
23237eeb815SJakub Kruzik 
23337eeb815SJakub Kruzik /*
23437eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
23537eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
23637eeb815SJakub Kruzik 
23737eeb815SJakub Kruzik    Input Parameter:
23837eeb815SJakub Kruzik .  pc - the preconditioner context
23937eeb815SJakub Kruzik 
24037eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
24137eeb815SJakub Kruzik */
24237eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
24337eeb815SJakub Kruzik {
24437eeb815SJakub Kruzik   PetscErrorCode ierr;
24537eeb815SJakub Kruzik 
24637eeb815SJakub Kruzik   PetscFunctionBegin;
24737eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
24837eeb815SJakub Kruzik 
24937eeb815SJakub Kruzik   /*
25037eeb815SJakub Kruzik       Free the private data structure that was hanging off the PC
25137eeb815SJakub Kruzik   */
25237eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
25337eeb815SJakub Kruzik   PetscFunctionReturn(0);
25437eeb815SJakub Kruzik }
25537eeb815SJakub Kruzik 
25637eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
25737eeb815SJakub Kruzik {
25837eeb815SJakub Kruzik   PC_Deflation      *jac = (PC_Deflation*)pc->data;
25937eeb815SJakub Kruzik   PetscErrorCode ierr;
26037eeb815SJakub Kruzik   PetscBool      flg;
26137eeb815SJakub Kruzik   PCDeflationType   deflt,type;
26237eeb815SJakub Kruzik 
26337eeb815SJakub Kruzik   PetscFunctionBegin;
26437eeb815SJakub Kruzik   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
26537eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
26637eeb815SJakub Kruzik   ierr = PetscOptionsEnum("-pc_jacobi_type","How to construct diagonal matrix","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
26737eeb815SJakub Kruzik   if (flg) {
26837eeb815SJakub Kruzik     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
26937eeb815SJakub Kruzik   }
27037eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
27137eeb815SJakub Kruzik   PetscFunctionReturn(0);
27237eeb815SJakub Kruzik }
27337eeb815SJakub Kruzik 
27437eeb815SJakub Kruzik /* -------------------------------------------------------------------------- */
27537eeb815SJakub Kruzik /*
27637eeb815SJakub Kruzik    PCCreate_Deflation - Creates a Deflation preconditioner context, PC_DeflationPC_Deflation,
27737eeb815SJakub Kruzik    and sets this as the private data within the generic preconditioning
27837eeb815SJakub Kruzik    context, PC, that was created within PCCreate().
27937eeb815SJakub Kruzik 
28037eeb815SJakub Kruzik    Input Parameter:
28137eeb815SJakub Kruzik .  pc - the preconditioner context
28237eeb815SJakub Kruzik 
28337eeb815SJakub Kruzik    Application Interface Routine: PCCreate()
28437eeb815SJakub Kruzik */
28537eeb815SJakub Kruzik 
28637eeb815SJakub Kruzik /*MC
28737eeb815SJakub Kruzik      PCDEFLATION - Deflation (i.e. diagonal scaling preconditioning)
28837eeb815SJakub Kruzik 
28937eeb815SJakub Kruzik    Options Database Key:
29037eeb815SJakub Kruzik +    -pc_jacobi_type <diagonal,rowmax,rowsum> - approach for forming the preconditioner
29137eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
29237eeb815SJakub Kruzik 
29337eeb815SJakub Kruzik    Level: beginner
29437eeb815SJakub Kruzik 
29537eeb815SJakub Kruzik   Concepts: Deflation, diagonal scaling, preconditioners
29637eeb815SJakub Kruzik 
29737eeb815SJakub Kruzik   Notes:
29837eeb815SJakub Kruzik     By using KSPSetPCSide(ksp,PC_SYMMETRIC) or -ksp_pc_side symmetric
29937eeb815SJakub Kruzik          can scale each side of the matrix by the square root of the diagonal entries.
30037eeb815SJakub Kruzik 
30137eeb815SJakub Kruzik          Zero entries along the diagonal are replaced with the value 1.0
30237eeb815SJakub Kruzik 
30337eeb815SJakub Kruzik          See PCPBDEFLATION for a point-block Deflation preconditioner
30437eeb815SJakub Kruzik 
30537eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
306*e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
30737eeb815SJakub Kruzik M*/
30837eeb815SJakub Kruzik 
30937eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
31037eeb815SJakub Kruzik {
31137eeb815SJakub Kruzik   PC_Deflation   *def;
31237eeb815SJakub Kruzik   PetscErrorCode ierr;
31337eeb815SJakub Kruzik 
31437eeb815SJakub Kruzik   PetscFunctionBegin;
31537eeb815SJakub Kruzik   /*
31637eeb815SJakub Kruzik      Creates the private data structure for this preconditioner and
31737eeb815SJakub Kruzik      attach it to the PC object.
31837eeb815SJakub Kruzik   */
31937eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
32037eeb815SJakub Kruzik   pc->data = (void*)def;
32137eeb815SJakub Kruzik 
322*e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
323*e662bc50SJakub Kruzik   def->pre           = PETSC_TRUE;
324*e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
325*e662bc50SJakub Kruzik   def->truenorm      = PETSC_TRUE;
326*e662bc50SJakub Kruzik   def->reductionfact = -1;
327*e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
328*e662bc50SJakub Kruzik   def->spacesize     = 1;
329*e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
330*e662bc50SJakub Kruzik   def->nestedlvl     = 0;
331*e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
332*e662bc50SJakub Kruzik   def->adaptiveconv  = PETSC_FALSE;
333*e662bc50SJakub Kruzik   def->adaptiveconst = 1.0;
33437eeb815SJakub Kruzik 
33537eeb815SJakub Kruzik   /*
33637eeb815SJakub Kruzik       Set the pointers for the functions that are provided above.
33737eeb815SJakub Kruzik       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
33837eeb815SJakub Kruzik       are called, they will automatically call these functions.  Note we
33937eeb815SJakub Kruzik       choose not to provide a couple of these functions since they are
34037eeb815SJakub Kruzik       not needed.
34137eeb815SJakub Kruzik   */
34237eeb815SJakub Kruzik   pc->ops->apply               = PCApply_Deflation;
34337eeb815SJakub Kruzik   pc->ops->applytranspose      = PCApply_Deflation;
34437eeb815SJakub Kruzik   pc->ops->setup               = PCSetUp_Deflation;
34537eeb815SJakub Kruzik   pc->ops->reset               = PCReset_Deflation;
34637eeb815SJakub Kruzik   pc->ops->destroy             = PCDestroy_Deflation;
34737eeb815SJakub Kruzik   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
34837eeb815SJakub Kruzik   pc->ops->view                = 0;
34937eeb815SJakub Kruzik   pc->ops->applyrichardson     = 0;
35037eeb815SJakub Kruzik   pc->ops->applysymmetricleft  = 0;
35137eeb815SJakub Kruzik   pc->ops->applysymmetricright = 0;
35237eeb815SJakub Kruzik 
35337eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
35437eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
355*e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
35637eeb815SJakub Kruzik   PetscFunctionReturn(0);
35737eeb815SJakub Kruzik }
35837eeb815SJakub Kruzik 
359