xref: /petsc/src/ksp/pc/impls/deflation/deflation.c (revision 880d05ce2aa56e94179dba63c9fa9a8b0cd598c1)
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 
5337eeb815SJakub Kruzik const char *const PCDeflationTypes[]    = {"INIT","PRE","POST","PCDeflationType","PC_DEFLATION_",0};
5437eeb815SJakub Kruzik 
558513960bSJakub Kruzik const char *const PCDeflationSpaceTypes[] = {
568513960bSJakub Kruzik   "haar",
578513960bSJakub Kruzik   "jacket-haar",
588513960bSJakub Kruzik   "db2",
598513960bSJakub Kruzik   "db4",
608513960bSJakub Kruzik   "db8",
618513960bSJakub Kruzik   "db16",
628513960bSJakub Kruzik   "biorth22",
638513960bSJakub Kruzik   "meyer",
648513960bSJakub Kruzik   "aggregation",
658513960bSJakub Kruzik   "slepc",
668513960bSJakub Kruzik   "slepc-cheap",
678513960bSJakub Kruzik   "user",
688513960bSJakub Kruzik   "DdefSpaceType",
698513960bSJakub Kruzik   "Ddef_SPACE_",
708513960bSJakub Kruzik   0
718513960bSJakub Kruzik };
728513960bSJakub Kruzik 
738513960bSJakub Kruzik 
74e662bc50SJakub Kruzik static PetscErrorCode PCDeflationSetSpace_Deflation(PC pc,Mat W,PetscBool transpose)
75e662bc50SJakub Kruzik {
76e662bc50SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
77e662bc50SJakub Kruzik   PetscErrorCode ierr;
78e662bc50SJakub Kruzik 
79e662bc50SJakub Kruzik   PetscFunctionBegin;
80e662bc50SJakub Kruzik   if (transpose) {
81e662bc50SJakub Kruzik     def->Wt = W;
82e662bc50SJakub Kruzik     def->W = NULL;
83e662bc50SJakub Kruzik   } else {
84e662bc50SJakub Kruzik     def->W = W;
85e662bc50SJakub Kruzik   }
86e662bc50SJakub Kruzik   ierr = PetscObjectReference((PetscObject)W);CHKERRQ(ierr);
87e662bc50SJakub Kruzik   PetscFunctionReturn(0);
88e662bc50SJakub Kruzik }
89e662bc50SJakub Kruzik 
90e662bc50SJakub Kruzik /* TODO create PCDeflationSetSpaceTranspose? */
91e662bc50SJakub Kruzik /*@
92e662bc50SJakub Kruzik    PCDeflationSetSpace - Set deflation space matrix (or its transpose).
93e662bc50SJakub Kruzik 
94e662bc50SJakub Kruzik    Logically Collective on PC
95e662bc50SJakub Kruzik 
96e662bc50SJakub Kruzik    Input Parameters:
97e662bc50SJakub Kruzik +  pc - the preconditioner context
98e662bc50SJakub Kruzik .  W  - deflation matrix
99e662bc50SJakub Kruzik -  tranpose - indicates that W is an explicit transpose of the deflation matrix
100e662bc50SJakub Kruzik 
101e662bc50SJakub Kruzik    Level: intermediate
102e662bc50SJakub Kruzik 
103e662bc50SJakub Kruzik .seealso: PCDEFLATION
104e662bc50SJakub Kruzik @*/
105e662bc50SJakub Kruzik PetscErrorCode PCDeflationSetSpace(PC pc,Mat W,PetscBool transpose)
106e662bc50SJakub Kruzik {
107e662bc50SJakub Kruzik   PetscErrorCode ierr;
108e662bc50SJakub Kruzik 
109e662bc50SJakub Kruzik   PetscFunctionBegin;
110e662bc50SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
111e662bc50SJakub Kruzik   PetscValidHeaderSpecific(W,MAT_CLASSID,2);
112e662bc50SJakub Kruzik   PetscValidLogicalCollectiveBool(pc,transpose,3);
113e662bc50SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetSpace_C",(PC,Mat,PetscBool),(pc,W,transpose));CHKERRQ(ierr);
114e662bc50SJakub Kruzik   PetscFunctionReturn(0);
115e662bc50SJakub Kruzik }
116e662bc50SJakub Kruzik 
1175c0c31c5SJakub Kruzik static PetscErrorCode PCDeflationSetLvl_Deflation(PC pc,PetscInt current,PetscInt max)
1185c0c31c5SJakub Kruzik {
1195c0c31c5SJakub Kruzik   PC_Deflation   *def = (PC_Deflation*)pc->data;
1205c0c31c5SJakub Kruzik 
1215c0c31c5SJakub Kruzik   PetscFunctionBegin;
1225c0c31c5SJakub Kruzik   def->nestedlvl = current;
1235c0c31c5SJakub Kruzik   def->maxnestedlvl = max;
1245c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1255c0c31c5SJakub Kruzik }
1265c0c31c5SJakub Kruzik 
1275c0c31c5SJakub Kruzik /*@
1285c0c31c5SJakub Kruzik    PCDeflationSetMaxLvl - Set maximum level of deflation.
1295c0c31c5SJakub Kruzik 
1305c0c31c5SJakub Kruzik    Logically Collective on PC
1315c0c31c5SJakub Kruzik 
1325c0c31c5SJakub Kruzik    Input Parameters:
1335c0c31c5SJakub Kruzik +  pc  - the preconditioner context
1345c0c31c5SJakub Kruzik .  max - maximum deflation level
1355c0c31c5SJakub Kruzik 
1365c0c31c5SJakub Kruzik    Level: intermediate
1375c0c31c5SJakub Kruzik 
1385c0c31c5SJakub Kruzik .seealso: PCDEFLATION
1395c0c31c5SJakub Kruzik @*/
1405c0c31c5SJakub Kruzik PetscErrorCode PCDeflationSetMaxLvl(PC pc,PetscInt max)
1415c0c31c5SJakub Kruzik {
1425c0c31c5SJakub Kruzik   PetscErrorCode ierr;
1435c0c31c5SJakub Kruzik 
1445c0c31c5SJakub Kruzik   PetscFunctionBegin;
1455c0c31c5SJakub Kruzik   PetscValidLogicalCollectiveInt(pc,max,2);
1465c0c31c5SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetLvl_C",(PC,PetscInt,PetscInt),(pc,0,max));CHKERRQ(ierr);
1475c0c31c5SJakub Kruzik   PetscFunctionReturn(0);
1485c0c31c5SJakub Kruzik }
149e662bc50SJakub Kruzik 
15037eeb815SJakub Kruzik /*
151b27c8092SJakub Kruzik   TODO CP corection?
152b27c8092SJakub Kruzik   x <- x + W*(W'*A*W)^{-1}*W'*r  = x + Q*r
153b27c8092SJakub Kruzik */
154b27c8092SJakub Kruzik static PetscErrorCode PCPreSolve_Deflation(PC pc,KSP ksp,Vec b, Vec x)
155b27c8092SJakub Kruzik {
156b27c8092SJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
157b27c8092SJakub Kruzik   Mat              A;
158b27c8092SJakub Kruzik   Vec              r,w1,w2;
159b27c8092SJakub Kruzik   PetscBool        nonzero;
160b27c8092SJakub Kruzik   PetscErrorCode   ierr;
16137eeb815SJakub Kruzik 
162b27c8092SJakub Kruzik   PetscFunctionBegin;
163b27c8092SJakub Kruzik   w1 = def->workcoarse[0];
164b27c8092SJakub Kruzik   w2 = def->workcoarse[1];
165b27c8092SJakub Kruzik   r  = def->work;
166b27c8092SJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
16737eeb815SJakub Kruzik 
168b27c8092SJakub Kruzik   ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr);
169b27c8092SJakub Kruzik   ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
170b27c8092SJakub Kruzik   if (nonzero) {
171b27c8092SJakub Kruzik     ierr = MatMult(A,x,r);CHKERRQ(ierr);               /*    r  <- b - Ax              */
172b27c8092SJakub Kruzik     ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr);
173b27c8092SJakub Kruzik   } else {
174b27c8092SJakub Kruzik     ierr = VecCopy(b,r);CHKERRQ(ierr);                 /*    r  <- b (x is 0)          */
175b27c8092SJakub Kruzik   }
176b27c8092SJakub Kruzik 
177b27c8092SJakub Kruzik   ierr = MatMultTranspose(def->W,r,w1);CHKERRQ(ierr);  /*    w1 <- W'*r                */
178b27c8092SJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);   /*    w2 <- (W'*A*W)^{-1}*w1    */
179b27c8092SJakub Kruzik   ierr = MatMult(def->W,w2,r);CHKERRQ(ierr);           /*    r  <- W*w2                */
180b27c8092SJakub Kruzik   ierr = VecAYPX(x,1.0,r);CHKERRQ(ierr);
181b27c8092SJakub Kruzik   PetscFunctionReturn(0);
182b27c8092SJakub Kruzik }
18337eeb815SJakub Kruzik 
184f8bf229cSJakub Kruzik /*
185f8bf229cSJakub Kruzik   if (def->correct) {
186f8bf229cSJakub Kruzik     z <- r - W*(W'*A*W)^{-1}*W'*(A*r -r) = (P-Q)*r
187f8bf229cSJakub Kruzik   } else {
188f8bf229cSJakub Kruzik     z <- r - W*(W'*A*W)^{-1}*W'*A*r = P*r
189f8bf229cSJakub Kruzik   }
19037eeb815SJakub Kruzik */
191f8bf229cSJakub Kruzik static PetscErrorCode PCApply_Deflation(PC pc,Vec r, Vec z)
192f8bf229cSJakub Kruzik {
193f8bf229cSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
194f8bf229cSJakub Kruzik   Mat              A;
195f8bf229cSJakub Kruzik   Vec              u,w1,w2;
196f8bf229cSJakub Kruzik   PetscErrorCode   ierr;
197f8bf229cSJakub Kruzik 
198f8bf229cSJakub Kruzik   PetscFunctionBegin;
199f8bf229cSJakub Kruzik   w1 = def->workcoarse[0];
200f8bf229cSJakub Kruzik   w2 = def->workcoarse[1];
201f8bf229cSJakub Kruzik   u  = def->work;
202f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr);
203f8bf229cSJakub Kruzik 
204f8bf229cSJakub Kruzik   if (!def->AW) {
205f8bf229cSJakub Kruzik     ierr = MatMult(A,r,u);CHKERRQ(ierr);                       /*    u  <- A*r                 */
206f8bf229cSJakub Kruzik     if (def->correct) ierr = VecAXPY(u,-1.0,r);CHKERRQ(ierr);  /*    u  <- A*r -r              */
207f8bf229cSJakub Kruzik     ierr = MatMultTranspose(def->W,u,w1);CHKERRQ(ierr);        /*    w1 <- W'*u                */
208f8bf229cSJakub Kruzik   } else {
209f8bf229cSJakub Kruzik     ierr = MatMultTranspose(def->AW,u,w1);CHKERRQ(ierr);       /*    u  <- A*r                 */
210f8bf229cSJakub Kruzik     if (def->correct) {
211f8bf229cSJakub Kruzik       ierr = MatMultTranspose(def->W,r,w2);CHKERRQ(ierr);      /*    w2 <- W'*u                */
212f8bf229cSJakub Kruzik       ierr = VecAXPY(w1,-1.0,w2);CHKERRQ(ierr);                /*    w1 <- w1 - w2             */
213f8bf229cSJakub Kruzik     }
214f8bf229cSJakub Kruzik   }
215f8bf229cSJakub Kruzik   ierr = KSPSolve(def->WtAWinv,w1,w2);CHKERRQ(ierr);           /*    w2 <- (W'*A*W)^{-1}*w1    */
216f8bf229cSJakub Kruzik   ierr = MatMult(def->W,w2,u);CHKERRQ(ierr);                   /*    u  <- W*w2                */
217f8bf229cSJakub Kruzik   ierr = VecWAXPY(z,-1.0,u,r);CHKERRQ(ierr);                   /*    z  <- r - u               */
218f8bf229cSJakub Kruzik   PetscFunctionReturn(0);
219f8bf229cSJakub Kruzik }
220f8bf229cSJakub Kruzik 
221b27c8092SJakub Kruzik static PetscErrorCode  PCDeflationSetType_Deflation(PC pc,PCDeflationType type)
222b27c8092SJakub Kruzik {
223b27c8092SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
224b27c8092SJakub Kruzik 
225b27c8092SJakub Kruzik   PetscFunctionBegin;
226b27c8092SJakub Kruzik   def->init = PETSC_FALSE;
227b27c8092SJakub Kruzik   def->pre = PETSC_FALSE;
228b27c8092SJakub Kruzik   if (type == PC_DEFLATION_POST) {
229b27c8092SJakub Kruzik     //pc->ops->postsolve = PCPostSolve_Deflation;
230b27c8092SJakub Kruzik     pc->ops->presolve = 0;
231b27c8092SJakub Kruzik   } else {
232b27c8092SJakub Kruzik     pc->ops->presolve = PCPreSolve_Deflation;
233b27c8092SJakub Kruzik     pc->ops->postsolve = 0;
234b27c8092SJakub Kruzik     if (type == PC_DEFLATION_INIT) {
235b27c8092SJakub Kruzik       def->init = PETSC_TRUE;
236b27c8092SJakub Kruzik       pc->ops->apply = 0;
237b27c8092SJakub Kruzik     } else {
238b27c8092SJakub Kruzik       def->pre  = PETSC_TRUE;
239b27c8092SJakub Kruzik     }
240b27c8092SJakub Kruzik   }
241b27c8092SJakub Kruzik   PetscFunctionReturn(0);
242b27c8092SJakub Kruzik }
243b27c8092SJakub Kruzik 
244b27c8092SJakub Kruzik /*@
245b27c8092SJakub Kruzik    PCDeflationSetType - Causes the deflation preconditioner to use only a special
246b27c8092SJakub Kruzik     initial gues or pre/post solve solution update
247b27c8092SJakub Kruzik 
248b27c8092SJakub Kruzik    Logically Collective on PC
249b27c8092SJakub Kruzik 
250b27c8092SJakub Kruzik    Input Parameters:
251b27c8092SJakub Kruzik +  pc - the preconditioner context
252b27c8092SJakub Kruzik -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
253b27c8092SJakub Kruzik 
254b27c8092SJakub Kruzik    Options Database Key:
255b27c8092SJakub Kruzik .  -pc_deflation_type <pre,init,post>
256b27c8092SJakub Kruzik 
257b27c8092SJakub Kruzik    Level: intermediate
258b27c8092SJakub Kruzik 
259b27c8092SJakub Kruzik    Concepts: Deflation preconditioner
260b27c8092SJakub Kruzik 
261b27c8092SJakub Kruzik .seealso: PCDeflationGetType()
262b27c8092SJakub Kruzik @*/
263b27c8092SJakub Kruzik PetscErrorCode  PCDeflationSetType(PC pc,PCDeflationType type)
264b27c8092SJakub Kruzik {
265b27c8092SJakub Kruzik   PetscErrorCode ierr;
266b27c8092SJakub Kruzik 
267b27c8092SJakub Kruzik   PetscFunctionBegin;
268b27c8092SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
269b27c8092SJakub Kruzik   ierr = PetscTryMethod(pc,"PCDeflationSetType_C",(PC,PCDeflationType),(pc,type));CHKERRQ(ierr);
270b27c8092SJakub Kruzik   PetscFunctionReturn(0);
271b27c8092SJakub Kruzik }
272b27c8092SJakub Kruzik 
273b27c8092SJakub Kruzik static PetscErrorCode  PCDeflationGetType_Deflation(PC pc,PCDeflationType *type)
274b27c8092SJakub Kruzik {
275b27c8092SJakub Kruzik   PC_Deflation *def = (PC_Deflation*)pc->data;
276b27c8092SJakub Kruzik 
277b27c8092SJakub Kruzik   PetscFunctionBegin;
278b27c8092SJakub Kruzik   if (def->init) {
279b27c8092SJakub Kruzik     *type = PC_DEFLATION_INIT;
280b27c8092SJakub Kruzik   } else if (def->pre) {
281b27c8092SJakub Kruzik     *type = PC_DEFLATION_PRE;
282b27c8092SJakub Kruzik   } else {
283b27c8092SJakub Kruzik     *type = PC_DEFLATION_POST;
284b27c8092SJakub Kruzik   }
285b27c8092SJakub Kruzik   PetscFunctionReturn(0);
286b27c8092SJakub Kruzik }
287b27c8092SJakub Kruzik 
288b27c8092SJakub Kruzik /*@
289b27c8092SJakub Kruzik    PCDeflationGetType - Gets how the diagonal matrix is produced for the preconditioner
290b27c8092SJakub Kruzik 
291b27c8092SJakub Kruzik    Not Collective on PC
292b27c8092SJakub Kruzik 
293b27c8092SJakub Kruzik    Input Parameter:
294b27c8092SJakub Kruzik .  pc - the preconditioner context
295b27c8092SJakub Kruzik 
296b27c8092SJakub Kruzik    Output Parameter:
297b27c8092SJakub Kruzik -  type - PC_DEFLATION_PRE, PC_DEFLATION_INIT, PC_DEFLATION_POST
298b27c8092SJakub Kruzik 
299b27c8092SJakub Kruzik    Level: intermediate
300b27c8092SJakub Kruzik 
301b27c8092SJakub Kruzik    Concepts: Deflation preconditioner
302b27c8092SJakub Kruzik 
303b27c8092SJakub Kruzik .seealso: PCDeflationSetType()
304b27c8092SJakub Kruzik @*/
305b27c8092SJakub Kruzik PetscErrorCode  PCDeflationGetType(PC pc,PCDeflationType *type)
306b27c8092SJakub Kruzik {
307b27c8092SJakub Kruzik   PetscErrorCode ierr;
308b27c8092SJakub Kruzik 
309b27c8092SJakub Kruzik   PetscFunctionBegin;
310b27c8092SJakub Kruzik   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
311b27c8092SJakub Kruzik   ierr = PetscUseMethod(pc,"PCDeflationGetType_C",(PC,PCDeflationType*),(pc,type));CHKERRQ(ierr);
312b27c8092SJakub Kruzik   PetscFunctionReturn(0);
313b27c8092SJakub Kruzik }
314b27c8092SJakub Kruzik 
31537eeb815SJakub Kruzik static PetscErrorCode PCSetUp_Deflation(PC pc)
31637eeb815SJakub Kruzik {
317409a357bSJakub Kruzik   PC_Deflation     *def = (PC_Deflation*)pc->data;
318409a357bSJakub Kruzik   KSP              innerksp;
319409a357bSJakub Kruzik   PC               pcinner;
320409a357bSJakub Kruzik   Mat              Amat,nextDef=NULL,*mats;
321409a357bSJakub Kruzik   PetscInt         i,m,red,size,commsize;
322409a357bSJakub Kruzik   PetscBool        match,flgspd,transp=PETSC_FALSE;
323409a357bSJakub Kruzik   MatCompositeType ctype;
324409a357bSJakub Kruzik   MPI_Comm         comm;
325409a357bSJakub Kruzik   const char       *prefix;
32637eeb815SJakub Kruzik   PetscErrorCode   ierr;
32737eeb815SJakub Kruzik 
32837eeb815SJakub Kruzik   PetscFunctionBegin;
329409a357bSJakub Kruzik   if (pc->setupcalled) PetscFunctionReturn(0);
330409a357bSJakub Kruzik   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
331f8bf229cSJakub Kruzik   ierr = PCGetOperators(pc,NULL,&Amat);CHKERRQ(ierr);
332f8bf229cSJakub Kruzik 
333f8bf229cSJakub Kruzik   /* compute a deflation space */
334409a357bSJakub Kruzik   if (def->W || def->Wt) {
335409a357bSJakub Kruzik     def->spacetype = PC_DEFLATION_SPACE_USER;
336409a357bSJakub Kruzik   } else {
337e53e0a0dSJakub Kruzik     ierr = PCDeflationComputeSpace(pc);CHKERRQ(ierr);
33837eeb815SJakub Kruzik   }
33937eeb815SJakub Kruzik 
340409a357bSJakub Kruzik   /* nested deflation */
341409a357bSJakub Kruzik   if (def->W) {
342409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->W,MATCOMPOSITE,&match);CHKERRQ(ierr);
343409a357bSJakub Kruzik     if (match) {
344409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->W,&ctype);CHKERRQ(ierr);
345409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->W,&size);CHKERRQ(ierr);
34637eeb815SJakub Kruzik     }
347409a357bSJakub Kruzik   } else {
348409a357bSJakub Kruzik     ierr = MatCreateTranspose(def->Wt,&def->W);CHKERRQ(ierr);
349409a357bSJakub Kruzik     ierr = PetscObjectTypeCompare((PetscObject)def->Wt,MATCOMPOSITE,&match);CHKERRQ(ierr);
350409a357bSJakub Kruzik     if (match) {
351409a357bSJakub Kruzik       ierr = MatCompositeGetType(def->Wt,&ctype);CHKERRQ(ierr);
352409a357bSJakub Kruzik       ierr = MatCompositeGetNumberMat(def->Wt,&size);CHKERRQ(ierr);
353409a357bSJakub Kruzik     }
354409a357bSJakub Kruzik     transp = PETSC_TRUE;
355409a357bSJakub Kruzik   }
356409a357bSJakub Kruzik   if (match && ctype == MAT_COMPOSITE_MULTIPLICATIVE) {
357409a357bSJakub Kruzik     if (!transp) {
35828da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
35928da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
360409a357bSJakub Kruzik         for (i=0; i<size; i++) {
361409a357bSJakub Kruzik           ierr = MatCompositeGetMat(def->W,i,&mats[i]);CHKERRQ(ierr);
362409a357bSJakub Kruzik         }
363409a357bSJakub Kruzik         size -= 1;
364409a357bSJakub Kruzik         ierr = MatDestroy(&def->W);CHKERRQ(ierr);
365409a357bSJakub Kruzik         def->W = mats[size];
366409a357bSJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[size]);CHKERRQ(ierr);
367409a357bSJakub Kruzik         if (size > 1) {
368409a357bSJakub Kruzik           ierr = MatCreateComposite(comm,size,mats,&nextDef);CHKERRQ(ierr);
369409a357bSJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
370409a357bSJakub Kruzik         } else {
371409a357bSJakub Kruzik           nextDef = mats[0];
372409a357bSJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
373409a357bSJakub Kruzik         }
37428da0170SJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
375409a357bSJakub Kruzik       } else {
376409a357bSJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
377409a357bSJakub Kruzik         ierr = MatCompositeMerge(def->W);CHKERRQ(ierr);
378409a357bSJakub Kruzik       }
37928da0170SJakub Kruzik     } else {
38028da0170SJakub Kruzik       if (def->nestedlvl < def->maxnestedlvl) {
38128da0170SJakub Kruzik         ierr = PetscMalloc1(size,&mats);CHKERRQ(ierr);
38228da0170SJakub Kruzik         for (i=0; i<size; i++) {
38328da0170SJakub Kruzik           ierr = MatCompositeGetMat(def->Wt,i,&mats[i]);CHKERRQ(ierr);
38428da0170SJakub Kruzik         }
38528da0170SJakub Kruzik         size -= 1;
38628da0170SJakub Kruzik         ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
38728da0170SJakub Kruzik         def->Wt = mats[0];
38828da0170SJakub Kruzik         ierr = PetscObjectReference((PetscObject)mats[0]);CHKERRQ(ierr);
38928da0170SJakub Kruzik         if (size > 1) {
39028da0170SJakub Kruzik           ierr = MatCreateComposite(comm,size,&mats[1],&nextDef);CHKERRQ(ierr);
39128da0170SJakub Kruzik           ierr = MatCompositeSetType(nextDef,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
39228da0170SJakub Kruzik         } else {
39328da0170SJakub Kruzik           nextDef = mats[1];
39428da0170SJakub Kruzik           ierr = PetscObjectReference((PetscObject)mats[1]);CHKERRQ(ierr);
395409a357bSJakub Kruzik         }
396409a357bSJakub Kruzik         ierr = PetscFree(mats);CHKERRQ(ierr);
39728da0170SJakub Kruzik       } else {
39828da0170SJakub Kruzik         /* ierr = MatCompositeSetMergeType(def->W,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr); */
39928da0170SJakub Kruzik         ierr = MatCompositeMerge(def->Wt);CHKERRQ(ierr);
40028da0170SJakub Kruzik       }
40128da0170SJakub Kruzik     }
40228da0170SJakub Kruzik   }
40328da0170SJakub Kruzik 
40428da0170SJakub Kruzik   if (transp) {
40528da0170SJakub Kruzik     ierr = MatDestroy(&def->W);CHKERRQ(ierr);
40628da0170SJakub Kruzik     ierr = MatTranspose(def->Wt,MAT_INITIAL_MATRIX,&def->W);CHKERRQ(ierr);
407409a357bSJakub Kruzik   }
408409a357bSJakub Kruzik 
409409a357bSJakub Kruzik   /* setup coarse problem */
410409a357bSJakub Kruzik   if (!def->WtAWinv) {
41128da0170SJakub Kruzik     ierr = MatGetSize(def->W,NULL,&m);CHKERRQ(ierr);
412409a357bSJakub Kruzik     if (!def->WtAW) {
413409a357bSJakub Kruzik       /* TODO add implicit product version ? */
414409a357bSJakub Kruzik       if (!def->AW) {
415409a357bSJakub Kruzik         ierr = MatPtAP(Amat,def->W,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
416409a357bSJakub Kruzik       } else {
417409a357bSJakub Kruzik         ierr = MatTransposeMatMult(def->W,def->AW,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&def->WtAW);CHKERRQ(ierr);
418409a357bSJakub Kruzik       }
419409a357bSJakub Kruzik       /* TODO create MatInheritOption(Mat,MatOption) */
420409a357bSJakub Kruzik       ierr = MatGetOption(Amat,MAT_SPD,&flgspd);CHKERRQ(ierr);
421409a357bSJakub Kruzik       ierr = MatSetOption(def->WtAW,MAT_SPD,flgspd);CHKERRQ(ierr);
422409a357bSJakub Kruzik #if defined(PETSC_USE_DEBUG)
423409a357bSJakub Kruzik       /* Check WtAW is not sigular */
424409a357bSJakub Kruzik       PetscReal *norms;
425409a357bSJakub Kruzik       ierr = PetscMalloc1(m,&norms);CHKERRQ(ierr);
426409a357bSJakub Kruzik       ierr = MatGetColumnNorms(def->WtAW,NORM_INFINITY,norms);CHKERRQ(ierr);
427409a357bSJakub Kruzik       for (i=0; i<m; i++) {
428409a357bSJakub Kruzik         if (norms[i] < 100*PETSC_MACHINE_EPSILON) {
429409a357bSJakub Kruzik           SETERRQ1(comm,PETSC_ERR_SUP,"Column %D of W is in kernel of A.",i);
430409a357bSJakub Kruzik         }
431409a357bSJakub Kruzik       }
432409a357bSJakub Kruzik       ierr = PetscFree(norms);CHKERRQ(ierr);
433409a357bSJakub Kruzik #endif
434409a357bSJakub Kruzik     } else {
435409a357bSJakub Kruzik       ierr = MatGetOption(def->WtAW,MAT_SPD,&flgspd);CHKERRQ(ierr);
436409a357bSJakub Kruzik     }
437409a357bSJakub Kruzik     /* TODO use MATINV */
438ac66f006SJakub Kruzik     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
439409a357bSJakub Kruzik     ierr = KSPCreate(comm,&def->WtAWinv);CHKERRQ(ierr);
440409a357bSJakub Kruzik     ierr = KSPSetOperators(def->WtAWinv,def->WtAW,def->WtAW);CHKERRQ(ierr);
441409a357bSJakub Kruzik     ierr = KSPGetPC(def->WtAWinv,&pcinner);CHKERRQ(ierr);
442557a2f7dSJakub Kruzik     /* Setup KSP and PC */
443557a2f7dSJakub Kruzik     if (nextDef) { /* next level for multilevel deflation */
444557a2f7dSJakub Kruzik       innerksp = def->WtAWinv;
445557a2f7dSJakub Kruzik       /* set default KSPtype */
446557a2f7dSJakub Kruzik       if (!def->ksptype) {
447557a2f7dSJakub Kruzik         def->ksptype = KSPFGMRES;
448557a2f7dSJakub Kruzik         if (flgspd) { /* SPD system */
449557a2f7dSJakub Kruzik           def->ksptype = KSPFCG;
450557a2f7dSJakub Kruzik         }
451557a2f7dSJakub Kruzik       }
452557a2f7dSJakub Kruzik       ierr = KSPSetType(innerksp,def->ksptype);CHKERRQ(ierr); /* TODO iherit from KSP */
453557a2f7dSJakub Kruzik       ierr = PCSetType(pcinner,PCDEFLATION);CHKERRQ(ierr); /* TODO create coarse preconditinoner M_c = WtMW ? */
454557a2f7dSJakub Kruzik       ierr = PCDeflationSetSpace(pcinner,nextDef,transp);CHKERRQ(ierr);
455557a2f7dSJakub Kruzik       ierr = PCDeflationSetLvl_Deflation(pcinner,def->nestedlvl+1,def->maxnestedlvl);CHKERRQ(ierr);
456557a2f7dSJakub Kruzik       /* inherit options TODO if not set */
457557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->ksptype = def->ksptype;
458557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->correct = def->correct;
459557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->adaptiveconv = def->adaptiveconv;
460557a2f7dSJakub Kruzik       ((PC_Deflation*)(pcinner->data))->adaptiveconst = def->adaptiveconst;
461557a2f7dSJakub Kruzik       ierr = MatDestroy(&nextDef);CHKERRQ(ierr);
462557a2f7dSJakub Kruzik     } else { /* the last level */
463557a2f7dSJakub Kruzik       ierr = KSPSetType(def->WtAWinv,KSPPREONLY);CHKERRQ(ierr);
464409a357bSJakub Kruzik       ierr = PCSetType(pcinner,PCTELESCOPE);CHKERRQ(ierr);
465ac66f006SJakub Kruzik       /* ugly hack to not have overwritten PCTELESCOPE */
466ac66f006SJakub Kruzik       if (prefix) {
467ac66f006SJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,prefix);CHKERRQ(ierr);
468ac66f006SJakub Kruzik         ierr = KSPAppendOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
469ac66f006SJakub Kruzik       } else {
470ac66f006SJakub Kruzik         ierr = KSPSetOptionsPrefix(def->WtAWinv,"tel_");CHKERRQ(ierr);
471ac66f006SJakub Kruzik       }
472ac66f006SJakub Kruzik       ierr = PCSetFromOptions(pcinner);CHKERRQ(ierr);
473409a357bSJakub Kruzik       /* Reduction factor choice */
474409a357bSJakub Kruzik       red = def->reductionfact;
475409a357bSJakub Kruzik       if (red < 0) {
476409a357bSJakub Kruzik         ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr);
477409a357bSJakub Kruzik         red  = ceil((float)commsize/ceil((float)m/commsize));
478409a357bSJakub Kruzik         ierr = PetscObjectTypeCompareAny((PetscObject)(def->WtAW),&match,MATSEQDENSE,MATMPIDENSE,MATDENSE,"");CHKERRQ(ierr);
479409a357bSJakub Kruzik         if (match) red = commsize;
480409a357bSJakub Kruzik         ierr = PetscInfo1(pc,"Auto choosing reduction factor %D\n",red);CHKERRQ(ierr); /* TODO add level? */
481409a357bSJakub Kruzik       }
482409a357bSJakub Kruzik       ierr = PCTelescopeSetReductionFactor(pcinner,red);CHKERRQ(ierr);
483ac66f006SJakub Kruzik       ierr = PCSetUp(pcinner);CHKERRQ(ierr);
484409a357bSJakub Kruzik       ierr = PCTelescopeGetKSP(pcinner,&innerksp);CHKERRQ(ierr);
485ac66f006SJakub Kruzik       if (innerksp) {
486409a357bSJakub Kruzik         ierr = KSPGetPC(innerksp,&pcinner);CHKERRQ(ierr);
487409a357bSJakub Kruzik         /* TODO Cholesky if flgspd? */
488ac66f006SJakub Kruzik         ierr = PCSetType(pcinner,PCLU);CHKERRQ(ierr);
489409a357bSJakub Kruzik         //TODO remove explicit matSolverPackage
490409a357bSJakub Kruzik         if (commsize == red) {
491ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU);CHKERRQ(ierr);
492409a357bSJakub Kruzik         } else {
493ac66f006SJakub Kruzik           ierr = PCFactorSetMatSolverType(pcinner,MATSOLVERSUPERLU_DIST);CHKERRQ(ierr);
494409a357bSJakub Kruzik         }
495409a357bSJakub Kruzik       }
496557a2f7dSJakub Kruzik     }
497557a2f7dSJakub Kruzik 
498557a2f7dSJakub Kruzik     if (innerksp) {
499409a357bSJakub Kruzik       /* TODO use def_[lvl]_ if lvl > 0? */
500409a357bSJakub Kruzik       if (prefix) {
501409a357bSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,prefix);CHKERRQ(ierr);
502409a357bSJakub Kruzik         ierr = KSPAppendOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
503409a357bSJakub Kruzik       } else {
504409a357bSJakub Kruzik         ierr = KSPSetOptionsPrefix(innerksp,"def_");CHKERRQ(ierr);
505409a357bSJakub Kruzik       }
506557a2f7dSJakub Kruzik       ierr = KSPSetFromOptions(innerksp);CHKERRQ(ierr);
507557a2f7dSJakub Kruzik       ierr = KSPSetUp(innerksp);CHKERRQ(ierr);
508ac66f006SJakub Kruzik     }
509409a357bSJakub Kruzik     /* TODO: check if WtAWinv is KSP and move following from this if */
510409a357bSJakub Kruzik     //if (def->adaptiveconv) {
511409a357bSJakub Kruzik     //  PetscReal *rnorm;
512409a357bSJakub Kruzik     //  PetscNew(&rnorm);
513409a357bSJakub Kruzik     //  ierr = KSPSetConvergenceTest(def->WtAWinv,KSPDCGConvergedAdaptive_DCG,rnorm,NULL);CHKERRQ(ierr);
514409a357bSJakub Kruzik     //}
515409a357bSJakub Kruzik   }
516557a2f7dSJakub Kruzik   ierr = KSPSetFromOptions(def->WtAWinv);CHKERRQ(ierr);
517557a2f7dSJakub Kruzik   ierr = KSPSetUp(def->WtAWinv);CHKERRQ(ierr);
518409a357bSJakub Kruzik 
519f8bf229cSJakub Kruzik   /* create work vecs */
520f8bf229cSJakub Kruzik   ierr = MatCreateVecs(Amat,NULL,&def->work);CHKERRQ(ierr);
521f8bf229cSJakub Kruzik   ierr = KSPCreateVecs(def->WtAWinv,2,&def->workcoarse,0,NULL);CHKERRQ(ierr);
52237eeb815SJakub Kruzik   PetscFunctionReturn(0);
52337eeb815SJakub Kruzik }
524b30d4775SJakub Kruzik 
52537eeb815SJakub Kruzik static PetscErrorCode PCReset_Deflation(PC pc)
52637eeb815SJakub Kruzik {
527b30d4775SJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
52837eeb815SJakub Kruzik   PetscErrorCode    ierr;
52937eeb815SJakub Kruzik 
53037eeb815SJakub Kruzik   PetscFunctionBegin;
531b30d4775SJakub Kruzik   ierr = VecDestroy(&def->work);CHKERRQ(ierr);
532b30d4775SJakub Kruzik   ierr = VecDestroyVecs(2,&def->workcoarse);CHKERRQ(ierr);
533b30d4775SJakub Kruzik   ierr = MatDestroy(&def->W);CHKERRQ(ierr);
534b30d4775SJakub Kruzik   ierr = MatDestroy(&def->Wt);CHKERRQ(ierr);
535b30d4775SJakub Kruzik   ierr = MatDestroy(&def->AW);CHKERRQ(ierr);
536b30d4775SJakub Kruzik   ierr = MatDestroy(&def->WtAW);CHKERRQ(ierr);
537b30d4775SJakub Kruzik   ierr = KSPDestroy(&def->WtAWinv);CHKERRQ(ierr);
53837eeb815SJakub Kruzik   PetscFunctionReturn(0);
53937eeb815SJakub Kruzik }
54037eeb815SJakub Kruzik 
54137eeb815SJakub Kruzik /*
54237eeb815SJakub Kruzik    PCDestroy_Deflation - Destroys the private context for the Deflation preconditioner
54337eeb815SJakub Kruzik    that was created with PCCreate_Deflation().
54437eeb815SJakub Kruzik 
54537eeb815SJakub Kruzik    Input Parameter:
54637eeb815SJakub Kruzik .  pc - the preconditioner context
54737eeb815SJakub Kruzik 
54837eeb815SJakub Kruzik    Application Interface Routine: PCDestroy()
54937eeb815SJakub Kruzik */
55037eeb815SJakub Kruzik static PetscErrorCode PCDestroy_Deflation(PC pc)
55137eeb815SJakub Kruzik {
55237eeb815SJakub Kruzik   PetscErrorCode ierr;
55337eeb815SJakub Kruzik 
55437eeb815SJakub Kruzik   PetscFunctionBegin;
55537eeb815SJakub Kruzik   ierr = PCReset_Deflation(pc);CHKERRQ(ierr);
55637eeb815SJakub Kruzik   ierr = PetscFree(pc->data);CHKERRQ(ierr);
55737eeb815SJakub Kruzik   PetscFunctionReturn(0);
55837eeb815SJakub Kruzik }
55937eeb815SJakub Kruzik 
5608513960bSJakub Kruzik static PetscErrorCode PCView_Deflation(PC pc,PetscViewer viewer)
5618513960bSJakub Kruzik {
5628513960bSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
5638513960bSJakub Kruzik   PetscBool         iascii;
5648513960bSJakub Kruzik   PetscErrorCode    ierr;
5658513960bSJakub Kruzik 
5668513960bSJakub Kruzik   PetscFunctionBegin;
5678513960bSJakub Kruzik   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5688513960bSJakub Kruzik   if (iascii) {
5698513960bSJakub Kruzik     //if (cg->adaptiveconv) {
5708513960bSJakub Kruzik     //  ierr = PetscViewerASCIIPrintf(viewer,"  DCG: using adaptive precision for inner solve with C=%.1e\n",cg->adaptiveconst);CHKERRQ(ierr);
5718513960bSJakub Kruzik     //}
5728513960bSJakub Kruzik     if (def->correct) {
5738513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Using CP correction\n");CHKERRQ(ierr);
5748513960bSJakub Kruzik     }
5758513960bSJakub Kruzik     if (!def->nestedlvl) {
5768513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  Deflation space type: %s\n",PCDeflationSpaceTypes[def->spacetype]);CHKERRQ(ierr);
5778513960bSJakub Kruzik       ierr = PetscViewerASCIIPrintf(viewer,"  DCG %s\n",def->extendsp ? "extended" : "truncated");CHKERRQ(ierr);
5788513960bSJakub Kruzik     }
5798513960bSJakub Kruzik 
5808513960bSJakub Kruzik     ierr = PetscViewerASCIIPrintf(viewer,"--- Coarse solver\n");CHKERRQ(ierr);
5818513960bSJakub Kruzik     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
5828513960bSJakub Kruzik     ierr = KSPView(def->WtAWinv,viewer);CHKERRQ(ierr);
5838513960bSJakub Kruzik     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
5848513960bSJakub Kruzik   }
5858513960bSJakub Kruzik   PetscFunctionReturn(0);
5868513960bSJakub Kruzik }
5878513960bSJakub Kruzik 
58837eeb815SJakub Kruzik static PetscErrorCode PCSetFromOptions_Deflation(PetscOptionItems *PetscOptionsObject,PC pc)
58937eeb815SJakub Kruzik {
590*880d05ceSJakub Kruzik   PC_Deflation      *def = (PC_Deflation*)pc->data;
59137eeb815SJakub Kruzik   PetscErrorCode    ierr;
59237eeb815SJakub Kruzik   PetscBool         flg;
59337eeb815SJakub Kruzik   PCDeflationType   deflt,type;
59437eeb815SJakub Kruzik 
59537eeb815SJakub Kruzik   PetscFunctionBegin;
59637eeb815SJakub Kruzik   ierr = PCDeflationGetType(pc,&deflt);CHKERRQ(ierr);
59737eeb815SJakub Kruzik   ierr = PetscOptionsHead(PetscOptionsObject,"Deflation options");CHKERRQ(ierr);
598*880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_type","Determine type of deflation","PCDeflationSetType",PCDeflationTypes,(PetscEnum)deflt,(PetscEnum*)&type,&flg);CHKERRQ(ierr);
59937eeb815SJakub Kruzik   if (flg) {
60037eeb815SJakub Kruzik     ierr = PCDeflationSetType(pc,type);CHKERRQ(ierr);
60137eeb815SJakub Kruzik   }
602*880d05ceSJakub Kruzik   ierr = PetscOptionsEnum("-pc_deflation_compute_space","Compute deflation space","PCDeflationSetSpace",PCDeflationSpaceTypes,(PetscEnum)def->spacetype,(PetscEnum*)&def->spacetype,NULL);CHKERRQ(ierr);
603*880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_compute_space_size","Set size of the deflation space to compute","PCDeflationSetSpace",def->spacesize,&def->spacesize,NULL);CHKERRQ(ierr);
604*880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_space_extend","Extend deflation space instead of truncating (wavelets)","PCDeflation",def->extendsp,&def->extendsp,NULL);CHKERRQ(ierr);
605*880d05ceSJakub Kruzik //TODO add set function and fix manpages
606*880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_initdef","Use only initialization step - Initdef","PCDeflation",def->init,&def->init,NULL);CHKERRQ(ierr);
607*880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_correct","Add Qr to descent direction","PCDeflation",def->correct,&def->correct,NULL);CHKERRQ(ierr);
608*880d05ceSJakub Kruzik   ierr = PetscOptionsBool("-pc_deflation_adaptive","Adaptive stopping criteria","PCDeflation",def->adaptiveconv,&def->adaptiveconv,NULL);CHKERRQ(ierr);
609*880d05ceSJakub Kruzik   ierr = PetscOptionsReal("-pc_deflation_adaptive_const","Adaptive stopping criteria constant","PCDeflation",def->adaptiveconst,&def->adaptiveconst,NULL);CHKERRQ(ierr);
610*880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_redfact","Reduction factor for coarse problem solution","PCDeflation",def->reductionfact,&def->reductionfact,NULL);CHKERRQ(ierr);
611*880d05ceSJakub Kruzik   ierr = PetscOptionsInt("-pc_deflation_max_nested_lvl","Maximum of nested deflation levels","PCDeflation",def->maxnestedlvl,&def->maxnestedlvl,NULL);CHKERRQ(ierr);
61237eeb815SJakub Kruzik   ierr = PetscOptionsTail();CHKERRQ(ierr);
61337eeb815SJakub Kruzik   PetscFunctionReturn(0);
61437eeb815SJakub Kruzik }
61537eeb815SJakub Kruzik 
61637eeb815SJakub Kruzik /*MC
617e361f8b5SJakub Kruzik      PCDEFLATION - Deflation preconditioner shifts part of the spectrum to zero (deflates)
618e361f8b5SJakub Kruzik      or to a predefined value
61937eeb815SJakub Kruzik 
62037eeb815SJakub Kruzik    Options Database Key:
621e361f8b5SJakub Kruzik +    -pc_deflation_type <init,pre,post> - selects approach to deflation (default: pre)
62237eeb815SJakub Kruzik -    -pc_jacobi_abs - use the absolute value of the diagonal entry
62337eeb815SJakub Kruzik 
62437eeb815SJakub Kruzik    Level: beginner
62537eeb815SJakub Kruzik 
62637eeb815SJakub Kruzik   Notes:
627e361f8b5SJakub Kruzik     todo
62837eeb815SJakub Kruzik 
62937eeb815SJakub Kruzik .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
630e662bc50SJakub Kruzik            PCDeflationSetType(), PCDeflationSetSpace()
63137eeb815SJakub Kruzik M*/
63237eeb815SJakub Kruzik 
63337eeb815SJakub Kruzik PETSC_EXTERN PetscErrorCode PCCreate_Deflation(PC pc)
63437eeb815SJakub Kruzik {
63537eeb815SJakub Kruzik   PC_Deflation   *def;
63637eeb815SJakub Kruzik   PetscErrorCode ierr;
63737eeb815SJakub Kruzik 
63837eeb815SJakub Kruzik   PetscFunctionBegin;
63937eeb815SJakub Kruzik   ierr     = PetscNewLog(pc,&def);CHKERRQ(ierr);
64037eeb815SJakub Kruzik   pc->data = (void*)def;
64137eeb815SJakub Kruzik 
642e662bc50SJakub Kruzik   def->init          = PETSC_FALSE;
643e662bc50SJakub Kruzik   def->pre           = PETSC_TRUE;
644e662bc50SJakub Kruzik   def->correct       = PETSC_FALSE;
645e662bc50SJakub Kruzik   def->truenorm      = PETSC_TRUE;
646e662bc50SJakub Kruzik   def->reductionfact = -1;
647e662bc50SJakub Kruzik   def->spacetype     = PC_DEFLATION_SPACE_HAAR;
648e662bc50SJakub Kruzik   def->spacesize     = 1;
649e662bc50SJakub Kruzik   def->extendsp      = PETSC_FALSE;
650e662bc50SJakub Kruzik   def->nestedlvl     = 0;
651e662bc50SJakub Kruzik   def->maxnestedlvl  = 0;
652e662bc50SJakub Kruzik   def->adaptiveconv  = PETSC_FALSE;
653e662bc50SJakub Kruzik   def->adaptiveconst = 1.0;
65437eeb815SJakub Kruzik 
65537eeb815SJakub Kruzik   /*
65637eeb815SJakub Kruzik       Set the pointers for the functions that are provided above.
65737eeb815SJakub Kruzik       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
65837eeb815SJakub Kruzik       are called, they will automatically call these functions.  Note we
65937eeb815SJakub Kruzik       choose not to provide a couple of these functions since they are
66037eeb815SJakub Kruzik       not needed.
66137eeb815SJakub Kruzik   */
66237eeb815SJakub Kruzik   pc->ops->apply               = PCApply_Deflation;
66337eeb815SJakub Kruzik   pc->ops->applytranspose      = PCApply_Deflation;
664b27c8092SJakub Kruzik   pc->ops->presolve            = PCPreSolve_Deflation;
665b27c8092SJakub Kruzik   pc->ops->postsolve           = 0;
66637eeb815SJakub Kruzik   pc->ops->setup               = PCSetUp_Deflation;
66737eeb815SJakub Kruzik   pc->ops->reset               = PCReset_Deflation;
66837eeb815SJakub Kruzik   pc->ops->destroy             = PCDestroy_Deflation;
66937eeb815SJakub Kruzik   pc->ops->setfromoptions      = PCSetFromOptions_Deflation;
6708513960bSJakub Kruzik   pc->ops->view                = PCView_Deflation;
67137eeb815SJakub Kruzik   pc->ops->applyrichardson     = 0;
67237eeb815SJakub Kruzik   pc->ops->applysymmetricleft  = 0;
67337eeb815SJakub Kruzik   pc->ops->applysymmetricright = 0;
67437eeb815SJakub Kruzik 
67537eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetType_C",PCDeflationSetType_Deflation);CHKERRQ(ierr);
67637eeb815SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationGetType_C",PCDeflationGetType_Deflation);CHKERRQ(ierr);
677e662bc50SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetSpace_C",PCDeflationSetSpace_Deflation);CHKERRQ(ierr);
6785c0c31c5SJakub Kruzik   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCDeflationSetLvl_C",PCDeflationSetLvl_Deflation);CHKERRQ(ierr);
67937eeb815SJakub Kruzik   PetscFunctionReturn(0);
68037eeb815SJakub Kruzik }
68137eeb815SJakub Kruzik 
682