xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
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 
9 #include "../src/ksp/pc/impls/factor/lu/lu.h"  /*I "petscpc.h" I*/
10 
11 EXTERN_C_BEGIN
12 #undef __FUNCT__
13 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
14 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
15 {
16   PC_LU *lu = (PC_LU*)pc->data;
17 
18   PetscFunctionBegin;
19   lu->nonzerosalongdiagonal = PETSC_TRUE;
20   if (z == PETSC_DECIDE) {
21     lu->nonzerosalongdiagonaltol = 1.e-10;
22   } else {
23     lu->nonzerosalongdiagonaltol = z;
24   }
25   PetscFunctionReturn(0);
26 }
27 EXTERN_C_END
28 
29 EXTERN_C_BEGIN
30 #undef __FUNCT__
31 #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
32 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
33 {
34   PC_LU *lu = (PC_LU*)pc->data;
35 
36   PetscFunctionBegin;
37   lu->reuseordering = flag;
38   PetscFunctionReturn(0);
39 }
40 EXTERN_C_END
41 
42 EXTERN_C_BEGIN
43 #undef __FUNCT__
44 #define __FUNCT__ "PCFactorSetReuseFill_LU"
45 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
46 {
47   PC_LU *lu = (PC_LU*)pc->data;
48 
49   PetscFunctionBegin;
50   lu->reusefill = flag;
51   PetscFunctionReturn(0);
52 }
53 EXTERN_C_END
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "PCSetFromOptions_LU"
57 static PetscErrorCode PCSetFromOptions_LU(PC pc)
58 {
59   PC_LU           *lu = (PC_LU*)pc->data;
60   PetscErrorCode  ierr;
61   PetscTruth      flg,set;
62   char            tname[256], solvertype[64];
63   PetscFList      ordlist;
64   PetscReal       tol;
65 
66   PetscFunctionBegin;
67   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
68   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
69     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
70     if (flg) {
71       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
72     }
73     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr);
74 
75     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
76     if (flg) {
77         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
78     }
79     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);CHKERRQ(ierr);
80     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
81     if (flg) {
82       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
83     }
84     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);
85 
86     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
87     if (flg) {
88       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
89     }
90     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
91     if (flg) {
92       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
93     }
94 
95     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
96     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr);
97     if (flg) {
98       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
99     }
100 
101     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
102     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
103     if (flg) {
104       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
105     }
106 
107     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
108     if (flg) {
109       tol = PETSC_DECIDE;
110       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
111       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
112     }
113 
114     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",((PC_Factor*)lu)->info.dtcol,&((PC_Factor*)lu)->info.dtcol,&flg);CHKERRQ(ierr);
115 
116     flg = ((PC_Factor*)lu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
117     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
118     if (set) {
119       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
120     }
121   ierr = PetscOptionsTail();CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PCView_LU"
127 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
128 {
129   PC_LU          *lu = (PC_LU*)pc->data;
130   PetscErrorCode ierr;
131   PetscTruth     iascii,isstring;
132 
133   PetscFunctionBegin;
134   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
135   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
136   if (iascii) {
137 
138     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
139     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
140     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);
141     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr);
142     if (((PC_Factor*)lu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
143     if (((PC_Factor*)lu)->fact) {
144       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
145       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
146       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
147       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
148       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
149       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
150       ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr);
151       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
152       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
153       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
154       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
155     }
156     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
157     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
158   } else if (isstring) {
159     ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
160   } else {
161     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
162   }
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "PCSetUp_LU"
168 static PetscErrorCode PCSetUp_LU(PC pc)
169 {
170   PetscErrorCode ierr;
171   PC_LU          *dir = (PC_LU*)pc->data;
172 
173   PetscFunctionBegin;
174   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
175 
176   if (dir->inplace) {
177     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
178     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
179     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
180     if (dir->row) {
181       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
182       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
183     }
184     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
185     ((PC_Factor*)dir)->fact = pc->pmat;
186   } else {
187     MatInfo info;
188     if (!pc->setupcalled) {
189       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
190       if (dir->nonzerosalongdiagonal) {
191         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
192       }
193       if (dir->row) {
194         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
195         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
196       }
197       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
198       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
199       ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
200       dir->actualfill = info.fill_ratio_needed;
201       ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr);
202     } else if (pc->flag != SAME_NONZERO_PATTERN) {
203       if (!dir->reuseordering) {
204         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
205         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
206         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
207         if (dir->nonzerosalongdiagonal) {
208           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
209         }
210         if (dir->row) {
211           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
212           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
213         }
214       }
215       ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
216       ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
217       ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((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 = MatLUFactorNumeric(((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_LU"
229 static PetscErrorCode PCDestroy_LU(PC pc)
230 {
231   PC_LU          *dir = (PC_LU*)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 && dir->col && dir->row != dir->col) {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_LU"
246 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
247 {
248   PC_LU          *dir = (PC_LU*)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_LU"
259 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
260 {
261   PC_LU          *dir = (PC_LU*)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_LU"
275 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
276 {
277   PC_LU *dir = (PC_LU*)pc->data;
278 
279   PetscFunctionBegin;
280   dir->inplace = PETSC_TRUE;
281   PetscFunctionReturn(0);
282 }
283 EXTERN_C_END
284 
285 /* ------------------------------------------------------------------------ */
286 
287 /*MC
288    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
289 
290    Options Database Keys:
291 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
292 .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
293 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
294 .  -pc_factor_fill <fill> - Sets fill amount
295 .  -pc_factor_in_place - Activates in-place factorization
296 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
297 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
298                                          stability of factorization.
299 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
300 .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
301         is optional with PETSC_TRUE being the default
302 -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
303         diagonal.
304 
305    Notes: Not all options work for all matrix formats
306           Run with -help to see additional options for particular matrix formats or factorization
307           algorithms
308 
309    Level: beginner
310 
311    Concepts: LU factorization, direct solver
312 
313    Notes: Usually this will compute an "exact" solution in one iteration and does
314           not need a Krylov method (i.e. you can use -ksp_type preonly, or
315           KSPSetType(ksp,KSPPREONLY) for the Krylov method
316 
317 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
318            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
319            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
320            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
321 M*/
322 
323 EXTERN_C_BEGIN
324 #undef __FUNCT__
325 #define __FUNCT__ "PCCreate_LU"
326 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
327 {
328   PetscErrorCode ierr;
329   PetscMPIInt    size;
330   PC_LU          *dir;
331 
332   PetscFunctionBegin;
333   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
334 
335   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
336   ((PC_Factor*)dir)->fact                  = 0;
337   dir->inplace               = PETSC_FALSE;
338   dir->nonzerosalongdiagonal = PETSC_FALSE;
339 
340   ((PC_Factor*)dir)->info.fill           = 5.0;
341   ((PC_Factor*)dir)->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
342   ((PC_Factor*)dir)->info.shiftnz        = 0.0;
343   ((PC_Factor*)dir)->info.zeropivot      = 1.e-12;
344   ((PC_Factor*)dir)->info.pivotinblocks  = 1.0;
345   ((PC_Factor*)dir)->info.shiftpd        = 0.0; /* false */
346   dir->col                 = 0;
347   dir->row                 = 0;
348 
349   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
350   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
351   if (size == 1) {
352     ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
353   } else {
354     ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
355   }
356   dir->reusefill        = PETSC_FALSE;
357   dir->reuseordering    = PETSC_FALSE;
358   pc->data              = (void*)dir;
359 
360   pc->ops->destroy           = PCDestroy_LU;
361   pc->ops->apply             = PCApply_LU;
362   pc->ops->applytranspose    = PCApplyTranspose_LU;
363   pc->ops->setup             = PCSetUp_LU;
364   pc->ops->setfromoptions    = PCSetFromOptions_LU;
365   pc->ops->view              = PCView_LU;
366   pc->ops->applyrichardson   = 0;
367   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
368 
369   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
370                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
371   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
372                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
373   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
374                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
375   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
376                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
377   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
378                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
379 
380   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
381                     PCFactorSetFill_Factor);CHKERRQ(ierr);
382   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
383                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
384   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
385                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
386   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
387                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
388   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
389                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
390   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_Factor",
391                     PCFactorSetPivoting_Factor);CHKERRQ(ierr);
392   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
393                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
394   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
395                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
396   PetscFunctionReturn(0);
397 }
398 EXTERN_C_END
399