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