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