xref: /petsc/src/ksp/pc/impls/factor/factor.c (revision e1222f9f1d0bf3ef979e3eaa00088ee6427d6243)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
25e8efad8SHong Zhang 
3*e1222f9fSBarry Smith #include "src/ksp/pc/impls/factor/factor.h"  /*I "petscpc.h" I*/
45e8efad8SHong Zhang 
5ee45ca4aSHong Zhang #undef __FUNCT__
6ee45ca4aSHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot"
7ee45ca4aSHong Zhang /*@
8ee45ca4aSHong Zhang    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
9ee45ca4aSHong Zhang 
10ee45ca4aSHong Zhang    Collective on PC
11ee45ca4aSHong Zhang 
12ee45ca4aSHong Zhang    Input Parameters:
13afaefe49SHong Zhang +  pc - the preconditioner context
14afaefe49SHong Zhang -  zero - all pivots smaller than this will be considered zero
15ee45ca4aSHong Zhang 
16ee45ca4aSHong Zhang    Options Database Key:
17ee45ca4aSHong Zhang .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
18ee45ca4aSHong Zhang 
19ee45ca4aSHong Zhang    Level: intermediate
20ee45ca4aSHong Zhang 
21ee45ca4aSHong Zhang .keywords: PC, set, factorization, direct, fill
22ee45ca4aSHong Zhang 
23ee45ca4aSHong Zhang .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd()
24ee45ca4aSHong Zhang @*/
25dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero)
26ee45ca4aSHong Zhang {
27afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
28afaefe49SHong Zhang 
29ee45ca4aSHong Zhang   PetscFunctionBegin;
30afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
31afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
32afaefe49SHong Zhang   if (f) {
33afaefe49SHong Zhang     ierr = (*f)(pc,zero);CHKERRQ(ierr);
34afaefe49SHong Zhang   }
35ee45ca4aSHong Zhang   PetscFunctionReturn(0);
36ee45ca4aSHong Zhang }
37ee45ca4aSHong Zhang 
385e8efad8SHong Zhang #undef __FUNCT__
395e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero"
405e8efad8SHong Zhang /*@
415e8efad8SHong Zhang    PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during
425e8efad8SHong Zhang      numerical factorization, thus the matrix has nonzero pivots
435e8efad8SHong Zhang 
445e8efad8SHong Zhang    Collective on PC
455e8efad8SHong Zhang 
465e8efad8SHong Zhang    Input Parameters:
47afaefe49SHong Zhang +  pc - the preconditioner context
48afaefe49SHong Zhang -  shift - amount of shift
495e8efad8SHong Zhang 
505e8efad8SHong Zhang    Options Database Key:
519f95998fSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
525e8efad8SHong Zhang 
535e8efad8SHong Zhang    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
545e8efad8SHong Zhang          pivot, then the shift is doubled until this is alleviated.
555e8efad8SHong Zhang 
565e8efad8SHong Zhang    Level: intermediate
575e8efad8SHong Zhang 
585e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill
595e8efad8SHong Zhang 
60ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd()
615e8efad8SHong Zhang @*/
62dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift)
635e8efad8SHong Zhang {
64afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
65afaefe49SHong Zhang 
665e8efad8SHong Zhang   PetscFunctionBegin;
67afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
68afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr);
69afaefe49SHong Zhang   if (f) {
70afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
715e8efad8SHong Zhang   }
725e8efad8SHong Zhang   PetscFunctionReturn(0);
735e8efad8SHong Zhang }
740a29876aSHong Zhang 
750a29876aSHong Zhang #undef __FUNCT__
760a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd"
770a29876aSHong Zhang /*@
780a29876aSHong Zhang    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
790a29876aSHong Zhang    If a matrix factorisation breaks down because of nonpositive pivots,
800a29876aSHong Zhang    adding sufficient identity to the diagonal will remedy this.
810a29876aSHong Zhang    Setting this causes a bisection method to find the minimum shift that
820a29876aSHong Zhang    will lead to a well-defined matrix factor.
830a29876aSHong Zhang 
840a29876aSHong Zhang    Collective on PC
850a29876aSHong Zhang 
860a29876aSHong Zhang    Input parameters:
87afaefe49SHong Zhang +  pc - the preconditioner context
88afaefe49SHong Zhang -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
890a29876aSHong Zhang 
900a29876aSHong Zhang    Options Database Key:
919f95998fSHong Zhang .  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
929f95998fSHong Zhang    is optional with PETSC_TRUE being the default
930a29876aSHong Zhang 
940a29876aSHong Zhang    Level: intermediate
950a29876aSHong Zhang 
960a29876aSHong Zhang .keywords: PC, indefinite, factorization
970a29876aSHong Zhang 
98ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero()
990a29876aSHong Zhang @*/
100dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift)
1010a29876aSHong Zhang {
102afaefe49SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
103afaefe49SHong Zhang 
1040a29876aSHong Zhang   PetscFunctionBegin;
105afaefe49SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
106afaefe49SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr);
107afaefe49SHong Zhang   if (f) {
108afaefe49SHong Zhang     ierr = (*f)(pc,shift);CHKERRQ(ierr);
109afaefe49SHong Zhang   }
1100a29876aSHong Zhang   PetscFunctionReturn(0);
1110a29876aSHong Zhang }
112