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