1*dba47a55SKris Buschelman #define PETSCKSP_DLL 25e8efad8SHong Zhang 35e8efad8SHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 45e8efad8SHong Zhang 5ee45ca4aSHong Zhang /* Options Database Keys: ??? 6ee45ca4aSHong Zhang . -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots 7ee45ca4aSHong Zhang . -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner 8ee45ca4aSHong Zhang . -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot 9ee45ca4aSHong Zhang */ 10ee45ca4aSHong Zhang 11ee45ca4aSHong Zhang #undef __FUNCT__ 12ee45ca4aSHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot" 13ee45ca4aSHong Zhang /*@ 14ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 15ee45ca4aSHong Zhang 16ee45ca4aSHong Zhang Collective on PC 17ee45ca4aSHong Zhang 18ee45ca4aSHong Zhang Input Parameters: 19afaefe49SHong Zhang + pc - the preconditioner context 20afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 21ee45ca4aSHong Zhang 22ee45ca4aSHong Zhang Options Database Key: 23ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 24ee45ca4aSHong Zhang 25ee45ca4aSHong Zhang Level: intermediate 26ee45ca4aSHong Zhang 27ee45ca4aSHong Zhang .keywords: PC, set, factorization, direct, fill 28ee45ca4aSHong Zhang 29ee45ca4aSHong Zhang .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd() 30ee45ca4aSHong Zhang @*/ 31*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot(PC pc,PetscReal zero) 32ee45ca4aSHong Zhang { 33afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 34afaefe49SHong Zhang 35ee45ca4aSHong Zhang PetscFunctionBegin; 36afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 37afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 38afaefe49SHong Zhang if (f) { 39afaefe49SHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 40afaefe49SHong Zhang } 41ee45ca4aSHong Zhang PetscFunctionReturn(0); 42ee45ca4aSHong Zhang } 43ee45ca4aSHong Zhang 445e8efad8SHong Zhang #undef __FUNCT__ 455e8efad8SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero" 465e8efad8SHong Zhang /*@ 475e8efad8SHong Zhang PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during 485e8efad8SHong Zhang numerical factorization, thus the matrix has nonzero pivots 495e8efad8SHong Zhang 505e8efad8SHong Zhang Collective on PC 515e8efad8SHong Zhang 525e8efad8SHong Zhang Input Parameters: 53afaefe49SHong Zhang + pc - the preconditioner context 54afaefe49SHong Zhang - shift - amount of shift 555e8efad8SHong Zhang 565e8efad8SHong Zhang Options Database Key: 579f95998fSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 585e8efad8SHong Zhang 595e8efad8SHong Zhang Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero 605e8efad8SHong Zhang pivot, then the shift is doubled until this is alleviated. 615e8efad8SHong Zhang 625e8efad8SHong Zhang Level: intermediate 635e8efad8SHong Zhang 645e8efad8SHong Zhang .keywords: PC, set, factorization, direct, fill 655e8efad8SHong Zhang 66ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd() 675e8efad8SHong Zhang @*/ 68*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero(PC pc,PetscReal shift) 695e8efad8SHong Zhang { 70afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 71afaefe49SHong Zhang 725e8efad8SHong Zhang PetscFunctionBegin; 73afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 74afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);CHKERRQ(ierr); 75afaefe49SHong Zhang if (f) { 76afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 775e8efad8SHong Zhang } 785e8efad8SHong Zhang PetscFunctionReturn(0); 795e8efad8SHong Zhang } 800a29876aSHong Zhang 810a29876aSHong Zhang #undef __FUNCT__ 820a29876aSHong Zhang #define __FUNCT__ "PCFactorSetShiftPd" 830a29876aSHong Zhang /*@ 840a29876aSHong Zhang PCFactorSetShiftPd - specify whether to use Manteuffel shifting. 850a29876aSHong Zhang If a matrix factorisation breaks down because of nonpositive pivots, 860a29876aSHong Zhang adding sufficient identity to the diagonal will remedy this. 870a29876aSHong Zhang Setting this causes a bisection method to find the minimum shift that 880a29876aSHong Zhang will lead to a well-defined matrix factor. 890a29876aSHong Zhang 900a29876aSHong Zhang Collective on PC 910a29876aSHong Zhang 920a29876aSHong Zhang Input parameters: 93afaefe49SHong Zhang + pc - the preconditioner context 94afaefe49SHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 950a29876aSHong Zhang 960a29876aSHong Zhang Options Database Key: 979f95998fSHong Zhang . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 989f95998fSHong Zhang is optional with PETSC_TRUE being the default 990a29876aSHong Zhang 1000a29876aSHong Zhang Level: intermediate 1010a29876aSHong Zhang 1020a29876aSHong Zhang .keywords: PC, indefinite, factorization 1030a29876aSHong Zhang 104ee45ca4aSHong Zhang .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero() 1050a29876aSHong Zhang @*/ 106*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd(PC pc,PetscTruth shift) 1070a29876aSHong Zhang { 108afaefe49SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 109afaefe49SHong Zhang 1100a29876aSHong Zhang PetscFunctionBegin; 111afaefe49SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 112afaefe49SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);CHKERRQ(ierr); 113afaefe49SHong Zhang if (f) { 114afaefe49SHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 115afaefe49SHong Zhang } 1160a29876aSHong Zhang PetscFunctionReturn(0); 1170a29876aSHong Zhang } 118