xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 289f3aa11b3ac459c99a6895ea3b9e0a0e6ded3a)
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: PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks()
24 @*/
25 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero)
26 {
27   PetscErrorCode ierr,(*f)(PC,PetscReal);
28 
29   PetscFunctionBegin;
30   PetscValidHeaderSpecific(pc,PC_COOKIE,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__ "PCFactorSetShiftNonzero"
40 /*@
41    PCFactorSetShiftNonzero - adds this 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 -  shift - amount of shift
49 
50    Options Database Key:
51 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
52 
53    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
54          pivot, then the shift is doubled until this is alleviated.
55 
56    Level: intermediate
57 
58 .keywords: PC, set, factorization, direct, fill
59 
60 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd(), PCFactorSetShiftInBlocks()
61 @*/
62 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift)
63 {
64   PetscErrorCode ierr,(*f)(PC,PetscReal);
65 
66   PetscFunctionBegin;
67   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
68   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr);
69   if (f) {
70     ierr = (*f)(pc,shift);CHKERRQ(ierr);
71   }
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "PCFactorSetShiftPd"
77 /*@
78    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
79    If a matrix factorisation breaks down because of nonpositive pivots,
80    adding sufficient identity to the diagonal will remedy this.
81    Setting this causes a bisection method to find the minimum shift that
82    will lead to a well-defined matrix factor.
83 
84    Collective on PC
85 
86    Input parameters:
87 +  pc - the preconditioner context
88 -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
89 
90    Options Database Key:
91 .  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
92    is optional with false being the default
93 
94    Level: intermediate
95 
96 .keywords: PC, indefinite, factorization
97 
98 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftInBlocks()
99 @*/
100 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift)
101 {
102   PetscErrorCode ierr,(*f)(PC,PetscTruth);
103 
104   PetscFunctionBegin;
105   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
106   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr);
107   if (f) {
108     ierr = (*f)(pc,shift);CHKERRQ(ierr);
109   }
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNCT__
114 #define __FUNCT__ "PCFactorSetDropTolerance"
115 /*@
116    PCFactorSetDropTolerance - The preconditioner will use an ILU
117    based on a drop tolerance.
118 
119    Collective on PC
120 
121    Input Parameters:
122 +  pc - the preconditioner context
123 .  dt - the drop tolerance, try from 1.e-10 to .1
124 .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
125 -  maxrowcount - the max number of nonzeros allowed in a row, best value
126                  depends on the number of nonzeros in row of original matrix
127 
128    Options Database Key:
129 .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
130 
131    Level: intermediate
132 
133       There are NO default values for the 3 parameters, you must set them with reasonable values for your
134       matrix. We don't know how to compute reasonable values.
135 
136 .keywords: PC, levels, reordering, factorization, incomplete, ILU
137 @*/
138 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
139 {
140   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
141 
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
144   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
145   if (f) {
146     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
147   }
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "PCFactorSetLevels"
153 /*@
154    PCFactorSetLevels - Sets the number of levels of fill to use.
155 
156    Collective on PC
157 
158    Input Parameters:
159 +  pc - the preconditioner context
160 -  levels - number of levels of fill
161 
162    Options Database Key:
163 .  -pc_factor_levels <levels> - Sets fill level
164 
165    Level: intermediate
166 
167 .keywords: PC, levels, fill, factorization, incomplete, ILU
168 @*/
169 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
170 {
171   PetscErrorCode ierr,(*f)(PC,PetscInt);
172 
173   PetscFunctionBegin;
174   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
175   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
176   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
177   if (f) {
178     ierr = (*f)(pc,levels);CHKERRQ(ierr);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
185 /*@
186    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
187    treated as level 0 fill even if there is no non-zero location.
188 
189    Collective on PC
190 
191    Input Parameters:
192 +  pc - the preconditioner context
193 
194    Options Database Key:
195 .  -pc_factor_diagonal_fill
196 
197    Notes:
198    Does not apply with 0 fill.
199 
200    Level: intermediate
201 
202 .keywords: PC, levels, fill, factorization, incomplete, ILU
203 @*/
204 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
205 {
206   PetscErrorCode ierr,(*f)(PC);
207 
208   PetscFunctionBegin;
209   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
210   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
211   if (f) {
212     ierr = (*f)(pc);CHKERRQ(ierr);
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "PCFactorSetShiftInBlocks"
219 /*@
220     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if matrix factor on that block detects a zero pivot.
221 
222     Collective on PC
223 
224     Input Parameters:
225 +   pc - the preconditioner context
226 -   shift - amount to shift or PETSC_DECIDE
227 
228     Options Database Key:
229 .   -pc_factor_shift_in_blocks <shift>
230 
231     Level: intermediate
232 
233 .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot(), PCFactorSetShiftInBlocks()
234 @*/
235 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
236 {
237   PetscErrorCode ierr,(*f)(PC,PetscReal);
238 
239   PetscFunctionBegin;
240   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
241   if (f) {
242     ierr = (*f)(pc,shift);CHKERRQ(ierr);
243   }
244   PetscFunctionReturn(0);
245 }
246 
247 #undef __FUNCT__
248 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal"
249 /*@
250    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
251 
252    Collective on PC
253 
254    Input Parameters:
255 +  pc - the preconditioner context
256 -  tol - diagonal entries smaller than this in absolute value are considered zero
257 
258    Options Database Key:
259 .  -pc_factor_nonzeros_along_diagonal
260 
261    Level: intermediate
262 
263 .keywords: PC, set, factorization, direct, fill
264 
265 .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
266 @*/
267 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
268 {
269   PetscErrorCode ierr,(*f)(PC,PetscReal);
270 
271   PetscFunctionBegin;
272   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
273   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
274   if (f) {
275     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
276   }
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "PCFactorSetMatSolverPackage"
282 /*@C
283    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
284 
285    Collective on PC
286 
287    Input Parameters:
288 +  pc - the preconditioner context
289 -  stype - for example, spooles, superlu, superlu_dist
290 
291    Options Database Key:
292 .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
293 
294    Level: intermediate
295 
296    Note:
297      By default this will use the PETSc factorization if it exists
298 
299 
300 .keywords: PC, set, factorization, direct, fill
301 
302 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
303 
304 @*/
305 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
306 {
307   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
311   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
312   if (f) {
313     ierr = (*f)(pc,stype);CHKERRQ(ierr);
314   }
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "PCFactorGetMatSolverPackage"
320 /*@C
321    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
322 
323    Collective on PC
324 
325    Input Parameter:
326 .  pc - the preconditioner context
327 
328    Output Parameter:
329 .   stype - for example, spooles, superlu, superlu_dist
330 
331    Level: intermediate
332 
333 
334 .keywords: PC, set, factorization, direct, fill
335 
336 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
337 
338 @*/
339 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
340 {
341   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
342 
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
345   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
346   if (f) {
347     ierr = (*f)(pc,stype);CHKERRQ(ierr);
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    Collective on PC
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 .keywords: PC, set, factorization, direct, fill
376 
377 @*/
378 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
379 {
380   PetscErrorCode ierr,(*f)(PC,PetscReal);
381 
382   PetscFunctionBegin;
383   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
384   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
385   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
386   if (f) {
387     ierr = (*f)(pc,fill);CHKERRQ(ierr);
388   }
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    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 PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
425 {
426   PetscErrorCode ierr,(*f)(PC);
427 
428   PetscFunctionBegin;
429   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
430   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
431   if (f) {
432     ierr = (*f)(pc);CHKERRQ(ierr);
433   }
434   PetscFunctionReturn(0);
435 }
436 
437 #undef __FUNCT__
438 #define __FUNCT__ "PCFactorSetMatOrderingType"
439 /*@C
440     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
441     be used in the LU factorization.
442 
443     Collective on PC
444 
445     Input Parameters:
446 +   pc - the preconditioner context
447 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
448 
449     Options Database Key:
450 .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
451 
452     Level: intermediate
453 
454     Notes: nested dissection is used by default
455 
456     For Cholesky and ICC and the SBAIJ format reorderings are not available,
457     since only the upper triangular part of the matrix is stored. You can use the
458     SeqAIJ format in this case to get reorderings.
459 
460 @*/
461 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
462 {
463   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
464 
465   PetscFunctionBegin;
466   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
467   if (f) {
468     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
469   }
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "PCFactorSetColumnPivot"
475 /*@
476     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
477       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
478       it is never done. For the Matlab and SuperLU factorization this is used.
479 
480     Collective on PC
481 
482     Input Parameters:
483 +   pc - the preconditioner context
484 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
485 
486     Options Database Key:
487 .   -pc_factor_pivoting <dtcol>
488 
489     Level: intermediate
490 
491 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
492 @*/
493 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
494 {
495   PetscErrorCode ierr,(*f)(PC,PetscReal);
496 
497   PetscFunctionBegin;
498   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
499   if (f) {
500     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
501   }
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "PCFactorSetPivotInBlocks"
507 /*@
508     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
509       with BAIJ or SBAIJ matrices
510 
511     Collective on PC
512 
513     Input Parameters:
514 +   pc - the preconditioner context
515 -   pivot - PETSC_TRUE or PETSC_FALSE
516 
517     Options Database Key:
518 .   -pc_factor_pivot_in_blocks <true,false>
519 
520     Level: intermediate
521 
522 .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
523 @*/
524 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
525 {
526   PetscErrorCode ierr,(*f)(PC,PetscTruth);
527 
528   PetscFunctionBegin;
529   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
530   if (f) {
531     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
532   }
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "PCFactorSetReuseFill"
538 /*@
539    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
540    this causes later ones to use the fill ratio computed in the initial factorization.
541 
542    Collective on PC
543 
544    Input Parameters:
545 +  pc - the preconditioner context
546 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
547 
548    Options Database Key:
549 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
550 
551    Level: intermediate
552 
553 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
554 
555 .seealso: PCFactorSetReuseOrdering()
556 @*/
557 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
558 {
559   PetscErrorCode ierr,(*f)(PC,PetscTruth);
560 
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
563   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
564   if (f) {
565     ierr = (*f)(pc,flag);CHKERRQ(ierr);
566   }
567   PetscFunctionReturn(0);
568 }
569