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