xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
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 "private/pcimpl.h"                /*I "petscpc.h" I*/
10 #include "src/ksp/pc/impls/factor/lu/lu.h"
11 
12 EXTERN_C_BEGIN
13 #undef __FUNCT__
14 #define __FUNCT__ "PCFactorSetMatSolverPackage_LU"
15 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_LU(PC pc,const MatSolverPackage stype)
16 {
17   PetscErrorCode ierr;
18   PC_LU          *lu = (PC_LU*)pc->data;
19 
20   PetscFunctionBegin;
21   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
22   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 EXTERN_C_END
26 
27 EXTERN_C_BEGIN
28 #undef __FUNCT__
29 #define __FUNCT__ "PCFactorSetZeroPivot_LU"
30 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z)
31 {
32   PC_LU *lu = (PC_LU*)pc->data;
33 
34   PetscFunctionBegin;
35   lu->info.zeropivot = z;
36   PetscFunctionReturn(0);
37 }
38 EXTERN_C_END
39 
40 EXTERN_C_BEGIN
41 #undef __FUNCT__
42 #define __FUNCT__ "PCFactorSetShiftNonzero_LU"
43 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift)
44 {
45   PC_LU *dir = (PC_LU*)pc->data;
46 
47   PetscFunctionBegin;
48   if (shift == (PetscReal) PETSC_DECIDE) {
49     dir->info.shiftnz = 1.e-12;
50   } else {
51     dir->info.shiftnz = shift;
52   }
53   PetscFunctionReturn(0);
54 }
55 EXTERN_C_END
56 
57 EXTERN_C_BEGIN
58 #undef __FUNCT__
59 #define __FUNCT__ "PCFactorSetShiftPd_LU"
60 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift)
61 {
62   PC_LU *dir = (PC_LU*)pc->data;
63 
64   PetscFunctionBegin;
65   if (shift) {
66     dir->info.shift_fraction = 0.0;
67     dir->info.shiftpd = 1.0;
68   } else {
69     dir->info.shiftpd = 0.0;
70   }
71   PetscFunctionReturn(0);
72 }
73 EXTERN_C_END
74 
75 EXTERN_C_BEGIN
76 #undef __FUNCT__
77 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
78 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
79 {
80   PC_LU *lu = (PC_LU*)pc->data;
81 
82   PetscFunctionBegin;
83   lu->nonzerosalongdiagonal = PETSC_TRUE;
84   if (z == PETSC_DECIDE) {
85     lu->nonzerosalongdiagonaltol = 1.e-10;
86   } else {
87     lu->nonzerosalongdiagonaltol = z;
88   }
89   PetscFunctionReturn(0);
90 }
91 EXTERN_C_END
92 
93 EXTERN_C_BEGIN
94 #undef __FUNCT__
95 #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
96 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag)
97 {
98   PC_LU *lu = (PC_LU*)pc->data;
99 
100   PetscFunctionBegin;
101   lu->reuseordering = flag;
102   PetscFunctionReturn(0);
103 }
104 EXTERN_C_END
105 
106 EXTERN_C_BEGIN
107 #undef __FUNCT__
108 #define __FUNCT__ "PCFactorSetReuseFill_LU"
109 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag)
110 {
111   PC_LU *lu = (PC_LU*)pc->data;
112 
113   PetscFunctionBegin;
114   lu->reusefill = flag;
115   PetscFunctionReturn(0);
116 }
117 EXTERN_C_END
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "PCSetFromOptions_LU"
121 static PetscErrorCode PCSetFromOptions_LU(PC pc)
122 {
123   PC_LU           *lu = (PC_LU*)pc->data;
124   PetscErrorCode  ierr;
125   PetscTruth      flg,set;
126   char            tname[256], solvertype[64];
127   PetscFList      ordlist;
128   PetscReal       tol;
129 
130   PetscFunctionBegin;
131   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
132   ierr = PetscOptionsHead("LU options");CHKERRQ(ierr);
133     ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr);
134     if (flg) {
135       ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
136     }
137     ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr);
138 
139     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
140     if (flg) {
141         ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
142     }
143     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr);
144     ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr);
145     if (flg) {
146       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
147     }
148     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr);
149 
150     ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr);
151     if (flg) {
152       ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
153     }
154     ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr);
155     if (flg) {
156       ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
157     }
158 
159     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
160     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr);
161     if (flg) {
162       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
163     }
164 
165     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
166     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
167     if (flg) {
168       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
169     }
170 
171     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
172     if (flg) {
173       tol = PETSC_DECIDE;
174       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
175       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
176     }
177 
178     ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr);
179 
180     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
181     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
182     if (set) {
183       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
184     }
185   ierr = PetscOptionsTail();CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "PCView_LU"
191 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
192 {
193   PC_LU          *lu = (PC_LU*)pc->data;
194   PetscErrorCode ierr;
195   PetscTruth     iascii,isstring;
196 
197   PetscFunctionBegin;
198   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
199   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
200   if (iascii) {
201 
202     if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);}
203     else             {ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);}
204     ierr = PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr);
205     ierr = PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr);
206     if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");CHKERRQ(ierr);}
207     if (lu->fact) {
208       ierr = PetscViewerASCIIPrintf(viewer,"  LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr);
209       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
210       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
211       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
212       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
213       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
214       ierr = MatView(lu->fact,viewer);CHKERRQ(ierr);
215       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
216       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
217       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
218       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
219     }
220     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
221     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
222   } else if (isstring) {
223     ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
224   } else {
225     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
226   }
227   PetscFunctionReturn(0);
228 }
229 
230 #undef __FUNCT__
231 #define __FUNCT__ "PCFactorGetMatrix_LU"
232 static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat)
233 {
234   PC_LU *dir = (PC_LU*)pc->data;
235 
236   PetscFunctionBegin;
237   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
238   *mat = dir->fact;
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "PCSetUp_LU"
244 static PetscErrorCode PCSetUp_LU(PC pc)
245 {
246   PetscErrorCode ierr;
247   PC_LU          *dir = (PC_LU*)pc->data;
248 
249   PetscFunctionBegin;
250   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
251 
252   if (dir->inplace) {
253     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
254     if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
255     ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
256     if (dir->row) {
257       ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
258       ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
259     }
260     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
261     dir->fact = pc->pmat;
262   } else {
263     MatInfo info;
264     if (!pc->setupcalled) {
265       ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
266       if (dir->nonzerosalongdiagonal) {
267         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
268       }
269       if (dir->row) {
270         ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
271         ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
272       }
273       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
274       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
275       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
276       dir->actualfill = info.fill_ratio_needed;
277       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
278     } else if (pc->flag != SAME_NONZERO_PATTERN) {
279       if (!dir->reuseordering) {
280         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
281         if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
282         ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
283         if (dir->nonzerosalongdiagonal) {
284           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
285         }
286         if (dir->row) {
287           ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);
288           ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr);
289         }
290       }
291       ierr = MatDestroy(dir->fact);CHKERRQ(ierr);
292       ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr);
293       ierr = MatLUFactorSymbolic(dir->fact,pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr);
294       ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
295       dir->actualfill = info.fill_ratio_needed;
296       ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr);
297     }
298     ierr = MatLUFactorNumeric(dir->fact,pc->pmat,&dir->info);CHKERRQ(ierr);
299   }
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "PCDestroy_LU"
305 static PetscErrorCode PCDestroy_LU(PC pc)
306 {
307   PC_LU          *dir = (PC_LU*)pc->data;
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);}
312   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);}
313   if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);}
314   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
315   ierr = PetscStrfree(dir->solvertype);CHKERRQ(ierr);
316   ierr = PetscFree(dir);CHKERRQ(ierr);
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "PCApply_LU"
322 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
323 {
324   PC_LU          *dir = (PC_LU*)pc->data;
325   PetscErrorCode ierr;
326 
327   PetscFunctionBegin;
328   if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);}
329   else              {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);}
330   PetscFunctionReturn(0);
331 }
332 
333 #undef __FUNCT__
334 #define __FUNCT__ "PCApplyTranspose_LU"
335 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
336 {
337   PC_LU          *dir = (PC_LU*)pc->data;
338   PetscErrorCode ierr;
339 
340   PetscFunctionBegin;
341   if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);}
342   else              {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);}
343   PetscFunctionReturn(0);
344 }
345 
346 /* -----------------------------------------------------------------------------------*/
347 
348 EXTERN_C_BEGIN
349 #undef __FUNCT__
350 #define __FUNCT__ "PCFactorSetFill_LU"
351 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill)
352 {
353   PC_LU *dir = (PC_LU*)pc->data;
354 
355   PetscFunctionBegin;
356   dir->info.fill = fill;
357   PetscFunctionReturn(0);
358 }
359 EXTERN_C_END
360 
361 EXTERN_C_BEGIN
362 #undef __FUNCT__
363 #define __FUNCT__ "PCFactorSetUseInPlace_LU"
364 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc)
365 {
366   PC_LU *dir = (PC_LU*)pc->data;
367 
368   PetscFunctionBegin;
369   dir->inplace = PETSC_TRUE;
370   PetscFunctionReturn(0);
371 }
372 EXTERN_C_END
373 
374 EXTERN_C_BEGIN
375 #undef __FUNCT__
376 #define __FUNCT__ "PCFactorSetMatOrderingType_LU"
377 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering)
378 {
379   PC_LU          *dir = (PC_LU*)pc->data;
380   PetscErrorCode ierr;
381 
382   PetscFunctionBegin;
383   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
384   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
385   PetscFunctionReturn(0);
386 }
387 EXTERN_C_END
388 
389 EXTERN_C_BEGIN
390 #undef __FUNCT__
391 #define __FUNCT__ "PCFactorSetPivoting_LU"
392 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol)
393 {
394   PC_LU *dir = (PC_LU*)pc->data;
395 
396   PetscFunctionBegin;
397   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
398   dir->info.dtcol = dtcol;
399   PetscFunctionReturn(0);
400 }
401 EXTERN_C_END
402 
403 EXTERN_C_BEGIN
404 #undef __FUNCT__
405 #define __FUNCT__ "PCFactorSetPivotInBlocks_LU"
406 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
407 {
408   PC_LU *dir = (PC_LU*)pc->data;
409 
410   PetscFunctionBegin;
411   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
412   PetscFunctionReturn(0);
413 }
414 EXTERN_C_END
415 
416 /* -----------------------------------------------------------------------------------*/
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
420 /*@
421    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
422 
423    Collective on PC
424 
425    Input Parameters:
426 +  pc - the preconditioner context
427 -  tol - diagonal entries smaller than this in absolute value are considered zero
428 
429    Options Database Key:
430 .  -pc_factor_nonzeros_along_diagonal
431 
432    Level: intermediate
433 
434 .keywords: PC, set, factorization, direct, fill
435 
436 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
437 @*/
438 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
439 {
440   PetscErrorCode ierr,(*f)(PC,PetscReal);
441 
442   PetscFunctionBegin;
443   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
444   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
445   if (f) {
446     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
447   }
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "PCFactorSetMatSolverPackage"
453 /*@
454    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
455 
456    Collective on PC
457 
458    Input Parameters:
459 +  pc - the preconditioner context
460 -  stype - for example, spooles, superlu, superlu_d
461 
462    Options Database Key:
463 .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
464 
465    Level: intermediate
466 
467    Note:
468      By default this will use the PETSc factorization if it exists
469 
470 .keywords: PC, set, factorization, direct, fill
471 
472 .seealso: MatGetFactor(), MatSolverPackage
473 
474 @*/
475 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
476 {
477   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
478 
479   PetscFunctionBegin;
480   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
481   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
482   if (f) {
483     ierr = (*f)(pc,stype);CHKERRQ(ierr);
484   }
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "PCFactorSetFill"
490 /*@
491    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
492    fill = number nonzeros in factor/number nonzeros in original matrix.
493 
494    Collective on PC
495 
496    Input Parameters:
497 +  pc - the preconditioner context
498 -  fill - amount of expected fill
499 
500    Options Database Key:
501 .  -pc_factor_fill <fill> - Sets fill amount
502 
503    Level: intermediate
504 
505    Note:
506    For sparse matrix factorizations it is difficult to predict how much
507    fill to expect. By running with the option -info PETSc will print the
508    actual amount of fill used; allowing you to set the value accurately for
509    future runs. Default PETSc uses a value of 5.0
510 
511 .keywords: PC, set, factorization, direct, fill
512 
513 @*/
514 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
515 {
516   PetscErrorCode ierr,(*f)(PC,PetscReal);
517 
518   PetscFunctionBegin;
519   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
520   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
521   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
522   if (f) {
523     ierr = (*f)(pc,fill);CHKERRQ(ierr);
524   }
525   PetscFunctionReturn(0);
526 }
527 
528 #undef __FUNCT__
529 #define __FUNCT__ "PCFactorSetUseInPlace"
530 /*@
531    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
532    For dense matrices, this enables the solution of much larger problems.
533    For sparse matrices the factorization cannot be done truly in-place
534    so this does not save memory during the factorization, but after the matrix
535    is factored, the original unfactored matrix is freed, thus recovering that
536    space.
537 
538    Collective on PC
539 
540    Input Parameters:
541 .  pc - the preconditioner context
542 
543    Options Database Key:
544 .  -pc_factor_in_place - Activates in-place factorization
545 
546    Notes:
547    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
548    a different matrix is provided for the multiply and the preconditioner in
549    a call to KSPSetOperators().
550    This is because the Krylov space methods require an application of the
551    matrix multiplication, which is not possible here because the matrix has
552    been factored in-place, replacing the original matrix.
553 
554    Level: intermediate
555 
556 .keywords: PC, set, factorization, direct, inplace, in-place, LU
557 
558 .seealso: PCILUSetUseInPlace()
559 @*/
560 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
561 {
562   PetscErrorCode ierr,(*f)(PC);
563 
564   PetscFunctionBegin;
565   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
566   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
567   if (f) {
568     ierr = (*f)(pc);CHKERRQ(ierr);
569   }
570   PetscFunctionReturn(0);
571 }
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "PCFactorSetMatOrderingType"
575 /*@C
576     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
577     be used in the LU factorization.
578 
579     Collective on PC
580 
581     Input Parameters:
582 +   pc - the preconditioner context
583 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
584 
585     Options Database Key:
586 .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
587 
588     Level: intermediate
589 
590     Notes: nested dissection is used by default
591 
592     For Cholesky and ICC and the SBAIJ format reorderings are not available,
593     since only the upper triangular part of the matrix is stored. You can use the
594     SeqAIJ format in this case to get reorderings.
595 
596 @*/
597 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
598 {
599   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
600 
601   PetscFunctionBegin;
602   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
603   if (f) {
604     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
605   }
606   PetscFunctionReturn(0);
607 }
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "PCFactorSetPivoting"
611 /*@
612     PCFactorSetPivoting - Determines when pivoting is done during LU.
613       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
614       it is never done. For the Matlab and SuperLU factorization this is used.
615 
616     Collective on PC
617 
618     Input Parameters:
619 +   pc - the preconditioner context
620 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
621 
622     Options Database Key:
623 .   -pc_factor_pivoting <dtcol>
624 
625     Level: intermediate
626 
627 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
628 @*/
629 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
630 {
631   PetscErrorCode ierr,(*f)(PC,PetscReal);
632 
633   PetscFunctionBegin;
634   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
635   if (f) {
636     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
637   }
638   PetscFunctionReturn(0);
639 }
640 
641 #undef __FUNCT__
642 #define __FUNCT__ "PCFactorSetPivotInBlocks"
643 /*@
644     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
645       with BAIJ or SBAIJ matrices
646 
647     Collective on PC
648 
649     Input Parameters:
650 +   pc - the preconditioner context
651 -   pivot - PETSC_TRUE or PETSC_FALSE
652 
653     Options Database Key:
654 .   -pc_factor_pivot_in_blocks <true,false>
655 
656     Level: intermediate
657 
658 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
659 @*/
660 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
661 {
662   PetscErrorCode ierr,(*f)(PC,PetscTruth);
663 
664   PetscFunctionBegin;
665   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
666   if (f) {
667     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
668   }
669   PetscFunctionReturn(0);
670 }
671 
672 /* ------------------------------------------------------------------------ */
673 
674 /*MC
675    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
676 
677    Options Database Keys:
678 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
679 .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
680 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
681 .  -pc_factor_fill <fill> - Sets fill amount
682 .  -pc_factor_in_place - Activates in-place factorization
683 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
684 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
685                                          stability of factorization.
686 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
687 .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
688         is optional with PETSC_TRUE being the default
689 -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
690         diagonal.
691 
692    Notes: Not all options work for all matrix formats
693           Run with -help to see additional options for particular matrix formats or factorization
694           algorithms
695 
696    Level: beginner
697 
698    Concepts: LU factorization, direct solver
699 
700    Notes: Usually this will compute an "exact" solution in one iteration and does
701           not need a Krylov method (i.e. you can use -ksp_type preonly, or
702           KSPSetType(ksp,KSPPREONLY) for the Krylov method
703 
704 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
705            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
706            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(),
707            PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal()
708 M*/
709 
710 EXTERN_C_BEGIN
711 #undef __FUNCT__
712 #define __FUNCT__ "PCCreate_LU"
713 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc)
714 {
715   PetscErrorCode ierr;
716   PetscMPIInt    size;
717   PC_LU          *dir;
718 
719   PetscFunctionBegin;
720   ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr);
721 
722   ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr);
723   dir->fact                  = 0;
724   dir->inplace               = PETSC_FALSE;
725   dir->nonzerosalongdiagonal = PETSC_FALSE;
726 
727   dir->info.fill           = 5.0;
728   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
729   dir->info.shiftnz        = 0.0;
730   dir->info.zeropivot      = 1.e-12;
731   dir->info.pivotinblocks  = 1.0;
732   dir->info.shiftpd        = 0.0; /* false */
733   dir->info.shift_fraction = 0.0;
734   dir->col                 = 0;
735   dir->row                 = 0;
736 
737   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&dir->solvertype);CHKERRQ(ierr);
738   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
739   if (size == 1) {
740     ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr);
741   } else {
742     ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr);
743   }
744   dir->reusefill        = PETSC_FALSE;
745   dir->reuseordering    = PETSC_FALSE;
746   pc->data              = (void*)dir;
747 
748   pc->ops->destroy           = PCDestroy_LU;
749   pc->ops->apply             = PCApply_LU;
750   pc->ops->applytranspose    = PCApplyTranspose_LU;
751   pc->ops->setup             = PCSetUp_LU;
752   pc->ops->setfromoptions    = PCSetFromOptions_LU;
753   pc->ops->view              = PCView_LU;
754   pc->ops->applyrichardson   = 0;
755   pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU;
756 
757   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_LU",
758                     PCFactorSetMatSolverPackage_LU);CHKERRQ(ierr);
759   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU",
760                     PCFactorSetZeroPivot_LU);CHKERRQ(ierr);
761   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU",
762                     PCFactorSetShiftNonzero_LU);CHKERRQ(ierr);
763   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU",
764                     PCFactorSetShiftPd_LU);CHKERRQ(ierr);
765 
766   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU",
767                     PCFactorSetFill_LU);CHKERRQ(ierr);
768   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU",
769                     PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
770   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU",
771                     PCFactorSetMatOrderingType_LU);CHKERRQ(ierr);
772   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU",
773                     PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
774   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU",
775                     PCFactorSetReuseFill_LU);CHKERRQ(ierr);
776   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU",
777                     PCFactorSetPivoting_LU);CHKERRQ(ierr);
778   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU",
779                     PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr);
780   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU",
781                     PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
782   PetscFunctionReturn(0);
783 }
784 EXTERN_C_END
785