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