xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision 7d6bfa3b9d7db0ccd4cc481237114ca8dbb0dbff)
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 "../src/ksp/pc/impls/factor/factor.h"         /*I "petscpc.h" I*/
9 
10 typedef struct {
11   PC_Factor        hdr;
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   PetscTruth       reuseordering;    /* reuses previous reordering computed */
16   PetscTruth       reusefill;        /* reuse fill from previous Cholesky */
17 } PC_Cholesky;
18 
19 EXTERN_C_BEGIN
20 #undef __FUNCT__
21 #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky"
22 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
23 {
24   PC_Cholesky *lu;
25 
26   PetscFunctionBegin;
27   lu               = (PC_Cholesky*)pc->data;
28   lu->reuseordering = flag;
29   PetscFunctionReturn(0);
30 }
31 EXTERN_C_END
32 
33 EXTERN_C_BEGIN
34 #undef __FUNCT__
35 #define __FUNCT__ "PCFactorSetReuseFill_Cholesky"
36 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag)
37 {
38   PC_Cholesky *lu;
39 
40   PetscFunctionBegin;
41   lu = (PC_Cholesky*)pc->data;
42   lu->reusefill = flag;
43   PetscFunctionReturn(0);
44 }
45 EXTERN_C_END
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "PCSetFromOptions_Cholesky"
49 static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
50 {
51   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
52   PetscErrorCode ierr;
53   PetscTruth     flg;
54   char           tname[256], solvertype[64];
55   PetscFList     ordlist;
56 
57   PetscFunctionBegin;
58   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
59   ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr);
60     ierr = PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
61     if (flg) {
62       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
63     }
64     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr);
65 
66     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
67     if (flg) {
68       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
69     }
70     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
71     if (flg) {
72       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
73     }
74 
75     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
76     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr);
77     if (flg) {
78       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
79     }
80 
81     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
82     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific Cholesky solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
83     if (flg) {
84       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
85     }
86 
87     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
88     if (flg) {
89       ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
90     }
91     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);CHKERRQ(ierr);
92     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
93     if (flg) {
94       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
95     }
96     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)lu)->info.zeropivot,&((PC_Factor*)lu)->info.zeropivot,0);CHKERRQ(ierr);
97 
98   ierr = PetscOptionsTail();CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 #undef __FUNCT__
103 #define __FUNCT__ "PCView_Cholesky"
104 static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
105 {
106   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
107   PetscErrorCode ierr;
108   PetscTruth     iascii,isstring;
109 
110   PetscFunctionBegin;
111   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
112   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
113   if (iascii) {
114 
115     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");CHKERRQ(ierr);}
116     else             {ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");CHKERRQ(ierr);}
117     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);
118     if (((PC_Factor*)lu)->fact) {
119       ierr = PetscViewerASCIIPrintf(viewer,"  Cholesky: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
120       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
121       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
122       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
123       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
124       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
125       ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr);
126       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
127       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
128       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
129       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
130     }
131     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
132     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
133   } else if (isstring) {
134     ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
135   } else {
136     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name);
137   }
138   PetscFunctionReturn(0);
139 }
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) ((PC_Factor*)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,((PC_Factor*)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) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
168     ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
169     ((PC_Factor*)dir)->fact = pc->pmat;
170   } else {
171     MatInfo info;
172     if (!pc->setupcalled) {
173       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
174       /* check if dir->row == dir->col */
175       ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr);
176       if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
177       ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */
178       dir->col=0;
179 
180       ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
181       if (flg) {
182         PetscReal tol = 1.e-10;
183         ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
184         ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
185       }
186       if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
187       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
188       ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
189       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
190       dir->actualfill = info.fill_ratio_needed;
191       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
192     } else if (pc->flag != SAME_NONZERO_PATTERN) {
193       if (!dir->reuseordering) {
194         if (dir->row && dir->col && (dir->row != dir->col)) {
195           ierr = ISDestroy(dir->row);CHKERRQ(ierr);
196           dir->row = 0;
197         }
198         if (dir->col) {
199           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
200           dir->col =0;
201         }
202         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
203         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
204           ierr = ISDestroy(dir->col);CHKERRQ(ierr);
205           dir->col=0;
206         }
207         ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr);
208         if (flg) {
209           PetscReal tol = 1.e-10;
210           ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr);
211           ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr);
212         }
213         if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);}
214       }
215       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
216       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
217       ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
218       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
219       dir->actualfill = info.fill_ratio_needed;
220       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
221     }
222     ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "PCDestroy_Cholesky"
229 static PetscErrorCode PCDestroy_Cholesky(PC pc)
230 {
231   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
236   if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
237   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
238   ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
239   ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
240   ierr = PetscFree(dir);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "PCApply_Cholesky"
246 static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
247 {
248   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
249   PetscErrorCode ierr;
250 
251   PetscFunctionBegin;
252   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
253   else              {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
254   PetscFunctionReturn(0);
255 }
256 
257 #undef __FUNCT__
258 #define __FUNCT__ "PCApplyTranspose_Cholesky"
259 static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
260 {
261   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
262   PetscErrorCode ierr;
263 
264   PetscFunctionBegin;
265   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
266   else              {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);}
267   PetscFunctionReturn(0);
268 }
269 
270 /* -----------------------------------------------------------------------------------*/
271 
272 EXTERN_C_BEGIN
273 #undef __FUNCT__
274 #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky"
275 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc)
276 {
277   PC_Cholesky *dir;
278 
279   PetscFunctionBegin;
280   dir = (PC_Cholesky*)pc->data;
281   dir->inplace = PETSC_TRUE;
282   PetscFunctionReturn(0);
283 }
284 EXTERN_C_END
285 
286 /* -----------------------------------------------------------------------------------*/
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "PCFactorSetReuseOrdering"
290 /*@
291    PCFactorSetReuseOrdering - When similar matrices are factored, this
292    causes the ordering computed in the first factor to be used for all
293    following factors.
294 
295    Collective on PC
296 
297    Input Parameters:
298 +  pc - the preconditioner context
299 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
300 
301    Options Database Key:
302 .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
303 
304    Level: intermediate
305 
306 .keywords: PC, levels, reordering, factorization, incomplete, LU
307 
308 .seealso: PCFactorSetReuseFill()
309 @*/
310 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
311 {
312   PetscErrorCode ierr,(*f)(PC,PetscTruth);
313 
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
316   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
317   if (f) {
318     ierr = (*f)(pc,flag);CHKERRQ(ierr);
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 
324 /*MC
325    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
326 
327    Options Database Keys:
328 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
329 .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
330 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
331 .  -pc_factor_fill <fill> - Sets fill amount
332 .  -pc_factor_in_place - Activates in-place factorization
333 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
334 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
335 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
336    is optional with PETSC_TRUE being the default
337 
338    Notes: Not all options work for all matrix formats
339 
340    Level: beginner
341 
342    Concepts: Cholesky factorization, direct solver
343 
344    Notes: Usually this will compute an "exact" solution in one iteration and does
345           not need a Krylov method (i.e. you can use -ksp_type preonly, or
346           KSPSetType(ksp,KSPPREONLY) for the Krylov method
347 
348 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
349            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
350            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
351 	   PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
352 
353 M*/
354 
355 EXTERN_C_BEGIN
356 #undef __FUNCT__
357 #define __FUNCT__ "PCCreate_Cholesky"
358 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
359 {
360   PetscErrorCode ierr;
361   PC_Cholesky    *dir;
362 
363   PetscFunctionBegin;
364   ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr);
365 
366   ((PC_Factor*)dir)->fact                   = 0;
367   dir->inplace                = PETSC_FALSE;
368   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
369   ((PC_Factor*)dir)->info.fill              = 5.0;
370   ((PC_Factor*)dir)->info.shiftnz           = 0.0;
371   ((PC_Factor*)dir)->info.shiftpd           = 0.0; /* false */
372   ((PC_Factor*)dir)->info.pivotinblocks     = 1.0;
373   dir->col                    = 0;
374   dir->row                    = 0;
375   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
376   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
377   dir->reusefill        = PETSC_FALSE;
378   dir->reuseordering    = PETSC_FALSE;
379   pc->data              = (void*)dir;
380 
381   pc->ops->destroy           = PCDestroy_Cholesky;
382   pc->ops->apply             = PCApply_Cholesky;
383   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
384   pc->ops->setup             = PCSetUp_Cholesky;
385   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
386   pc->ops->view              = PCView_Cholesky;
387   pc->ops->applyrichardson   = 0;
388   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
389 
390   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
391                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
392   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
393                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
394   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
395                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
396   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
397                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
398   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
399                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
400 
401   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
402                     PCFactorSetFill_Factor);CHKERRQ(ierr);
403   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
404                     PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr);
405   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
406                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
407   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
408                     PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr);
409   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
410                     PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 EXTERN_C_END
414