xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 7a7894de0ae7ea508ad2b43522ac5a3ba7d8845d)
1 
2 /*
3    Defines a direct factorization preconditioner for any Mat implementation
4    Note: this need not be consided a preconditioner since it supplies
5          a direct solver.
6 */
7 #include "src/ksp/pc/pcimpl.h"                /*I "petscpc.h" I*/
8 
9 typedef struct {
10   Mat             fact;             /* factored matrix */
11   PetscReal       actualfill;       /* actual fill in factor */
12   PetscTruth      inplace;          /* flag indicating in-place factorization */
13   IS              row,col;          /* index sets used for reordering */
14   MatOrderingType ordering;         /* matrix ordering */
15   PetscTruth      reuseordering;    /* reuses previous reordering computed */
16   PetscTruth      reusefill;        /* reuse fill from previous Cholesky */
17   MatFactorInfo   info;
18 } PC_Cholesky;
19 
20 EXTERN_C_BEGIN
21 #undef __FUNCT__
22 #define __FUNCT__ "PCCholeskySetReuseOrdering_Cholesky"
23 PetscErrorCode PCCholeskySetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
24 {
25   PC_Cholesky *lu;
26 
27   PetscFunctionBegin;
28   lu               = (PC_Cholesky*)pc->data;
29   lu->reuseordering = flag;
30   PetscFunctionReturn(0);
31 }
32 EXTERN_C_END
33 
34 EXTERN_C_BEGIN
35 #undef __FUNCT__
36 #define __FUNCT__ "PCCholeskySetReuseFill_Cholesky"
37 PetscErrorCode PCCholeskySetReuseFill_Cholesky(PC pc,PetscTruth flag)
38 {
39   PC_Cholesky *lu;
40 
41   PetscFunctionBegin;
42   lu = (PC_Cholesky*)pc->data;
43   lu->reusefill = flag;
44   PetscFunctionReturn(0);
45 }
46 EXTERN_C_END
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "PCSetFromOptions_Cholesky"
50 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
51 {
52   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
53   PetscErrorCode ierr;
54   PetscTruth     flg;
55   char           tname[256];
56   PetscFList     ordlist;
57 
58   PetscFunctionBegin;
59   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
60   ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr);
61   ierr = PetscOptionsName("-pc_cholesky_in_place","Form Cholesky in the same memory as the matrix","PCCholeskySetUseInPlace",&flg);CHKERRQ(ierr);
62   if (flg) {
63     ierr = PCCholeskySetUseInPlace(pc);CHKERRQ(ierr);
64   }
65   ierr = PetscOptionsReal("-pc_cholesky_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCCholeskySetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
66 
67   ierr = PetscOptionsName("-pc_cholesky_reuse_fill","Use fill from previous factorization","PCCholeskySetReuseFill",&flg);CHKERRQ(ierr);
68   if (flg) {
69     ierr = PCCholeskySetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
70   }
71   ierr = PetscOptionsName("-pc_cholesky_reuse_ordering","Reuse ordering from previous factorization","PCCholeskySetReuseOrdering",&flg);CHKERRQ(ierr);
72   if (flg) {
73     ierr = PCCholeskySetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
74   }
75 
76   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
77   ierr = PetscOptionsList("-pc_cholesky_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCCholeskySetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
78   if (flg) {
79     ierr = PCCholeskySetMatOrdering(pc,tname);CHKERRQ(ierr);
80   }
81   ierr = PetscOptionsName("-pc_cholesky_damping","Damping added to diagonal","PCCholestkySetDamping",&flg);CHKERRQ(ierr);
82   if (flg) {
83     ierr = PCCholeskySetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
84   }
85   ierr = PetscOptionsReal("-pc_cholesky_damping","Damping added to diagonal","PCCholeskySetDamping",lu->info.damping,&lu->info.damping,0);CHKERRQ(ierr);
86   ierr = PetscOptionsName("-pc_cholesky_shift","Manteuffel shift applied to diagonal","PCCholeskySetShift",&flg);CHKERRQ(ierr);
87   if (flg) {
88     ierr = PCCholeskySetShift(pc,PETSC_TRUE);CHKERRQ(ierr);
89   }
90 
91   ierr = PetscOptionsTail();CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PCView_Cholesky"
97 static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
98 {
99   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
100   PetscErrorCode ierr;
101   PetscTruth     iascii,isstring;
102 
103   PetscFunctionBegin;
104   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
105   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
106   if (iascii) {
107     MatInfo info;
108 
109     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");CHKERRQ(ierr);}
110     else             {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");CHKERRQ(ierr);}
111     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
112     if (lu->fact) {
113       ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
114       ierr = PetscViewerASCIIPrintf(viewer,"    Cholesky nonzeros %g\n",info.nz_used);CHKERRQ(ierr);
115       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr);
116       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
117       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
118     }
119     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
120     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
121   } else if (isstring) {
122     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
123   } else {
124     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
125   }
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "PCGetFactoredMatrix_Cholesky"
131 static PetscErrorCode PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
132 {
133   PC_Cholesky *dir = (PC_Cholesky*)pc->data;
134 
135   PetscFunctionBegin;
136   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
137   *mat = dir->fact;
138   PetscFunctionReturn(0);
139 }
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "PCSetUp_Cholesky"
143 static PetscErrorCode PCSetUp_Cholesky(PC pc)
144 {
145   PetscErrorCode ierr;
146   PetscTruth     flg;
147   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
148 
149   PetscFunctionBegin;
150   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
151 
152   if (dir->inplace) {
153     if (dir->row && dir->col && (dir->row != dir->col)) {
154       ierr = ISDestroy(dir->row);CHKERRQ(ierr);
155       dir->row = 0;
156     }
157     if (dir->col) {
158       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
159       dir->col = 0;
160     }
161     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
162     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
163       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
164       dir->col=0;
165     }
166     if (dir->row) {PetscLogObjectParent(pc,dir->row);}
167     ierr = MatCholeskyFactor(pc->pmat,dir->row,&dir->info);CHKERRQ(ierr);
168     dir->fact = pc->pmat;
169   } else {
170     MatInfo info;
171     if (!pc->setupcalled) {
172       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
173       if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
174         ierr = ISDestroy(dir->col);CHKERRQ(ierr);
175         dir->col=0;
176       }
177       ierr = PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
178       if (flg) {
179         PetscReal tol = 1.e-10;
180         ierr = PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
181         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
182       }
183       if (dir->row) {PetscLogObjectParent(pc,dir->row);}
184       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
185       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
186       dir->actualfill = info.fill_ratio_needed;
187       PetscLogObjectParent(pc,dir->fact);
188     } else if (pc->flag != SAME_NONZERO_PATTERN) {
189       if (!dir->reuseordering) {
190         if (dir->row && dir->col && (dir->row != dir->col)) {
191           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
192           dir->row = 0;
193         }
194         if (dir->col) {
195           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
196           dir->col =0;
197         }
198         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
199         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
200           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
201           dir->col=0;
202         }
203         ierr = PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
204         if (flg) {
205           PetscReal tol = 1.e-10;
206           ierr = PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
207           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
208         }
209         if (dir->row) {PetscLogObjectParent(pc,dir->row);}
210       }
211       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
212       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
213       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
214       dir->actualfill = info.fill_ratio_needed;
215       PetscLogObjectParent(pc,dir->fact);
216     }
217     ierr = MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
218   }
219   PetscFunctionReturn(0);
220 }
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "PCDestroy_Cholesky"
224 static PetscErrorCode PCDestroy_Cholesky(PC pc)
225 {
226   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
231   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
232   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
233   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
234   ierr = PetscFree(dir);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "PCApply_Cholesky"
240 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
241 {
242   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
247   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
248   PetscFunctionReturn(0);
249 }
250 
251 #undef __FUNCT__
252 #define __FUNCT__ "PCApplyTranspose_Cholesky"
253 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
254 {
255   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
260   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
261   PetscFunctionReturn(0);
262 }
263 
264 /* -----------------------------------------------------------------------------------*/
265 
266 EXTERN_C_BEGIN
267 #undef __FUNCT__
268 #define __FUNCT__ "PCCholeskySetFill_Cholesky"
269 PetscErrorCode PCCholeskySetFill_Cholesky(PC pc,PetscReal fill)
270 {
271   PC_Cholesky *dir;
272 
273   PetscFunctionBegin;
274   dir = (PC_Cholesky*)pc->data;
275   dir->info.fill = fill;
276   PetscFunctionReturn(0);
277 }
278 EXTERN_C_END
279 
280 EXTERN_C_BEGIN
281 #undef __FUNCT__
282 #define __FUNCT__ "PCCholeskySetDamping_Cholesky"
283 PetscErrorCode PCCholeskySetDamping_Cholesky(PC pc,PetscReal damping)
284 {
285   PC_Cholesky *dir;
286 
287   PetscFunctionBegin;
288   dir = (PC_Cholesky*)pc->data;
289   if (damping == (PetscReal) PETSC_DECIDE) {
290     dir->info.damping = 1.e-12;
291   } else {
292     dir->info.damping = damping;
293   }
294   PetscFunctionReturn(0);
295 }
296 EXTERN_C_END
297 
298 EXTERN_C_BEGIN
299 #undef __FUNCT__
300 #define __FUNCT__ "PCCholeskySetShift_Cholesky"
301 PetscErrorCode PCCholeskySetShift_Cholesky(PC pc,PetscTruth shift)
302 {
303   PC_Cholesky *dir;
304 
305   PetscFunctionBegin;
306   dir = (PC_Cholesky*)pc->data;
307   dir->info.shift = shift;
308   if (shift) dir->info.shift_fraction = 0.0;
309   PetscFunctionReturn(0);
310 }
311 EXTERN_C_END
312 
313 EXTERN_C_BEGIN
314 #undef __FUNCT__
315 #define __FUNCT__ "PCCholeskySetUseInPlace_Cholesky"
316 PetscErrorCode PCCholeskySetUseInPlace_Cholesky(PC pc)
317 {
318   PC_Cholesky *dir;
319 
320   PetscFunctionBegin;
321   dir = (PC_Cholesky*)pc->data;
322   dir->inplace = PETSC_TRUE;
323   PetscFunctionReturn(0);
324 }
325 EXTERN_C_END
326 
327 EXTERN_C_BEGIN
328 #undef __FUNCT__
329 #define __FUNCT__ "PCCholeskySetMatOrdering_Cholesky"
330 PetscErrorCode PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
331 {
332   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
333   PetscErrorCode ierr;
334 
335   PetscFunctionBegin;
336   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
337   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
338   PetscFunctionReturn(0);
339 }
340 EXTERN_C_END
341 
342 /* -----------------------------------------------------------------------------------*/
343 
344 #undef __FUNCT__
345 #define __FUNCT__ "PCCholeskySetReuseOrdering"
346 /*@
347    PCCholeskySetReuseOrdering - When similar matrices are factored, this
348    causes the ordering computed in the first factor to be used for all
349    following factors.
350 
351    Collective on PC
352 
353    Input Parameters:
354 +  pc - the preconditioner context
355 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
356 
357    Options Database Key:
358 .  -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
359 
360    Level: intermediate
361 
362 .keywords: PC, levels, reordering, factorization, incomplete, LU
363 
364 .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
365 @*/
366 PetscErrorCode PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
367 {
368   PetscErrorCode ierr,(*f)(PC,PetscTruth);
369 
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
372   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
373   if (f) {
374     ierr = (*f)(pc,flag);CHKERRQ(ierr);
375   }
376   PetscFunctionReturn(0);
377 }
378 
379 #undef __FUNCT__
380 #define __FUNCT__ "PCCholeskySetReuseFill"
381 /*@
382    PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
383    this causes later ones to use the fill computed in the initial factorization.
384 
385    Collective on PC
386 
387    Input Parameters:
388 +  pc - the preconditioner context
389 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
390 
391    Options Database Key:
392 .  -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
393 
394    Level: intermediate
395 
396 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
397 
398 .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
399 @*/
400 PetscErrorCode PCCholeskySetReuseFill(PC pc,PetscTruth flag)
401 {
402   PetscErrorCode ierr,(*f)(PC,PetscTruth);
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
406   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
407   if (f) {
408     ierr = (*f)(pc,flag);CHKERRQ(ierr);
409   }
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "PCCholeskySetFill"
415 /*@
416    PCCholeskySetFill - Indicates the amount of fill you expect in the factored matrix,
417    fill = number nonzeros in factor/number nonzeros in original matrix.
418 
419    Collective on PC
420 
421    Input Parameters:
422 +  pc - the preconditioner context
423 -  fill - amount of expected fill
424 
425    Options Database Key:
426 .  -pc_cholesky_fill <fill> - Sets fill amount
427 
428    Level: intermediate
429 
430    Note:
431    For sparse matrix factorizations it is difficult to predict how much
432    fill to expect. By running with the option -log_info PETSc will print the
433    actual amount of fill used; allowing you to set the value accurately for
434    future runs. Default PETSc uses a value of 5.0
435 
436 .keywords: PC, set, factorization, direct, fill
437 
438 .seealso: PCILUSetFill()
439 @*/
440 PetscErrorCode PCCholeskySetFill(PC pc,PetscReal fill)
441 {
442   PetscErrorCode ierr,(*f)(PC,PetscReal);
443 
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
446   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
447   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
448   if (f) {
449     ierr = (*f)(pc,fill);CHKERRQ(ierr);
450   }
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "PCCholeskySetDamping"
456 /*@
457    PCCholeskySetDamping - Adds this quantity to the diagonal of the matrix during the
458    Cholesky numerical factorization.
459 
460    Collective on PC
461 
462    Input Parameters:
463 +  pc - the preconditioner context
464 -  damping - amount of damping
465 
466    Options Database Key:
467 .  -pc_cholesky_damping <damping> - Sets damping amount
468 
469    Level: intermediate
470 
471 .keywords: PC, set, factorization, direct, fill
472 
473 .seealso: PCCholeskySetFill(), PCILUSetDamping()
474 @*/
475 PetscErrorCode PCCholeskySetDamping(PC pc,PetscReal damping)
476 {
477   PetscErrorCode ierr,(*f)(PC,PetscReal);
478 
479   PetscFunctionBegin;
480   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
481   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
482   if (f) {
483     ierr = (*f)(pc,damping);CHKERRQ(ierr);
484   }
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "PCCholeskySetShift"
490 /*@
491    PCCholeskySetShift - specify whether to use Manteuffel shifting of Cholesky.
492    If an Cholesky factorisation breaks down because of nonpositive pivots,
493    adding sufficient identity to the diagonal will remedy this.
494    Setting this causes a bisection method to find the minimum shift that
495    will lead to a well-defined Cholesky.
496 
497    Input parameters:
498 +  pc - the preconditioner context
499 -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
500 
501    Options Database Key:
502 .  -pc_ilu_shift - Activate PCCholeskySetShift()
503 
504    Level: intermediate
505 
506 .keywords: PC, indefinite, factorization, incomplete, Cholesky
507 
508 .seealso: PCILUSetShift()
509 @*/
510 PetscErrorCode PCCholeskySetShift(PC pc,PetscTruth shift)
511 {
512   PetscErrorCode ierr,(*f)(PC,PetscTruth);
513 
514   PetscFunctionBegin;
515   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
516   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetShift_C",(void (**)(void))&f);CHKERRQ(ierr);
517   if (f) {
518     ierr = (*f)(pc,shift);CHKERRQ(ierr);
519   }
520   PetscFunctionReturn(0);
521 }
522 
523 #undef __FUNCT__
524 #define __FUNCT__ "PCCholeskySetUseInPlace"
525 /*@
526    PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
527    For dense matrices, this enables the solution of much larger problems.
528    For sparse matrices the factorization cannot be done truly in-place
529    so this does not save memory during the factorization, but after the matrix
530    is factored, the original unfactored matrix is freed, thus recovering that
531    space.
532 
533    Collective on PC
534 
535    Input Parameters:
536 .  pc - the preconditioner context
537 
538    Options Database Key:
539 .  -pc_cholesky_in_place - Activates in-place factorization
540 
541    Notes:
542    PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
543    a different matrix is provided for the multiply and the preconditioner in
544    a call to KSPSetOperators().
545    This is because the Krylov space methods require an application of the
546    matrix multiplication, which is not possible here because the matrix has
547    been factored in-place, replacing the original matrix.
548 
549    Level: intermediate
550 
551 .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
552 
553 .seealso: PCICholeskySetUseInPlace()
554 @*/
555 PetscErrorCode PCCholeskySetUseInPlace(PC pc)
556 {
557   PetscErrorCode ierr,(*f)(PC);
558 
559   PetscFunctionBegin;
560   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
561   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
562   if (f) {
563     ierr = (*f)(pc);CHKERRQ(ierr);
564   }
565   PetscFunctionReturn(0);
566 }
567 
568 #undef __FUNCT__
569 #define __FUNCT__ "PCCholeskySetMatOrdering"
570 /*@
571     PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
572     be used it the Cholesky factorization.
573 
574     Collective on PC
575 
576     Input Parameters:
577 +   pc - the preconditioner context
578 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
579 
580     Options Database Key:
581 .   -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
582 
583     Level: intermediate
584 
585 .seealso: PCICholeskySetMatOrdering()
586 @*/
587 PetscErrorCode PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
588 {
589   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
590 
591   PetscFunctionBegin;
592   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
593   if (f) {
594     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
595   }
596   PetscFunctionReturn(0);
597 }
598 
599 /*MC
600    PCCholesky - Uses a direct solver, based on Cholesky factorization, as a preconditioner
601 
602    Options Database Keys:
603 +  -pc_cholesky_reuse_ordering - Activate PCLUSetReuseOrdering()
604 .  -pc_cholesky_reuse_fill - Activates PCLUSetReuseFill()
605 .  -pc_cholesky_fill <fill> - Sets fill amount
606 .  -pc_cholesky_damping <damping> - Sets damping amount
607 .  -pc_cholesky_shift - Activates Manteuffel shift
608 .  -pc_cholesky_in_place - Activates in-place factorization
609 -  -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
610 
611    Notes: Not all options work for all matrix formats
612 
613    Level: beginner
614 
615    Concepts: Cholesky factorization, direct solver
616 
617    Notes: Usually this will compute an "exact" solution in one iteration and does
618           not need a Krylov method (i.e. you can use -ksp_type preonly, or
619           KSPSetType(ksp,KSPPREONLY) for the Krylov method
620 
621 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
622            PCILU, PCLU, PCICC, PCCholeskySetReuseOrdering(), PCCholeskySetReuseFill(), PCGetFactoredMatrix(),
623            PCCholeskySetFill(), PCCholeskySetDamping(), PCCholeskySetShift(),
624 	   PCCholeskySetUseInPlace(), PCCholeskySetMatOrdering()
625 
626 M*/
627 
628 EXTERN_C_BEGIN
629 #undef __FUNCT__
630 #define __FUNCT__ "PCCreate_Cholesky"
631 PetscErrorCode PCCreate_Cholesky(PC pc)
632 {
633   PetscErrorCode ierr;
634   PC_Cholesky    *dir;
635 
636   PetscFunctionBegin;
637   ierr = PetscNew(PC_Cholesky,&dir);CHKERRQ(ierr);
638   PetscLogObjectMemory(pc,sizeof(PC_Cholesky));
639 
640   dir->fact                   = 0;
641   dir->inplace                = PETSC_FALSE;
642   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
643   dir->info.fill              = 5.0;
644   dir->info.damping           = 0.0;
645   dir->info.shift             = PETSC_FALSE;
646   dir->info.shift_fraction    = 0.0;
647   dir->info.pivotinblocks     = 1.0;
648   dir->col                    = 0;
649   dir->row                    = 0;
650   ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
651   dir->reusefill        = PETSC_FALSE;
652   dir->reuseordering    = PETSC_FALSE;
653   pc->data              = (void*)dir;
654 
655   pc->ops->destroy           = PCDestroy_Cholesky;
656   pc->ops->apply             = PCApply_Cholesky;
657   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
658   pc->ops->setup             = PCSetUp_Cholesky;
659   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
660   pc->ops->view              = PCView_Cholesky;
661   pc->ops->applyrichardson   = 0;
662   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;
663 
664   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetFill_C","PCCholeskySetFill_Cholesky",
665                     PCCholeskySetFill_Cholesky);CHKERRQ(ierr);
666   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetDamping_C","PCCholeskySetDamping_Cholesky",
667                     PCCholeskySetDamping_Cholesky);CHKERRQ(ierr);
668   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetShift_C","PCCholeskySetShift_Cholesky",
669                     PCCholeskySetShift_Cholesky);CHKERRQ(ierr);
670   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetUseInPlace_C","PCCholeskySetUseInPlace_Cholesky",
671                     PCCholeskySetUseInPlace_Cholesky);CHKERRQ(ierr);
672   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetMatOrdering_C","PCCholeskySetMatOrdering_Cholesky",
673                     PCCholeskySetMatOrdering_Cholesky);CHKERRQ(ierr);
674   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseOrdering_C","PCCholeskySetReuseOrdering_Cholesky",
675                     PCCholeskySetReuseOrdering_Cholesky);CHKERRQ(ierr);
676   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseFill_C","PCCholeskySetReuseFill_Cholesky",
677                     PCCholeskySetReuseFill_Cholesky);CHKERRQ(ierr);
678   PetscFunctionReturn(0);
679 }
680 EXTERN_C_END
681