xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 853170213946f745f70941b14c22455e5c54846f) !
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 [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
92    is optional with PETSC_TRUE 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 /*@
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      Use PCFactorGetMatrix(pc,&mat); MatFactorGetPackage(mat,&package); to see what package is being used
303 
304 
305 .keywords: PC, set, factorization, direct, fill
306 
307 .seealso: MatGetFactor(), MatSolverPackage
308 
309 @*/
310 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
311 {
312   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);
313 
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
316   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr);
317   if (f) {
318     ierr = (*f)(pc,stype);CHKERRQ(ierr);
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 #undef __FUNCT__
324 #define __FUNCT__ "PCFactorSetFill"
325 /*@
326    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
327    fill = number nonzeros in factor/number nonzeros in original matrix.
328 
329    Collective on PC
330 
331    Input Parameters:
332 +  pc - the preconditioner context
333 -  fill - amount of expected fill
334 
335    Options Database Key:
336 .  -pc_factor_fill <fill> - Sets fill amount
337 
338    Level: intermediate
339 
340    Note:
341    For sparse matrix factorizations it is difficult to predict how much
342    fill to expect. By running with the option -info PETSc will print the
343    actual amount of fill used; allowing you to set the value accurately for
344    future runs. Default PETSc uses a value of 5.0
345 
346 .keywords: PC, set, factorization, direct, fill
347 
348 @*/
349 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill)
350 {
351   PetscErrorCode ierr,(*f)(PC,PetscReal);
352 
353   PetscFunctionBegin;
354   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
355   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
356   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
357   if (f) {
358     ierr = (*f)(pc,fill);CHKERRQ(ierr);
359   }
360   PetscFunctionReturn(0);
361 }
362 
363 #undef __FUNCT__
364 #define __FUNCT__ "PCFactorSetUseInPlace"
365 /*@
366    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
367    For dense matrices, this enables the solution of much larger problems.
368    For sparse matrices the factorization cannot be done truly in-place
369    so this does not save memory during the factorization, but after the matrix
370    is factored, the original unfactored matrix is freed, thus recovering that
371    space.
372 
373    Collective on PC
374 
375    Input Parameters:
376 .  pc - the preconditioner context
377 
378    Options Database Key:
379 .  -pc_factor_in_place - Activates in-place factorization
380 
381    Notes:
382    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
383    a different matrix is provided for the multiply and the preconditioner in
384    a call to KSPSetOperators().
385    This is because the Krylov space methods require an application of the
386    matrix multiplication, which is not possible here because the matrix has
387    been factored in-place, replacing the original matrix.
388 
389    Level: intermediate
390 
391 .keywords: PC, set, factorization, direct, inplace, in-place, LU
392 
393 .seealso: PCILUSetUseInPlace()
394 @*/
395 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc)
396 {
397   PetscErrorCode ierr,(*f)(PC);
398 
399   PetscFunctionBegin;
400   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
401   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
402   if (f) {
403     ierr = (*f)(pc);CHKERRQ(ierr);
404   }
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCFactorSetMatOrderingType"
410 /*@C
411     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
412     be used in the LU factorization.
413 
414     Collective on PC
415 
416     Input Parameters:
417 +   pc - the preconditioner context
418 -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
419 
420     Options Database Key:
421 .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
422 
423     Level: intermediate
424 
425     Notes: nested dissection is used by default
426 
427     For Cholesky and ICC and the SBAIJ format reorderings are not available,
428     since only the upper triangular part of the matrix is stored. You can use the
429     SeqAIJ format in this case to get reorderings.
430 
431 @*/
432 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
433 {
434   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);
435 
436   PetscFunctionBegin;
437   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr);
438   if (f) {
439     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
440   }
441   PetscFunctionReturn(0);
442 }
443 
444 #undef __FUNCT__
445 #define __FUNCT__ "PCFactorSetPivoting"
446 /*@
447     PCFactorSetPivoting - Determines when pivoting is done during LU.
448       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
449       it is never done. For the Matlab and SuperLU factorization this is used.
450 
451     Collective on PC
452 
453     Input Parameters:
454 +   pc - the preconditioner context
455 -   dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
456 
457     Options Database Key:
458 .   -pc_factor_pivoting <dtcol>
459 
460     Level: intermediate
461 
462 .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
463 @*/
464 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol)
465 {
466   PetscErrorCode ierr,(*f)(PC,PetscReal);
467 
468   PetscFunctionBegin;
469   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr);
470   if (f) {
471     ierr = (*f)(pc,dtcol);CHKERRQ(ierr);
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "PCFactorSetPivotInBlocks"
478 /*@
479     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
480       with BAIJ or SBAIJ matrices
481 
482     Collective on PC
483 
484     Input Parameters:
485 +   pc - the preconditioner context
486 -   pivot - PETSC_TRUE or PETSC_FALSE
487 
488     Options Database Key:
489 .   -pc_factor_pivot_in_blocks <true,false>
490 
491     Level: intermediate
492 
493 .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
494 @*/
495 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
496 {
497   PetscErrorCode ierr,(*f)(PC,PetscTruth);
498 
499   PetscFunctionBegin;
500   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
501   if (f) {
502     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
503   }
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "PCFactorSetReuseFill"
509 /*@
510    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
511    this causes later ones to use the fill ratio computed in the initial factorization.
512 
513    Collective on PC
514 
515    Input Parameters:
516 +  pc - the preconditioner context
517 -  flag - PETSC_TRUE to reuse else PETSC_FALSE
518 
519    Options Database Key:
520 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
521 
522    Level: intermediate
523 
524 .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
525 
526 .seealso: PCFactorSetReuseOrdering()
527 @*/
528 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
529 {
530   PetscErrorCode ierr,(*f)(PC,PetscTruth);
531 
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(pc,PC_COOKIE,2);
534   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
535   if (f) {
536     ierr = (*f)(pc,flag);CHKERRQ(ierr);
537   }
538   PetscFunctionReturn(0);
539 }
540