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