xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision fcb31d998e9a4b4c0e7669003213f8c0db327f88)
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 
51292e2e67SJakub Kruzik #include <../src/ksp/pc/impls/deflation/deflation.h> /*I "petscpc.h" I*/  /* includes for fortran wrappers */
5237eeb815SJakub Kruzik 
538513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = {
548513960bSJakub Kruzik   "haar",
558513960bSJakub Kruzik   "jacket-haar",
568513960bSJakub Kruzik   "db2",
578513960bSJakub Kruzik   "db4",
588513960bSJakub Kruzik   "db8",
598513960bSJakub Kruzik   "db16",
608513960bSJakub Kruzik   "biorth22",
618513960bSJakub Kruzik   "meyer",
628513960bSJakub Kruzik   "aggregation",
638513960bSJakub Kruzik   "slepc",
648513960bSJakub Kruzik   "slepc-cheap",
658513960bSJakub Kruzik   "user",
668513960bSJakub Kruzik   "DdefSpaceType",
678513960bSJakub Kruzik   "Ddef_SPACE_",
688513960bSJakub Kruzik   0
698513960bSJakub Kruzik };
708513960bSJakub Kruzik 
718513960bSJakub Kruzik 
72e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
73e662bc50SJakub Kruzik {
74e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
75e662bc50SJakub Kruzik   PetscErrorCode ierr;
76e662bc50SJakub Kruzik 
77e662bc50SJakub Kruzik   PetscFunctionBegin;
78e662bc50SJakub Kruzik   if (transpose) {
79e662bc50SJakub Kruzik     def->Wt = W;
80e662bc50SJakub Kruzik     def->W = NULL;
81e662bc50SJakub Kruzik   } else {
82e662bc50SJakub Kruzik     def->W = W;
83e662bc50SJakub Kruzik   }
84e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
85e662bc50SJakub Kruzik   PetscFunctionReturn(0);
86e662bc50SJakub Kruzik }
87e662bc50SJakub Kruzik 
88e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */
89e662bc50SJakub Kruzik /*@
90e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
91e662bc50SJakub Kruzik 
92e662bc50SJakub Kruzik    Logically Collective on PC
93e662bc50SJakub Kruzik 
94e662bc50SJakub Kruzik    Input Parameters:
95e662bc50SJakub Kruzik +  pc - the preconditioner context
96e662bc50SJakub Kruzik .  W  - deflation matrix
97e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
98e662bc50SJakub Kruzik 
99e662bc50SJakub Kruzik    Level: intermediate
100e662bc50SJakub Kruzik 
101e662bc50SJakub Kruzik .seealso: PCDEFLATION
102e662bc50SJakub Kruzik @*/
103e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
104e662bc50SJakub Kruzik {
105e662bc50SJakub Kruzik   PetscErrorCode ierr;
106e662bc50SJakub Kruzik 
107e662bc50SJakub Kruzik   PetscFunctionBegin;
108e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
109e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
110e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
111e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
112e662bc50SJakub Kruzik   PetscFunctionReturn(0);
113e662bc50SJakub Kruzik }
114e662bc50SJakub Kruzik 
1155c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
1165c0c31c5SJakub Kruzik {
1175c0c31c5SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
1185c0c31c5SJakub Kruzik 
1195c0c31c5SJakub Kruzik   PetscFunctionBegin;
12081e2347eSJakub Kruzik   if (current) def->nestedlvl = current;
1215c0c31c5SJakub Kruzik   def->maxnestedlvl = max;
1225c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1235c0c31c5SJakub Kruzik }
1245c0c31c5SJakub Kruzik 
1255c0c31c5SJakub Kruzik /*@
1265c0c31c5SJakub Kruzik    PCDeflationSetMaxLvl - Set maximum level of deflation.
1275c0c31c5SJakub Kruzik 
1285c0c31c5SJakub Kruzik    Logically Collective on PC
1295c0c31c5SJakub Kruzik 
1305c0c31c5SJakub Kruzik    Input Parameters:
1315c0c31c5SJakub Kruzik +  pc  - the preconditioner context
13222b0793eSJakub Kruzik -  max - maximum deflation level
1335c0c31c5SJakub Kruzik 
1345c0c31c5SJakub Kruzik    Level: intermediate
1355c0c31c5SJakub Kruzik 
1365c0c31c5SJakub Kruzik .seealso: PCDEFLATION
1375c0c31c5SJakub Kruzik @*/
1385c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
1395c0c31c5SJakub Kruzik {
1405c0c31c5SJakub Kruzik   PetscErrorCode ierr;
1415c0c31c5SJakub Kruzik 
1425c0c31c5SJakub Kruzik   PetscFunctionBegin;
14322b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1445c0c31c5SJakub Kruzik   PetscValidLogicalCollectiveInt(pc,max,2);
1455c0c31c5SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
1465c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1475c0c31c5SJakub Kruzik }
148e662bc50SJakub Kruzik 
14922b0793eSJakub Kruzik static PetscErrorCode PCDeflationSetPC_Deflation(PC pc,PC apc)
15022b0793eSJakub Kruzik {
15122b0793eSJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
15222b0793eSJakub Kruzik   PetscErrorCode ierr;
15322b0793eSJakub Kruzik 
15422b0793eSJakub Kruzik   PetscFunctionBegin;
15522b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
15622b0793eSJakub Kruzik   def->pc = apc;
15722b0793eSJakub Kruzik   ierr = PetscObjectReference((PetscObject)apc);CHKERRQ(ierr);
15822b0793eSJakub Kruzik   ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)def->pc);CHKERRQ(ierr);
15922b0793eSJakub Kruzik   PetscFunctionReturn(0);
16022b0793eSJakub Kruzik }
16122b0793eSJakub Kruzik 
16222b0793eSJakub Kruzik /*@
16322b0793eSJakub Kruzik    PCDeflationSetPC - Set additional preconditioner.
16422b0793eSJakub Kruzik 
16522b0793eSJakub Kruzik    Logically Collective on PC
16622b0793eSJakub Kruzik 
16722b0793eSJakub Kruzik    Input Parameters:
16822b0793eSJakub Kruzik +  pc  - the preconditioner context
16922b0793eSJakub Kruzik -  apc - additional preconditioner
17022b0793eSJakub Kruzik 
17122b0793eSJakub Kruzik    Level: advanced
17222b0793eSJakub Kruzik 
17322b0793eSJakub Kruzik .seealso: PCDEFLATION
17422b0793eSJakub Kruzik @*/
17522b0793eSJakub Kruzik PetscErrorCode PCDeflationSetPC(PC pc,PC apc)
17622b0793eSJakub Kruzik {
17722b0793eSJakub Kruzik   PetscErrorCode ierr;
17822b0793eSJakub Kruzik 
17922b0793eSJakub Kruzik   PetscFunctionBegin;
18022b0793eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
18122b0793eSJakub Kruzik   PetscValidHeaderSpecific(apc,PC_CLASSID,2);
18222b0793eSJakub Kruzik   PetscCheckSameComm(pc,1,apc,2);
18322b0793eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetPC_C",(PC,PC),(pc,apc));CHKERRQ(ierr);
18422b0793eSJakub Kruzik   PetscFunctionReturn(0);
18522b0793eSJakub Kruzik }
18622b0793eSJakub Kruzik 
1874a99276eSJakub Kruzik static PetscErrorCode PCDeflationGetCoarseKSP_Deflation(PC pc,KSP *ksp)
1884a99276eSJakub Kruzik {
1894a99276eSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
1904a99276eSJakub Kruzik 
1914a99276eSJakub Kruzik   PetscFunctionBegin;
1924a99276eSJakub Kruzik   *ksp = def->WtAWinv;
1934a99276eSJakub Kruzik   PetscFunctionReturn(0);
1944a99276eSJakub Kruzik }
1954a99276eSJakub Kruzik 
1964a99276eSJakub Kruzik /*@
1974a99276eSJakub Kruzik    PCDeflationGetCoarseKSP - Returns a pointer to the coarse problem KSP.
1984a99276eSJakub Kruzik 
1994a99276eSJakub Kruzik    Not Collective
2004a99276eSJakub Kruzik 
2014a99276eSJakub Kruzik    Input Parameters:
2024a99276eSJakub Kruzik .  pc - preconditioner context
2034a99276eSJakub Kruzik 
2044a99276eSJakub Kruzik    Output Parameter:
2054a99276eSJakub Kruzik .  ksp - coarse problem KSP context
2064a99276eSJakub Kruzik 
2074a99276eSJakub Kruzik    Level: developer
2084a99276eSJakub Kruzik 
2094a99276eSJakub Kruzik .seealso: PCDEFLATION
2104a99276eSJakub Kruzik @*/
2114a99276eSJakub Kruzik PetscErrorCode  PCDeflationGetCoarseKSP(PC pc,KSP *ksp)
2124a99276eSJakub Kruzik {
2134a99276eSJakub Kruzik   PetscErrorCode ierr;
2144a99276eSJakub Kruzik 
2154a99276eSJakub Kruzik   PetscFunctionBegin;
2164a99276eSJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2174a99276eSJakub Kruzik   PetscValidPointer(ksp,2);
2184a99276eSJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationGetCoarseKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr);
2194a99276eSJakub Kruzik   PetscFunctionReturn(0);
2204a99276eSJakub Kruzik }
2214a99276eSJakub Kruzik 
22237eeb815SJakub Kruzik /*
223b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
224b27c8092SJakub Kruzik */
225b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
226b27c8092SJakub Kruzik {
227b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
228b27c8092SJakub Kruzik   Mat              A;
229b27c8092SJakub Kruzik   Vec              r,w1,w2;
230b27c8092SJakub Kruzik   PetscBool        nonzero;
231b27c8092SJakub Kruzik   PetscErrorCode   ierr;
23237eeb815SJakub Kruzik 
233b27c8092SJakub Kruzik   PetscFunctionBegin;
234b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
235b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
236b27c8092SJakub Kruzik   r  = def->work;
237b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
23837eeb815SJakub Kruzik 
239b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
240b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
241b27c8092SJakub Kruzik   if (nonzero) {
242b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);               /*    r  <- b - Ax              */
243b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
244b27c8092SJakub Kruzik   } else {
245b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                 /*    r  <- b (x is 0)          */
246b27c8092SJakub Kruzik   }
247b27c8092SJakub Kruzik 
248b27c8092SJakub Kruzik   ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
249b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);   /*    w2 <- (W'*A*W)^{-1}*w1    */
250b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);           /*    r  <- W*w2                */
251b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
252b27c8092SJakub Kruzik   PetscFunctionReturn(0);
253b27c8092SJakub Kruzik }
25437eeb815SJakub Kruzik 
255f8bf229cSJakub Kruzik /*
256f8bf229cSJakub Kruzik   if (def->correct) {
25722b0793eSJakub Kruzik     z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*M^{-1}*r
258f8bf229cSJakub Kruzik   } else {
25922b0793eSJakub Kruzik     z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*M^{-1}*r
260f8bf229cSJakub Kruzik   }
26137eeb815SJakub Kruzik */
262f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r,Vec z)
263f8bf229cSJakub Kruzik {
264f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
265f8bf229cSJakub Kruzik   Mat              A;
266f8bf229cSJakub Kruzik   Vec              u,w1,w2;
267f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
268f8bf229cSJakub Kruzik 
269f8bf229cSJakub Kruzik   PetscFunctionBegin;
270f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
271f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
272f8bf229cSJakub Kruzik   u  = def->work;
273f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
274f8bf229cSJakub Kruzik 
27522b0793eSJakub Kruzik   ierr = PCApply(def->pc,r,z);CHKERRQ(ierr);
27622b0793eSJakub Kruzik   if (!def->init) {
277f8bf229cSJakub Kruzik     if (!def->AW) {
27822b0793eSJakub Kruzik       ierr = MatMult(A,z,u);CHKERRQ(ierr);                                        /*    u  <- A*z                 */
279*fcb31d99SJakub Kruzik       if (def->correct) ierr = VecAXPY(u,-1.0*def->correctfact,z);CHKERRQ(ierr);  /*    u  <- A*z -z              */
280f8bf229cSJakub Kruzik       ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr);                         /*    w1 <- W'*u                */
281f8bf229cSJakub Kruzik     } else {
28222b0793eSJakub Kruzik       /* ONLY if A SYM */
28322b0793eSJakub Kruzik       ierr = MatMultTranspose(def->AW,z,w1);CHKERRQ(ierr);                        /*    w1  <- W'*A*r             */
284f8bf229cSJakub Kruzik       if (def->correct) {
28522b0793eSJakub Kruzik         ierr = MatMultTranspose(def->W,z,w2);CHKERRQ(ierr);                       /*    w2 <- W'*u                */
286*fcb31d99SJakub Kruzik         ierr = VecAXPY(w1,-1.0*def->correctfact,w2);CHKERRQ(ierr);                /*    w1 <- w1 - w2             */
287f8bf229cSJakub Kruzik       }
288f8bf229cSJakub Kruzik     }
289f8bf229cSJakub Kruzik     ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);                            /*    w2 <- (W'*A*W)^{-1}*w1    */
290f8bf229cSJakub Kruzik     ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                                    /*    u  <- W*w2                */
29122b0793eSJakub Kruzik     ierr = VecAXPY(z,-1.0,u);CHKERRQ(ierr);                                       /*    z  <- z - u               */
29222b0793eSJakub Kruzik   }
293f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
294f8bf229cSJakub Kruzik }
295f8bf229cSJakub Kruzik 
29637eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
29737eeb815SJakub Kruzik {
298409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
299409a357bSJakub Kruzik   KSP              innerksp;
300409a357bSJakub Kruzik   PC               pcinner;
301409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
302409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
303409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
304409a357bSJakub Kruzik   MatCompositeType ctype;
305409a357bSJakub Kruzik   MPI_Comm         comm;
306409a357bSJakub Kruzik   const char       *prefix;
30737eeb815SJakub Kruzik   PetscErrorCode   ierr;
30837eeb815SJakub Kruzik 
30937eeb815SJakub Kruzik   PetscFunctionBegin;
310409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
311409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
312f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
313f8bf229cSJakub Kruzik 
314f8bf229cSJakub Kruzik   /* compute a deflation space */
315409a357bSJakub Kruzik   if (def->W || def->Wt) {
316409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
317409a357bSJakub Kruzik   } else {
318e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
31937eeb815SJakub Kruzik   }
32037eeb815SJakub Kruzik 
321409a357bSJakub Kruzik   /* nested deflation */
322409a357bSJakub Kruzik   if (def->W) {
323409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
324409a357bSJakub Kruzik     if (match) {
325409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
326409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
32737eeb815SJakub Kruzik     }
328409a357bSJakub Kruzik   } else {
329409a357bSJakub Kruzik     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
330409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
331409a357bSJakub Kruzik     if (match) {
332409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
333409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
334409a357bSJakub Kruzik     }
335409a357bSJakub Kruzik     transp = PETSC_TRUE;
336409a357bSJakub Kruzik   }
337409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
338409a357bSJakub Kruzik     if (!transp) {
33928da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
34028da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
341409a357bSJakub Kruzik         for (i=0; i<size; i++) {
342409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
343409a357bSJakub Kruzik         }
344409a357bSJakub Kruzik         size -= 1;
345409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
346409a357bSJakub Kruzik         def->W = mats[size];
347409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
348409a357bSJakub Kruzik         if (size > 1) {
349409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
350409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
351409a357bSJakub Kruzik         } else {
352409a357bSJakub Kruzik           nextDef = mats[0];
353409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
354409a357bSJakub Kruzik         }
35528da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
356409a357bSJakub Kruzik       } else {
357409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
358409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
359409a357bSJakub Kruzik       }
36028da0170SJakub Kruzik     } else {
36128da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
36228da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
36328da0170SJakub Kruzik         for (i=0; i<size; i++) {
36428da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
36528da0170SJakub Kruzik         }
36628da0170SJakub Kruzik         size -= 1;
36728da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
36828da0170SJakub Kruzik         def->Wt = mats[0];
36928da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
37028da0170SJakub Kruzik         if (size > 1) {
37128da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
37228da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
37328da0170SJakub Kruzik         } else {
37428da0170SJakub Kruzik           nextDef = mats[1];
37528da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
376409a357bSJakub Kruzik         }
377409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
37828da0170SJakub Kruzik       } else {
37928da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
38028da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
38128da0170SJakub Kruzik       }
38228da0170SJakub Kruzik     }
38328da0170SJakub Kruzik   }
38428da0170SJakub Kruzik 
38528da0170SJakub Kruzik   if (transp) {
38628da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
38728da0170SJakub Kruzik     ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
388409a357bSJakub Kruzik   }
389409a357bSJakub Kruzik 
39022b0793eSJakub Kruzik 
39122b0793eSJakub Kruzik   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
39222b0793eSJakub Kruzik 
393409a357bSJakub Kruzik   /* setup coarse problem */
394409a357bSJakub Kruzik   if (!def->WtAWinv) {
39528da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
396409a357bSJakub Kruzik     if (!def->WtAW) {
397409a357bSJakub Kruzik       /* TODO add implicit product version ? */
398409a357bSJakub Kruzik       if (!def->AW) {
399409a357bSJakub Kruzik         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
400409a357bSJakub Kruzik       } else {
401409a357bSJakub Kruzik         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
402409a357bSJakub Kruzik       }
403409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
404409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
405409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
406409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
407409a357bSJakub Kruzik       /* Check WtAW is not sigular */
408409a357bSJakub Kruzik       PetscReal *norms;
409409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
410409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
411409a357bSJakub Kruzik       for (i=0; i<m; i++) {
412409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
413409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
414409a357bSJakub Kruzik         }
415409a357bSJakub Kruzik       }
416409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
417409a357bSJakub Kruzik #endif
418409a357bSJakub Kruzik     } else {
419409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
420409a357bSJakub Kruzik     }
421409a357bSJakub Kruzik     /* TODO use MATINV */
422409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
423409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
424409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
425557a2f7dSJakub Kruzik     /* Setup KSP and PC */
426557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
427557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
428557a2f7dSJakub Kruzik       /* set default KSPtype */
429557a2f7dSJakub Kruzik       if (!def->ksptype) {
430557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
431557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
432557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
433557a2f7dSJakub Kruzik         }
434557a2f7dSJakub Kruzik       }
435557a2f7dSJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
436557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
437557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
438557a2f7dSJakub Kruzik       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
439557a2f7dSJakub Kruzik       /* inherit options TODO if not set */
440557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
441557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
442557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
443557a2f7dSJakub Kruzik     } else { /* the last level */
444557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
445409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
446ac66f006SJakub Kruzik       /* ugly hack to not have overwritten PCTELESCOPE */
447ac66f006SJakub Kruzik       if (prefix) {
448ac66f006SJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
449ac66f006SJakub Kruzik       }
45022b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
451ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
452409a357bSJakub Kruzik       /* Reduction factor choice */
453409a357bSJakub Kruzik       red = def->reductionfact;
454409a357bSJakub Kruzik       if (red < 0) {
455409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
456409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
457409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
458409a357bSJakub Kruzik         if (match) red = commsize;
459409a357bSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
460409a357bSJakub Kruzik       }
461409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
462ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
463409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
464ac66f006SJakub Kruzik       if (innerksp) {
465409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
466409a357bSJakub Kruzik         /* TODO Cholesky if flgspd? */
467ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
468409a357bSJakub Kruzik         //TODO remove explicit matSolverPackage
469409a357bSJakub Kruzik         if (commsize == red) {
470ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
471409a357bSJakub Kruzik         } else {
472ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
473409a357bSJakub Kruzik         }
474409a357bSJakub Kruzik       }
475557a2f7dSJakub Kruzik     }
476557a2f7dSJakub Kruzik 
477557a2f7dSJakub Kruzik     if (innerksp) {
478409a357bSJakub Kruzik       /* TODO use def_[lvl]_ if lvl > 0? */
479409a357bSJakub Kruzik       if (prefix) {
480409a357bSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
481409a357bSJakub Kruzik       }
48222b0793eSJakub Kruzik       ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
483557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
484557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
485ac66f006SJakub Kruzik     }
486409a357bSJakub Kruzik   }
487557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
488557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
489409a357bSJakub Kruzik 
49022b0793eSJakub Kruzik   if (!def->pc) {
49122b0793eSJakub Kruzik     ierr = PCCreate(comm,&def->pc);CHKERRQ(ierr);
49222b0793eSJakub Kruzik     ierr = PCSetOperators(def->pc,Amat,Amat);CHKERRQ(ierr);
49322b0793eSJakub Kruzik     ierr = PCSetType(def->pc,PCNONE);CHKERRQ(ierr);
49422b0793eSJakub Kruzik     if (prefix) {
49522b0793eSJakub Kruzik       ierr = PCSetOptionsPrefix(def->pc,prefix);CHKERRQ(ierr);
49622b0793eSJakub Kruzik     }
49722b0793eSJakub Kruzik     ierr = PCAppendOptionsPrefix(def->pc,"def_pc_");CHKERRQ(ierr);
49822b0793eSJakub Kruzik     ierr = PCSetFromOptions(def->pc);CHKERRQ(ierr);
49922b0793eSJakub Kruzik     ierr = PCSetUp(def->pc);CHKERRQ(ierr);
50022b0793eSJakub Kruzik   }
50122b0793eSJakub Kruzik 
502f8bf229cSJakub Kruzik   /* create work vecs */
503f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
504f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
50537eeb815SJakub Kruzik   PetscFunctionReturn(0);
50637eeb815SJakub Kruzik }
507b30d4775SJakub Kruzik 
50837eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
50937eeb815SJakub Kruzik {
510b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
51137eeb815SJakub Kruzik   PetscErrorCode    ierr;
51237eeb815SJakub Kruzik 
51337eeb815SJakub Kruzik   PetscFunctionBegin;
514b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
515b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
516b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
517b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
518b30d4775SJakub Kruzik   ierr = MatDestroy(&def->AW);CHKERRQ(ierr);
519b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
520b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
52122b0793eSJakub Kruzik   ierr = PCDestroy(&def->pc);CHKERRQ(ierr);
52237eeb815SJakub Kruzik   PetscFunctionReturn(0);
52337eeb815SJakub Kruzik }
52437eeb815SJakub Kruzik 
52537eeb815SJakub Kruzik /*
52637eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
52737eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
52837eeb815SJakub Kruzik 
52937eeb815SJakub Kruzik    Input Parameter:
53037eeb815SJakub Kruzik .  pc - the preconditioner context
53137eeb815SJakub Kruzik 
53237eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
53337eeb815SJakub Kruzik */
53437eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
53537eeb815SJakub Kruzik {
53637eeb815SJakub Kruzik   PetscErrorCode ierr;
53737eeb815SJakub Kruzik 
53837eeb815SJakub Kruzik   PetscFunctionBegin;
53937eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
54037eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
54137eeb815SJakub Kruzik   PetscFunctionReturn(0);
54237eeb815SJakub Kruzik }
54337eeb815SJakub Kruzik 
5448513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
5458513960bSJakub Kruzik {
5468513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
5478513960bSJakub Kruzik   PetscBool         iascii;
5488513960bSJakub Kruzik   PetscErrorCode    ierr;
5498513960bSJakub Kruzik 
5508513960bSJakub Kruzik   PetscFunctionBegin;
5518513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5528513960bSJakub Kruzik   if (iascii) {
5538513960bSJakub Kruzik     //if (cg->adaptiveconv) {
5548513960bSJakub Kruzik     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
5558513960bSJakub Kruzik     //}
5568513960bSJakub Kruzik     if (def->correct) {
5578513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
5588513960bSJakub Kruzik     }
5598513960bSJakub Kruzik     if (!def->nestedlvl) {
5608513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
5618513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
5628513960bSJakub Kruzik     }
5638513960bSJakub Kruzik 
56422b0793eSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Additional PC\n");CHKERRQ(ierr);
56522b0793eSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
56622b0793eSJakub Kruzik     ierr = PCView(def->pc,viewer);CHKERRQ(ierr);
56722b0793eSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
56822b0793eSJakub Kruzik 
5698513960bSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
5708513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
5718513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
5728513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
5738513960bSJakub Kruzik   }
5748513960bSJakub Kruzik   PetscFunctionReturn(0);
5758513960bSJakub Kruzik }
5768513960bSJakub Kruzik 
57737eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
57837eeb815SJakub Kruzik {
579880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
58037eeb815SJakub Kruzik   PetscErrorCode    ierr;
58137eeb815SJakub Kruzik 
58237eeb815SJakub Kruzik   PetscFunctionBegin;
58337eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
584880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
585880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
586880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
587880d05ceSJakub Kruzik //TODO add set function and fix manpages
588880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr);
589*fcb31d99SJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correct","Add coarse problem correction Q to P","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr);
590*fcb31d99SJakub Kruzik   ierr = PetscOptionsReal("-pc_deflation_correct_val","Set multiple of Q to use as coarse problem correction","PCDeflation",def->correctfact,&def->correctfact,NULL);CHKERRQ(ierr);
591880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
592880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
59337eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
59437eeb815SJakub Kruzik   PetscFunctionReturn(0);
59537eeb815SJakub Kruzik }
59637eeb815SJakub Kruzik 
59737eeb815SJakub Kruzik /*MC
598e361f8b5SJakub Kruzik      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
599e361f8b5SJakub Kruzik      or to a predefined value
60037eeb815SJakub Kruzik 
60137eeb815SJakub Kruzik    Options Database Key:
602e361f8b5SJakub Kruzik +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
60337eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
60437eeb815SJakub Kruzik 
60537eeb815SJakub Kruzik    Level: beginner
60637eeb815SJakub Kruzik 
60737eeb815SJakub Kruzik   Notes:
608e361f8b5SJakub Kruzik     todo
60937eeb815SJakub Kruzik 
61037eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
611e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
61237eeb815SJakub Kruzik M*/
61337eeb815SJakub Kruzik 
61437eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
61537eeb815SJakub Kruzik {
61637eeb815SJakub Kruzik   PC_Deflation   *def;
61737eeb815SJakub Kruzik   PetscErrorCode ierr;
61837eeb815SJakub Kruzik 
61937eeb815SJakub Kruzik   PetscFunctionBegin;
62037eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
62137eeb815SJakub Kruzik   pc->data = (void*)def;
62237eeb815SJakub Kruzik 
623e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
624e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
625*fcb31d99SJakub Kruzik   def->correctfact   = 1.0;
626e662bc50SJakub Kruzik   def->reductionfact = -1;
627e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
628e662bc50SJakub Kruzik   def->spacesize     = 1;
629e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
630e662bc50SJakub Kruzik   def->nestedlvl     = 0;
631e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
63237eeb815SJakub Kruzik 
63337eeb815SJakub Kruzik   /*
63437eeb815SJakub Kruzik       Set the pointers for the functions that are provided above.
63537eeb815SJakub Kruzik       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
63637eeb815SJakub Kruzik       are called, they will automatically call these functions.  Note we
63737eeb815SJakub Kruzik       choose not to provide a couple of these functions since they are
63837eeb815SJakub Kruzik       not needed.
63937eeb815SJakub Kruzik   */
64037eeb815SJakub Kruzik   pc->ops->apply               = PCApply_Deflation;
64137eeb815SJakub Kruzik   pc->ops->applytranspose      = PCApply_Deflation;
642b27c8092SJakub Kruzik   pc->ops->presolve            = PCPreSolve_Deflation;
643b27c8092SJakub Kruzik   pc->ops->postsolve           = 0;
64437eeb815SJakub Kruzik   pc->ops->setup               = PCSetUp_Deflation;
64537eeb815SJakub Kruzik   pc->ops->reset               = PCReset_Deflation;
64637eeb815SJakub Kruzik   pc->ops->destroy             = PCDestroy_Deflation;
64737eeb815SJakub Kruzik   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
6488513960bSJakub Kruzik   pc->ops->view                = PCView_Deflation;
64937eeb815SJakub Kruzik   pc->ops->applyrichardson     = 0;
65037eeb815SJakub Kruzik   pc->ops->applysymmetricleft  = 0;
65137eeb815SJakub Kruzik   pc->ops->applysymmetricright = 0;
65237eeb815SJakub Kruzik 
653e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
6545c0c31c5SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
65522b0793eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetPC_C",PCDeflationSetPC_Deflation);CHKERRQ(ierr);
6564a99276eSJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetCoarseKSP_C",PCDeflationGetCoarseKSP_Deflation);CHKERRQ(ierr);
65737eeb815SJakub Kruzik   PetscFunctionReturn(0);
65837eeb815SJakub Kruzik }
65937eeb815SJakub Kruzik 
660