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