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