xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 7afc1a8b0cbfcb157f2be8aa81c6cdb9b187c271)
1 
2 #include <../src/ksp/pc/impls/factor/factor.h>  /*I "petscpc.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PCFactorSetUpMatSolverPackage"
6 /*@
7     PCFactorSetUpMatSolverPackage - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
8        set the options for that particular factorization object.
9 
10   Input Parameter:
11 .  pc  - the preconditioner context
12 
13   Notes: After you have called this function (which has to be after the KSPSetOperators() or PCSetOperators()) you can call PCFactorGetMatrix() and then set factor options on that matrix.
14 
15 .seealso: PCFactorSetMatSolverPackage(), PCFactorGetMatrix()
16 
17   Level: intermediate
18 
19 @*/
20 PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
26   ierr = PetscTryMethod(pc,"PCFactorSetUpMatSolverPackage_C",(PC),(pc));CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "PCFactorSetZeroPivot"
32 /*@
33    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
34 
35    Logically Collective on PC
36 
37    Input Parameters:
38 +  pc - the preconditioner context
39 -  zero - all pivots smaller than this will be considered zero
40 
41    Options Database Key:
42 .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
43 
44    Level: intermediate
45 
46 .keywords: PC, set, factorization, direct, fill
47 
48 .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
49 @*/
50 PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
51 {
52   PetscErrorCode ierr;
53 
54   PetscFunctionBegin;
55   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
56   PetscValidLogicalCollectiveReal(pc,zero,2);
57   ierr = PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));CHKERRQ(ierr);
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "PCFactorSetShiftType"
63 /*@
64    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
65      numerical factorization, thus the matrix has nonzero pivots
66 
67    Logically Collective on PC
68 
69    Input Parameters:
70 +  pc - the preconditioner context
71 -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
72 
73    Options Database Key:
74 .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
75 
76    Level: intermediate
77 
78 .keywords: PC, set, factorization,
79 
80 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
81 @*/
82 PetscErrorCode  PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
83 {
84   PetscErrorCode ierr;
85 
86   PetscFunctionBegin;
87   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
88   PetscValidLogicalCollectiveEnum(pc,shifttype,2);
89   ierr = PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "PCFactorSetShiftAmount"
95 /*@
96    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
97      numerical factorization, thus the matrix has nonzero pivots
98 
99    Logically Collective on PC
100 
101    Input Parameters:
102 +  pc - the preconditioner context
103 -  shiftamount - amount of shift
104 
105    Options Database Key:
106 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
107 
108    Level: intermediate
109 
110 .keywords: PC, set, factorization,
111 
112 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
113 @*/
114 PetscErrorCode  PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
115 {
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
120   PetscValidLogicalCollectiveReal(pc,shiftamount,2);
121   ierr = PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PCFactorSetDropTolerance"
127 /*
128    PCFactorSetDropTolerance - The preconditioner will use an ILU
129    based on a drop tolerance. (Under development)
130 
131    Logically Collective on PC
132 
133    Input Parameters:
134 +  pc - the preconditioner context
135 .  dt - the drop tolerance, try from 1.e-10 to .1
136 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
137 -  maxrowcount - the max number of nonzeros allowed in a row, best value
138                  depends on the number of nonzeros in row of original matrix
139 
140    Options Database Key:
141 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
142 
143    Level: intermediate
144 
145       There are NO default values for the 3 parameters, you must set them with reasonable values for your
146       matrix. We don't know how to compute reasonable values.
147 
148 .keywords: PC, levels, reordering, factorization, incomplete, ILU
149 */
150 PetscErrorCode  PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
151 {
152   PetscErrorCode ierr;
153 
154   PetscFunctionBegin;
155   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
156   PetscValidLogicalCollectiveReal(pc,dtcol,2);
157   PetscValidLogicalCollectiveInt(pc,maxrowcount,3);
158   ierr = PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "PCFactorGetLevels"
164 /*@
165    PCFactorGetLevels - Gets the number of levels of fill to use.
166 
167    Logically Collective on PC
168 
169    Input Parameters:
170 .  pc - the preconditioner context
171 
172    Output Parameter:
173 .  levels - number of levels of fill
174 
175    Level: intermediate
176 
177 .keywords: PC, levels, fill, factorization, incomplete, ILU
178 @*/
179 PetscErrorCode  PCFactorGetLevels(PC pc,PetscInt *levels)
180 {
181   PetscErrorCode ierr;
182 
183   PetscFunctionBegin;
184   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
185   ierr = PetscTryMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "PCFactorSetLevels"
191 /*@
192    PCFactorSetLevels - Sets the number of levels of fill to use.
193 
194    Logically Collective on PC
195 
196    Input Parameters:
197 +  pc - the preconditioner context
198 -  levels - number of levels of fill
199 
200    Options Database Key:
201 .  -pc_factor_levels <levels> - Sets fill level
202 
203    Level: intermediate
204 
205 .keywords: PC, levels, fill, factorization, incomplete, ILU
206 @*/
207 PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
208 {
209   PetscErrorCode ierr;
210 
211   PetscFunctionBegin;
212   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
213   if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
214   PetscValidLogicalCollectiveInt(pc,levels,2);
215   ierr = PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
221 /*@
222    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
223    treated as level 0 fill even if there is no non-zero location.
224 
225    Logically Collective on PC
226 
227    Input Parameters:
228 +  pc - the preconditioner context
229 
230    Options Database Key:
231 .  -pc_factor_diagonal_fill
232 
233    Notes:
234    Does not apply with 0 fill.
235 
236    Level: intermediate
237 
238 .keywords: PC, levels, fill, factorization, incomplete, ILU
239 @*/
240 PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc)
241 {
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
246   ierr = PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC),(pc));CHKERRQ(ierr);
247   PetscFunctionReturn(0);
248 }
249 
250 #undef __FUNCT__
251 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
252 /*@
253    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
254 
255    Logically Collective on PC
256 
257    Input Parameters:
258 +  pc - the preconditioner context
259 -  tol - diagonal entries smaller than this in absolute value are considered zero
260 
261    Options Database Key:
262 .  -pc_factor_nonzeros_along_diagonal
263 
264    Level: intermediate
265 
266 .keywords: PC, set, factorization, direct, fill
267 
268 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
269 @*/
270 PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
271 {
272   PetscErrorCode ierr;
273 
274   PetscFunctionBegin;
275   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
276   PetscValidLogicalCollectiveReal(pc,rtol,2);
277   ierr = PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "PCFactorSetMatSolverPackage"
283 /*@C
284    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
285 
286    Logically Collective on PC
287 
288    Input Parameters:
289 +  pc - the preconditioner context
290 -  stype - for example, superlu, superlu_dist
291 
292    Options Database Key:
293 .  -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse
294 
295    Level: intermediate
296 
297    Note:
298      By default this will use the PETSc factorization if it exists
299 
300 
301 .keywords: PC, set, factorization, direct, fill
302 
303 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
304 
305 @*/
306 PetscErrorCode  PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
307 {
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
312   ierr = PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));CHKERRQ(ierr);
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "PCFactorGetMatSolverPackage"
318 /*@C
319    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
320 
321    Not Collective
322 
323    Input Parameter:
324 .  pc - the preconditioner context
325 
326    Output Parameter:
327 .   stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
328 
329    Level: intermediate
330 
331 
332 .keywords: PC, set, factorization, direct, fill
333 
334 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
335 
336 @*/
337 PetscErrorCode  PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
338 {
339   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
340 
341   PetscFunctionBegin;
342   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
343   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);CHKERRQ(ierr);
344   if (f) {
345     ierr = (*f)(pc,stype);CHKERRQ(ierr);
346   } else {
347     *stype = NULL;
348   }
349   PetscFunctionReturn(0);
350 }
351 
352 #undef __FUNCT__
353 #define __FUNCT__ "PCFactorSetFill"
354 /*@
355    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
356    fill = number nonzeros in factor/number nonzeros in original matrix.
357 
358    Not Collective, each process can expect a different amount of fill
359 
360    Input Parameters:
361 +  pc - the preconditioner context
362 -  fill - amount of expected fill
363 
364    Options Database Key:
365 .  -pc_factor_fill <fill> - Sets fill amount
366 
367    Level: intermediate
368 
369    Note:
370    For sparse matrix factorizations it is difficult to predict how much
371    fill to expect. By running with the option -info PETSc will print the
372    actual amount of fill used; allowing you to set the value accurately for
373    future runs. Default PETSc uses a value of 5.0
374 
375    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
376 
377 
378 .keywords: PC, set, factorization, direct, fill
379 
380 @*/
381 PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387   if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
388   ierr = PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "PCFactorSetUseInPlace"
394 /*@
395    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
396    For dense matrices, this enables the solution of much larger problems.
397    For sparse matrices the factorization cannot be done truly in-place
398    so this does not save memory during the factorization, but after the matrix
399    is factored, the original unfactored matrix is freed, thus recovering that
400    space.
401 
402    Logically Collective on PC
403 
404    Input Parameters:
405 .  pc - the preconditioner context
406 
407    Options Database Key:
408 .  -pc_factor_in_place - Activates in-place factorization
409 
410    Notes:
411    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
412    a different matrix is provided for the multiply and the preconditioner in
413    a call to KSPSetOperators().
414    This is because the Krylov space methods require an application of the
415    matrix multiplication, which is not possible here because the matrix has
416    been factored in-place, replacing the original matrix.
417 
418    Level: intermediate
419 
420 .keywords: PC, set, factorization, direct, inplace, in-place, LU
421 
422 .seealso: PCILUSetUseInPlace()
423 @*/
424 PetscErrorCode  PCFactorSetUseInPlace(PC pc)
425 {
426   PetscErrorCode ierr;
427 
428   PetscFunctionBegin;
429   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
430   ierr = PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC),(pc));CHKERRQ(ierr);
431   PetscFunctionReturn(0);
432 }
433 
434 #undef __FUNCT__
435 #define __FUNCT__ "PCFactorSetMatOrderingType"
436 /*@C
437     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
438     be used in the LU factorization.
439 
440     Logically Collective on PC
441 
442     Input Parameters:
443 +   pc - the preconditioner context
444 -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
445 
446     Options Database Key:
447 .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
448 
449     Level: intermediate
450 
451     Notes: nested dissection is used by default
452 
453     For Cholesky and ICC and the SBAIJ format reorderings are not available,
454     since only the upper triangular part of the matrix is stored. You can use the
455     SeqAIJ format in this case to get reorderings.
456 
457 @*/
458 PetscErrorCode  PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
459 {
460   PetscErrorCode ierr;
461 
462   PetscFunctionBegin;
463   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
464   ierr = PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 #undef __FUNCT__
469 #define __FUNCT__ "PCFactorSetColumnPivot"
470 /*@
471     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
472       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
473       it is never done. For the MATLAB and SuperLU factorization this is used.
474 
475     Logically Collective on PC
476 
477     Input Parameters:
478 +   pc - the preconditioner context
479 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
480 
481     Options Database Key:
482 .   -pc_factor_pivoting <dtcol>
483 
484     Level: intermediate
485 
486 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
487 @*/
488 PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
489 {
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
494   PetscValidLogicalCollectiveReal(pc,dtcol,2);
495   ierr = PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));CHKERRQ(ierr);
496   PetscFunctionReturn(0);
497 }
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "PCFactorSetPivotInBlocks"
501 /*@
502     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
503       with BAIJ or SBAIJ matrices
504 
505     Logically Collective on PC
506 
507     Input Parameters:
508 +   pc - the preconditioner context
509 -   pivot - PETSC_TRUE or PETSC_FALSE
510 
511     Options Database Key:
512 .   -pc_factor_pivot_in_blocks <true,false>
513 
514     Level: intermediate
515 
516 .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
517 @*/
518 PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
519 {
520   PetscErrorCode ierr;
521 
522   PetscFunctionBegin;
523   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
524   PetscValidLogicalCollectiveBool(pc,pivot,2);
525   ierr = PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));CHKERRQ(ierr);
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "PCFactorSetReuseFill"
531 /*@
532    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
533    this causes later ones to use the fill ratio computed in the initial factorization.
534 
535    Logically Collective on PC
536 
537    Input Parameters:
538 +  pc - the preconditioner context
539 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
540 
541    Options Database Key:
542 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
543 
544    Level: intermediate
545 
546 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
547 
548 .seealso: PCFactorSetReuseOrdering()
549 @*/
550 PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool flag)
551 {
552   PetscErrorCode ierr;
553 
554   PetscFunctionBegin;
555   PetscValidHeaderSpecific(pc,PC_CLASSID,2);
556   PetscValidLogicalCollectiveBool(pc,flag,2);
557   ierr = PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr);
558   PetscFunctionReturn(0);
559 }
560