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