xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision ce0a2cd1da0658c2b28aad1be2e2c8e41567bece)
1 #define PETSCKSP_DLL
2 
3 /*
4    Defines a direct factorization preconditioner for any Mat implementation
5    Note: this need not be consided a preconditioner since it supplies
6          a direct solver.
7 */
8 #include "private/pcimpl.h"                /*I "petscpc.h" I*/
9 
10 typedef struct {
11   Mat             fact;             /* factored matrix */
12   PetscReal       actualfill;       /* actual fill in factor */
13   PetscTruth      inplace;          /* flag indicating in-place factorization */
14   IS              row,col;          /* index sets used for reordering */
15   MatOrderingType ordering;         /* matrix ordering */
16   PetscTruth      reuseordering;    /* reuses previous reordering computed */
17   PetscTruth      reusefill;        /* reuse fill from previous Cholesky */
18   MatFactorInfo   info;
19 } PC_Cholesky;
20 
21 EXTERN_C_BEGIN
22 #undef __FUNCT__
23 #define __FUNCT__ "PCFactorSetZeroPivot_Cholesky"
24 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Cholesky(PC pc,PetscReal z)
25 {
26   PC_Cholesky *ch;
27 
28   PetscFunctionBegin;
29   ch                 = (PC_Cholesky*)pc->data;
30   ch->info.zeropivot = z;
31   PetscFunctionReturn(0);
32 }
33 EXTERN_C_END
34 
35 EXTERN_C_BEGIN
36 #undef __FUNCT__
37 #define __FUNCT__ "PCFactorSetShiftNonzero_Cholesky"
38 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Cholesky(PC pc,PetscReal shift)
39 {
40   PC_Cholesky *dir;
41 
42   PetscFunctionBegin;
43   dir = (PC_Cholesky*)pc->data;
44   if (shift == (PetscReal) PETSC_DECIDE) {
45     dir->info.shiftnz = 1.e-12;
46   } else {
47     dir->info.shiftnz = shift;
48   }
49   PetscFunctionReturn(0);
50 }
51 EXTERN_C_END
52 
53 EXTERN_C_BEGIN
54 #undef __FUNCT__
55 #define __FUNCT__ "PCFactorSetShiftPd_Cholesky"
56 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Cholesky(PC pc,PetscTruth shift)
57 {
58   PC_Cholesky *dir;
59 
60   PetscFunctionBegin;
61   dir = (PC_Cholesky*)pc->data;
62   if (shift) {
63     dir->info.shift_fraction = 0.0;
64     dir->info.shiftpd = 1.0;
65   } else {
66     dir->info.shiftpd = 0.0;
67   }
68   PetscFunctionReturn(0);
69 }
70 EXTERN_C_END
71 
72 EXTERN_C_BEGIN
73 #undef __FUNCT__
74 #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky"
75 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
76 {
77   PC_Cholesky *lu;
78 
79   PetscFunctionBegin;
80   lu               = (PC_Cholesky*)pc->data;
81   lu->reuseordering = flag;
82   PetscFunctionReturn(0);
83 }
84 EXTERN_C_END
85 
86 EXTERN_C_BEGIN
87 #undef __FUNCT__
88 #define __FUNCT__ "PCFactorSetReuseFill_Cholesky"
89 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag)
90 {
91   PC_Cholesky *lu;
92 
93   PetscFunctionBegin;
94   lu = (PC_Cholesky*)pc->data;
95   lu->reusefill = flag;
96   PetscFunctionReturn(0);
97 }
98 EXTERN_C_END
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "PCSetFromOptions_Cholesky"
102 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
103 {
104   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
105   PetscErrorCode ierr;
106   PetscTruth     flg;
107   char           tname[256];
108   PetscFList     ordlist;
109 
110   PetscFunctionBegin;
111   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
112   ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr);
113   ierr = PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
114   if (flg) {
115     ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
116   }
117   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
118 
119   ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
120   if (flg) {
121     ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
122   }
123   ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
124   if (flg) {
125     ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
126   }
127 
128   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
129   ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
130   if (flg) {
131     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
132   }
133   ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
134   if (flg) {
135     ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
136   }
137   ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
138   ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
139   if (flg) {
140     ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
141   }
142   ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
143 
144   ierr = PetscOptionsTail();CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "PCView_Cholesky"
150 static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
151 {
152   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
153   PetscErrorCode ierr;
154   PetscTruth     iascii,isstring;
155 
156   PetscFunctionBegin;
157   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
158   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
159   if (iascii) {
160 
161     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");CHKERRQ(ierr);}
162     else             {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");CHKERRQ(ierr);}
163     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
164     if (lu->fact) {
165       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
166       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
167       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
168       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
169       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
170       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
171       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
172       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
173       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
174       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
175       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
176     }
177     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
178     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
179   } else if (isstring) {
180     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
181   } else {
182     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name);
183   }
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "PCFactorGetMatrix_Cholesky"
189 static PetscErrorCode PCFactorGetMatrix_Cholesky(PC pc,Mat *mat)
190 {
191   PC_Cholesky *dir = (PC_Cholesky*)pc->data;
192 
193   PetscFunctionBegin;
194   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
195   *mat = dir->fact;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "PCSetUp_Cholesky"
201 static PetscErrorCode PCSetUp_Cholesky(PC pc)
202 {
203   PetscErrorCode ierr;
204   PetscTruth     flg;
205   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
206 
207   PetscFunctionBegin;
208   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
209 
210   if (dir->inplace) {
211     if (dir->row && dir->col && (dir->row != dir->col)) {
212       ierr = ISDestroy(dir->row);CHKERRQ(ierr);
213       dir->row = 0;
214     }
215     if (dir->col) {
216       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
217       dir->col = 0;
218     }
219     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
220     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
221       ierr = ISDestroy(dir->col);CHKERRQ(ierr);
222       dir->col=0;
223     }
224     if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
225     ierr = MatCholeskyFactor(pc->pmat,dir->row,&dir->info);CHKERRQ(ierr);
226     dir->fact = pc->pmat;
227   } else {
228     MatInfo info;
229     if (!pc->setupcalled) {
230       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
231       /* check if dir->row == dir->col */
232       ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr);
233       if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
234       ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */
235       dir->col=0;
236 
237       ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
238       if (flg) {
239         PetscReal tol = 1.e-10;
240         ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
241         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
242       }
243       if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
244       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
245       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
246       dir->actualfill = info.fill_ratio_needed;
247       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
248     } else if (pc->flag != SAME_NONZERO_PATTERN) {
249       if (!dir->reuseordering) {
250         if (dir->row && dir->col && (dir->row != dir->col)) {
251           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
252           dir->row = 0;
253         }
254         if (dir->col) {
255           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
256           dir->col =0;
257         }
258         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
259         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
260           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
261           dir->col=0;
262         }
263         ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
264         if (flg) {
265           PetscReal tol = 1.e-10;
266           ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
267           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
268         }
269         if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
270       }
271       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
272       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
273       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
274       dir->actualfill = info.fill_ratio_needed;
275       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
276     }
277     ierr = MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
278   }
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "PCDestroy_Cholesky"
284 static PetscErrorCode PCDestroy_Cholesky(PC pc)
285 {
286   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
291   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
292   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
293   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
294   ierr = PetscFree(dir);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "PCApply_Cholesky"
300 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
301 {
302   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
303   PetscErrorCode ierr;
304 
305   PetscFunctionBegin;
306   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
307   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
308   PetscFunctionReturn(0);
309 }
310 
311 #undef __FUNCT__
312 #define __FUNCT__ "PCApplyTranspose_Cholesky"
313 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
314 {
315   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
316   PetscErrorCode ierr;
317 
318   PetscFunctionBegin;
319   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
320   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
321   PetscFunctionReturn(0);
322 }
323 
324 /* -----------------------------------------------------------------------------------*/
325 
326 EXTERN_C_BEGIN
327 #undef __FUNCT__
328 #define __FUNCT__ "PCFactorSetFill_Cholesky"
329 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill)
330 {
331   PC_Cholesky *dir;
332 
333   PetscFunctionBegin;
334   dir = (PC_Cholesky*)pc->data;
335   dir->info.fill = fill;
336   PetscFunctionReturn(0);
337 }
338 EXTERN_C_END
339 
340 EXTERN_C_BEGIN
341 #undef __FUNCT__
342 #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky"
343 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc)
344 {
345   PC_Cholesky *dir;
346 
347   PetscFunctionBegin;
348   dir = (PC_Cholesky*)pc->data;
349   dir->inplace = PETSC_TRUE;
350   PetscFunctionReturn(0);
351 }
352 EXTERN_C_END
353 
354 EXTERN_C_BEGIN
355 #undef __FUNCT__
356 #define __FUNCT__ "PCFactorSetMatOrderingType_Cholesky"
357 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Cholesky(PC pc,MatOrderingType ordering)
358 {
359   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
360   PetscErrorCode ierr;
361 
362   PetscFunctionBegin;
363   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
364   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
365   PetscFunctionReturn(0);
366 }
367 EXTERN_C_END
368 
369 /* -----------------------------------------------------------------------------------*/
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "PCFactorSetReuseOrdering"
373 /*@
374    PCFactorSetReuseOrdering - When similar matrices are factored, this
375    causes the ordering computed in the first factor to be used for all
376    following factors.
377 
378    Collective on PC
379 
380    Input Parameters:
381 +  pc - the preconditioner context
382 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
383 
384    Options Database Key:
385 .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
386 
387    Level: intermediate
388 
389 .keywords: PC, levels, reordering, factorization, incomplete, LU
390 
391 .seealso: PCFactorSetReuseFill()
392 @*/
393 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
394 {
395   PetscErrorCode ierr,(*f)(PC,PetscTruth);
396 
397   PetscFunctionBegin;
398   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
399   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
400   if (f) {
401     ierr = (*f)(pc,flag);CHKERRQ(ierr);
402   }
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "PCFactorSetReuseFill"
408 /*@
409    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
410    this causes later ones to use the fill ratio computed in the initial factorization.
411 
412    Collective on PC
413 
414    Input Parameters:
415 +  pc - the preconditioner context
416 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
417 
418    Options Database Key:
419 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
420 
421    Level: intermediate
422 
423 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
424 
425 .seealso: PCFactorSetReuseOrdering()
426 @*/
427 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
428 {
429   PetscErrorCode ierr,(*f)(PC,PetscTruth);
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
433   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
434   if (f) {
435     ierr = (*f)(pc,flag);CHKERRQ(ierr);
436   }
437   PetscFunctionReturn(0);
438 }
439 
440 /*MC
441    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
442 
443    Options Database Keys:
444 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
445 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
446 .  -pc_factor_fill <fill> - Sets fill amount
447 .  -pc_factor_in_place - Activates in-place factorization
448 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
449 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
450 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
451    is optional with PETSC_TRUE being the default
452 
453    Notes: Not all options work for all matrix formats
454 
455    Level: beginner
456 
457    Concepts: Cholesky factorization, direct solver
458 
459    Notes: Usually this will compute an "exact" solution in one iteration and does
460           not need a Krylov method (i.e. you can use -ksp_type preonly, or
461           KSPSetType(ksp,KSPPREONLY) for the Krylov method
462 
463 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
464            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
465            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
466 	   PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
467 
468 M*/
469 
470 EXTERN_C_BEGIN
471 #undef __FUNCT__
472 #define __FUNCT__ "PCCreate_Cholesky"
473 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
474 {
475   PetscErrorCode ierr;
476   PC_Cholesky    *dir;
477 
478   PetscFunctionBegin;
479   ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr);
480 
481   dir->fact                   = 0;
482   dir->inplace                = PETSC_FALSE;
483   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
484   dir->info.fill              = 5.0;
485   dir->info.shiftnz           = 0.0;
486   dir->info.shiftpd           = 0.0; /* false */
487   dir->info.shift_fraction    = 0.0;
488   dir->info.pivotinblocks     = 1.0;
489   dir->col                    = 0;
490   dir->row                    = 0;
491   ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
492   dir->reusefill        = PETSC_FALSE;
493   dir->reuseordering    = PETSC_FALSE;
494   pc->data              = (void*)dir;
495 
496   pc->ops->destroy           = PCDestroy_Cholesky;
497   pc->ops->apply             = PCApply_Cholesky;
498   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
499   pc->ops->setup             = PCSetUp_Cholesky;
500   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
501   pc->ops->view              = PCView_Cholesky;
502   pc->ops->applyrichardson   = 0;
503   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Cholesky;
504 
505   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky",
506                     PCFactorSetZeroPivot_Cholesky);CHKERRQ(ierr);
507   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky",
508                     PCFactorSetShiftNonzero_Cholesky);CHKERRQ(ierr);
509   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky",
510                     PCFactorSetShiftPd_Cholesky);CHKERRQ(ierr);
511 
512   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky",
513                     PCFactorSetFill_Cholesky);CHKERRQ(ierr);
514   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
515                     PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr);
516   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Cholesky",
517                     PCFactorSetMatOrderingType_Cholesky);CHKERRQ(ierr);
518   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
519                     PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr);
520   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
521                     PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr);
522   PetscFunctionReturn(0);
523 }
524 EXTERN_C_END
525