xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision 63dd3a1af947c5d07342e150a17d503f381a847e)
1 #define PETSCKSP_DLL
2 
3 #include "src/ksp/pc/pcimpl.h"                /*I "petscpc.h" I*/
4 
5 /*  Options Database Keys: ???
6 .  -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots
7 .  -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner
8 .  -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot
9  */
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "PCFactorSetZeroPivot"
13 /*@
14    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
15 
16    Collective on PC
17 
18    Input Parameters:
19 +  pc - the preconditioner context
20 -  zero - all pivots smaller than this will be considered zero
21 
22    Options Database Key:
23 .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
24 
25    Level: intermediate
26 
27 .keywords: PC, set, factorization, direct, fill
28 
29 .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd()
30 @*/
31 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero)
32 {
33   PetscErrorCode ierr,(*f)(PC,PetscReal);
34 
35   PetscFunctionBegin;
36   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
37   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
38   if (f) {
39     ierr = (*f)(pc,zero);CHKERRQ(ierr);
40   }
41   PetscFunctionReturn(0);
42 }
43 
44 #undef __FUNCT__
45 #define __FUNCT__ "PCFactorSetShiftNonzero"
46 /*@
47    PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during
48      numerical factorization, thus the matrix has nonzero pivots
49 
50    Collective on PC
51 
52    Input Parameters:
53 +  pc - the preconditioner context
54 -  shift - amount of shift
55 
56    Options Database Key:
57 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
58 
59    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
60          pivot, then the shift is doubled until this is alleviated.
61 
62    Level: intermediate
63 
64 .keywords: PC, set, factorization, direct, fill
65 
66 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd()
67 @*/
68 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift)
69 {
70   PetscErrorCode ierr,(*f)(PC,PetscReal);
71 
72   PetscFunctionBegin;
73   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
74   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr);
75   if (f) {
76     ierr = (*f)(pc,shift);CHKERRQ(ierr);
77   }
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "PCFactorSetShiftPd"
83 /*@
84    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
85    If a matrix factorisation breaks down because of nonpositive pivots,
86    adding sufficient identity to the diagonal will remedy this.
87    Setting this causes a bisection method to find the minimum shift that
88    will lead to a well-defined matrix factor.
89 
90    Collective on PC
91 
92    Input parameters:
93 +  pc - the preconditioner context
94 -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
95 
96    Options Database Key:
97 .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
98    is optional with PETSC_TRUE being the default
99 
100    Level: intermediate
101 
102 .keywords: PC, indefinite, factorization
103 
104 .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero()
105 @*/
106 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift)
107 {
108   PetscErrorCode ierr,(*f)(PC,PetscTruth);
109 
110   PetscFunctionBegin;
111   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
112   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr);
113   if (f) {
114     ierr = (*f)(pc,shift);CHKERRQ(ierr);
115   }
116   PetscFunctionReturn(0);
117 }
118