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