xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 9596e0b48258fba4fca4f68feb5185896facfe69)
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()
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()
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()
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__ "PCFactorSetUseDropTolerance"
115 /*@
116    PCFactorSetUseDropTolerance - 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_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
130 
131    Level: intermediate
132 
133     Notes:
134       This uses the iludt() code of Saad's SPARSKIT package
135 
136       There are NO default values for the 3 parameters, you must set them with reasonable values for your
137       matrix. We don't know how to compute reasonable values.
138 
139 .keywords: PC, levels, reordering, factorization, incomplete, ILU
140 @*/
141 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
142 {
143   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
144 
145   PetscFunctionBegin;
146   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
147   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
148   if (f) {
149     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
150   }
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNCT__
155 #define __FUNCT__ "PCFactorSetLevels"
156 /*@
157    PCFactorSetLevels - Sets the number of levels of fill to use.
158 
159    Collective on PC
160 
161    Input Parameters:
162 +  pc - the preconditioner context
163 -  levels - number of levels of fill
164 
165    Options Database Key:
166 .  -pc_factor_levels <levels> - Sets fill level
167 
168    Level: intermediate
169 
170 .keywords: PC, levels, fill, factorization, incomplete, ILU
171 @*/
172 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels)
173 {
174   PetscErrorCode ierr,(*f)(PC,PetscInt);
175 
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
178   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
179   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
180   if (f) {
181     ierr = (*f)(pc,levels);CHKERRQ(ierr);
182   }
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PCFactorSetAllowDiagonalFill"
188 /*@
189    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
190    treated as level 0 fill even if there is no non-zero location.
191 
192    Collective on PC
193 
194    Input Parameters:
195 +  pc - the preconditioner context
196 
197    Options Database Key:
198 .  -pc_factor_diagonal_fill
199 
200    Notes:
201    Does not apply with 0 fill.
202 
203    Level: intermediate
204 
205 .keywords: PC, levels, fill, factorization, incomplete, ILU
206 @*/
207 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc)
208 {
209   PetscErrorCode ierr,(*f)(PC);
210 
211   PetscFunctionBegin;
212   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
213   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
214   if (f) {
215     ierr = (*f)(pc);CHKERRQ(ierr);
216   }
217   PetscFunctionReturn(0);
218 }
219 
220 #undef __FUNCT__
221 #define __FUNCT__ "PCFactorSetShiftInBlocks"
222 /*@
223     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot.
224 
225     Collective on PC
226 
227     Input Parameters:
228 +   pc - the preconditioner context
229 -   shift - amount to shift or PETSC_DECIDE
230 
231     Options Database Key:
232 .   -pc_factor_shift_in_blocks <shift>
233 
234     Level: intermediate
235 
236 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks()
237 @*/
238 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
239 {
240   PetscErrorCode ierr,(*f)(PC,PetscReal);
241 
242   PetscFunctionBegin;
243   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
244   if (f) {
245     ierr = (*f)(pc,shift);CHKERRQ(ierr);
246   }
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    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 PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
271 {
272   PetscErrorCode ierr,(*f)(PC,PetscReal);
273 
274   PetscFunctionBegin;
275   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
276   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
277   if (f) {
278     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
279   }
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "PCFactorSetMatSolverPackage"
285 /*@C
286    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
287 
288    Collective on PC
289 
290    Input Parameters:
291 +  pc - the preconditioner context
292 -  stype - for example, spooles, superlu, superlu_dist
293 
294    Options Database Key:
295 .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps
296 
297    Level: intermediate
298 
299    Note:
300      By default this will use the PETSc factorization if it exists
301 
302 
303 .keywords: PC, set, factorization, direct, fill
304 
305 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
306 
307 @*/
308 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
309 {
310   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
311 
312   PetscFunctionBegin;
313   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
314   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
315   if (f) {
316     ierr = (*f)(pc,stype);CHKERRQ(ierr);
317   }
318   PetscFunctionReturn(0);
319 }
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "PCFactorGetMatSolverPackage"
323 /*@C
324    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
325 
326    Collective on PC
327 
328    Input Parameter:
329 .  pc - the preconditioner context
330 
331    Output Parameter:
332 .   stype - for example, spooles, superlu, superlu_dist
333 
334    Level: intermediate
335 
336 
337 .keywords: PC, set, factorization, direct, fill
338 
339 .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
340 
341 @*/
342 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
343 {
344   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
345 
346   PetscFunctionBegin;
347   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
348   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
349   if (f) {
350     ierr = (*f)(pc,stype);CHKERRQ(ierr);
351   }
352   PetscFunctionReturn(0);
353 }
354 
355 #undef __FUNCT__
356 #define __FUNCT__ "PCFactorSetFill"
357 /*@
358    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
359    fill = number nonzeros in factor/number nonzeros in original matrix.
360 
361    Collective on PC
362 
363    Input Parameters:
364 +  pc - the preconditioner context
365 -  fill - amount of expected fill
366 
367    Options Database Key:
368 .  -pc_factor_fill <fill> - Sets fill amount
369 
370    Level: intermediate
371 
372    Note:
373    For sparse matrix factorizations it is difficult to predict how much
374    fill to expect. By running with the option -info PETSc will print the
375    actual amount of fill used; allowing you to set the value accurately for
376    future runs. Default PETSc uses a value of 5.0
377 
378 .keywords: PC, set, factorization, direct, fill
379 
380 @*/
381 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
382 {
383   PetscErrorCode ierr,(*f)(PC,PetscReal);
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
387   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
388   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
389   if (f) {
390     ierr = (*f)(pc,fill);CHKERRQ(ierr);
391   }
392   PetscFunctionReturn(0);
393 }
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "PCFactorSetUseInPlace"
397 /*@
398    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
399    For dense matrices, this enables the solution of much larger problems.
400    For sparse matrices the factorization cannot be done truly in-place
401    so this does not save memory during the factorization, but after the matrix
402    is factored, the original unfactored matrix is freed, thus recovering that
403    space.
404 
405    Collective on PC
406 
407    Input Parameters:
408 .  pc - the preconditioner context
409 
410    Options Database Key:
411 .  -pc_factor_in_place - Activates in-place factorization
412 
413    Notes:
414    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
415    a different matrix is provided for the multiply and the preconditioner in
416    a call to KSPSetOperators().
417    This is because the Krylov space methods require an application of the
418    matrix multiplication, which is not possible here because the matrix has
419    been factored in-place, replacing the original matrix.
420 
421    Level: intermediate
422 
423 .keywords: PC, set, factorization, direct, inplace, in-place, LU
424 
425 .seealso: PCILUSetUseInPlace()
426 @*/
427 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
428 {
429   PetscErrorCode ierr,(*f)(PC);
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
433   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
434   if (f) {
435     ierr = (*f)(pc);CHKERRQ(ierr);
436   }
437   PetscFunctionReturn(0);
438 }
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "PCFactorSetMatOrderingType"
442 /*@C
443     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
444     be used in the LU factorization.
445 
446     Collective on PC
447 
448     Input Parameters:
449 +   pc - the preconditioner context
450 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
451 
452     Options Database Key:
453 .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
454 
455     Level: intermediate
456 
457     Notes: nested dissection is used by default
458 
459     For Cholesky and ICC and the SBAIJ format reorderings are not available,
460     since only the upper triangular part of the matrix is stored. You can use the
461     SeqAIJ format in this case to get reorderings.
462 
463 @*/
464 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
465 {
466   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
467 
468   PetscFunctionBegin;
469   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
470   if (f) {
471     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "PCFactorSetPivoting"
478 /*@
479     PCFactorSetPivoting - Determines when pivoting is done during LU.
480       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
481       it is never done. For the Matlab and SuperLU factorization this is used.
482 
483     Collective on PC
484 
485     Input Parameters:
486 +   pc - the preconditioner context
487 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
488 
489     Options Database Key:
490 .   -pc_factor_pivoting <dtcol>
491 
492     Level: intermediate
493 
494 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
495 @*/
496 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
497 {
498   PetscErrorCode ierr,(*f)(PC,PetscReal);
499 
500   PetscFunctionBegin;
501   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
502   if (f) {
503     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
504   }
505   PetscFunctionReturn(0);
506 }
507 
508 #undef __FUNCT__
509 #define __FUNCT__ "PCFactorSetPivotInBlocks"
510 /*@
511     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
512       with BAIJ or SBAIJ matrices
513 
514     Collective on PC
515 
516     Input Parameters:
517 +   pc - the preconditioner context
518 -   pivot - PETSC_TRUE or PETSC_FALSE
519 
520     Options Database Key:
521 .   -pc_factor_pivot_in_blocks <true,false>
522 
523     Level: intermediate
524 
525 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
526 @*/
527 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
528 {
529   PetscErrorCode ierr,(*f)(PC,PetscTruth);
530 
531   PetscFunctionBegin;
532   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
533   if (f) {
534     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
535   }
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "PCFactorSetReuseFill"
541 /*@
542    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
543    this causes later ones to use the fill ratio computed in the initial factorization.
544 
545    Collective on PC
546 
547    Input Parameters:
548 +  pc - the preconditioner context
549 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
550 
551    Options Database Key:
552 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
553 
554    Level: intermediate
555 
556 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
557 
558 .seealso: PCFactorSetReuseOrdering()
559 @*/
560 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
561 {
562   PetscErrorCode ierr,(*f)(PC,PetscTruth);
563 
564   PetscFunctionBegin;
565   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
566   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
567   if (f) {
568     ierr = (*f)(pc,flag);CHKERRQ(ierr);
569   }
570   PetscFunctionReturn(0);
571 }
572