xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision b24902e06ab141841ce6546a773244be2474cd18)
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 = MatGetFactor(pc->pmat,"petsc",MAT_FACTOR_CHOLESKY,&dir->fact);CHKERRQ(ierr);
245       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
246       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
247       dir->actualfill = info.fill_ratio_needed;
248       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
249     } else if (pc->flag != SAME_NONZERO_PATTERN) {
250       if (!dir->reuseordering) {
251         if (dir->row && dir->col && (dir->row != dir->col)) {
252           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
253           dir->row = 0;
254         }
255         if (dir->col) {
256           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
257           dir->col =0;
258         }
259         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
260         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
261           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
262           dir->col=0;
263         }
264         ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
265         if (flg) {
266           PetscReal tol = 1.e-10;
267           ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
268           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
269         }
270         if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
271       }
272       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
273       ierr = MatGetFactor(pc->pmat,"petsc",MAT_FACTOR_CHOLESKY,&dir->fact);CHKERRQ(ierr);
274       ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr);
275       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
276       dir->actualfill = info.fill_ratio_needed;
277       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
278     }
279     ierr = MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr);
280   }
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "PCDestroy_Cholesky"
286 static PetscErrorCode PCDestroy_Cholesky(PC pc)
287 {
288   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
289   PetscErrorCode ierr;
290 
291   PetscFunctionBegin;
292   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
293   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
294   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
295   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
296   ierr = PetscFree(dir);CHKERRQ(ierr);
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "PCApply_Cholesky"
302 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
303 {
304   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
309   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "PCApplyTranspose_Cholesky"
315 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
316 {
317   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
318   PetscErrorCode ierr;
319 
320   PetscFunctionBegin;
321   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
322   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
323   PetscFunctionReturn(0);
324 }
325 
326 /* -----------------------------------------------------------------------------------*/
327 
328 EXTERN_C_BEGIN
329 #undef __FUNCT__
330 #define __FUNCT__ "PCFactorSetFill_Cholesky"
331 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill)
332 {
333   PC_Cholesky *dir;
334 
335   PetscFunctionBegin;
336   dir = (PC_Cholesky*)pc->data;
337   dir->info.fill = fill;
338   PetscFunctionReturn(0);
339 }
340 EXTERN_C_END
341 
342 EXTERN_C_BEGIN
343 #undef __FUNCT__
344 #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky"
345 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc)
346 {
347   PC_Cholesky *dir;
348 
349   PetscFunctionBegin;
350   dir = (PC_Cholesky*)pc->data;
351   dir->inplace = PETSC_TRUE;
352   PetscFunctionReturn(0);
353 }
354 EXTERN_C_END
355 
356 EXTERN_C_BEGIN
357 #undef __FUNCT__
358 #define __FUNCT__ "PCFactorSetMatOrderingType_Cholesky"
359 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Cholesky(PC pc,MatOrderingType ordering)
360 {
361   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
366   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
367   PetscFunctionReturn(0);
368 }
369 EXTERN_C_END
370 
371 /* -----------------------------------------------------------------------------------*/
372 
373 #undef __FUNCT__
374 #define __FUNCT__ "PCFactorSetReuseOrdering"
375 /*@
376    PCFactorSetReuseOrdering - When similar matrices are factored, this
377    causes the ordering computed in the first factor to be used for all
378    following factors.
379 
380    Collective on PC
381 
382    Input Parameters:
383 +  pc - the preconditioner context
384 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
385 
386    Options Database Key:
387 .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
388 
389    Level: intermediate
390 
391 .keywords: PC, levels, reordering, factorization, incomplete, LU
392 
393 .seealso: PCFactorSetReuseFill()
394 @*/
395 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
396 {
397   PetscErrorCode ierr,(*f)(PC,PetscTruth);
398 
399   PetscFunctionBegin;
400   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
401   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
402   if (f) {
403     ierr = (*f)(pc,flag);CHKERRQ(ierr);
404   }
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCFactorSetReuseFill"
410 /*@
411    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
412    this causes later ones to use the fill ratio computed in the initial factorization.
413 
414    Collective on PC
415 
416    Input Parameters:
417 +  pc - the preconditioner context
418 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
419 
420    Options Database Key:
421 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
422 
423    Level: intermediate
424 
425 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
426 
427 .seealso: PCFactorSetReuseOrdering()
428 @*/
429 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
430 {
431   PetscErrorCode ierr,(*f)(PC,PetscTruth);
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
435   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
436   if (f) {
437     ierr = (*f)(pc,flag);CHKERRQ(ierr);
438   }
439   PetscFunctionReturn(0);
440 }
441 
442 /*MC
443    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
444 
445    Options Database Keys:
446 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
447 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
448 .  -pc_factor_fill <fill> - Sets fill amount
449 .  -pc_factor_in_place - Activates in-place factorization
450 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
451 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
452 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
453    is optional with PETSC_TRUE being the default
454 
455    Notes: Not all options work for all matrix formats
456 
457    Level: beginner
458 
459    Concepts: Cholesky factorization, direct solver
460 
461    Notes: Usually this will compute an "exact" solution in one iteration and does
462           not need a Krylov method (i.e. you can use -ksp_type preonly, or
463           KSPSetType(ksp,KSPPREONLY) for the Krylov method
464 
465 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
466            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
467            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
468 	   PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
469 
470 M*/
471 
472 EXTERN_C_BEGIN
473 #undef __FUNCT__
474 #define __FUNCT__ "PCCreate_Cholesky"
475 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
476 {
477   PetscErrorCode ierr;
478   PC_Cholesky    *dir;
479 
480   PetscFunctionBegin;
481   ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr);
482 
483   dir->fact                   = 0;
484   dir->inplace                = PETSC_FALSE;
485   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
486   dir->info.fill              = 5.0;
487   dir->info.shiftnz           = 0.0;
488   dir->info.shiftpd           = 0.0; /* false */
489   dir->info.shift_fraction    = 0.0;
490   dir->info.pivotinblocks     = 1.0;
491   dir->col                    = 0;
492   dir->row                    = 0;
493   ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
494   dir->reusefill        = PETSC_FALSE;
495   dir->reuseordering    = PETSC_FALSE;
496   pc->data              = (void*)dir;
497 
498   pc->ops->destroy           = PCDestroy_Cholesky;
499   pc->ops->apply             = PCApply_Cholesky;
500   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
501   pc->ops->setup             = PCSetUp_Cholesky;
502   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
503   pc->ops->view              = PCView_Cholesky;
504   pc->ops->applyrichardson   = 0;
505   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Cholesky;
506 
507   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky",
508                     PCFactorSetZeroPivot_Cholesky);CHKERRQ(ierr);
509   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky",
510                     PCFactorSetShiftNonzero_Cholesky);CHKERRQ(ierr);
511   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky",
512                     PCFactorSetShiftPd_Cholesky);CHKERRQ(ierr);
513 
514   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky",
515                     PCFactorSetFill_Cholesky);CHKERRQ(ierr);
516   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
517                     PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr);
518   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Cholesky",
519                     PCFactorSetMatOrderingType_Cholesky);CHKERRQ(ierr);
520   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
521                     PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr);
522   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
523                     PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr);
524   PetscFunctionReturn(0);
525 }
526 EXTERN_C_END
527