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