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