19b54502bSHong Zhang /* 29b54502bSHong Zhang Defines a Cholesky factorization preconditioner for any Mat implementation. 39b54502bSHong Zhang Presently only provided for MPIRowbs format (i.e. BlockSolve). 49b54502bSHong Zhang */ 59b54502bSHong Zhang 6*85b398ccSKris Buschelman #include "src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 79b54502bSHong Zhang 89b54502bSHong Zhang EXTERN_C_BEGIN 99b54502bSHong Zhang #undef __FUNCT__ 109b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering_ICC" 119b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering) 129b54502bSHong Zhang { 139b54502bSHong Zhang PC_ICC *dir = (PC_ICC*)pc->data; 149b54502bSHong Zhang PetscErrorCode ierr; 159b54502bSHong Zhang 169b54502bSHong Zhang PetscFunctionBegin; 179b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 189b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 199b54502bSHong Zhang PetscFunctionReturn(0); 209b54502bSHong Zhang } 219b54502bSHong Zhang EXTERN_C_END 229b54502bSHong Zhang 239b54502bSHong Zhang EXTERN_C_BEGIN 249b54502bSHong Zhang #undef __FUNCT__ 259b54502bSHong Zhang #define __FUNCT__ "PCICCSetDamping_ICC" 269b54502bSHong Zhang PetscErrorCode PCICCSetDamping_ICC(PC pc,PetscReal damping) 279b54502bSHong Zhang { 289b54502bSHong Zhang PC_ICC *dir; 299b54502bSHong Zhang 309b54502bSHong Zhang PetscFunctionBegin; 319b54502bSHong Zhang dir = (PC_ICC*)pc->data; 329b54502bSHong Zhang if (damping == (PetscReal) PETSC_DECIDE) { 339b54502bSHong Zhang dir->info.damping = 1.e-12; 349b54502bSHong Zhang } else { 359b54502bSHong Zhang dir->info.damping = damping; 369b54502bSHong Zhang } 379b54502bSHong Zhang PetscFunctionReturn(0); 389b54502bSHong Zhang } 399b54502bSHong Zhang EXTERN_C_END 409b54502bSHong Zhang 419b54502bSHong Zhang EXTERN_C_BEGIN 429b54502bSHong Zhang #undef __FUNCT__ 439b54502bSHong Zhang #define __FUNCT__ "PCICCSetShift_ICC" 449b54502bSHong Zhang PetscErrorCode PCICCSetShift_ICC(PC pc,PetscTruth shift) 459b54502bSHong Zhang { 469b54502bSHong Zhang PC_ICC *dir; 479b54502bSHong Zhang 489b54502bSHong Zhang PetscFunctionBegin; 499b54502bSHong Zhang dir = (PC_ICC*)pc->data; 509b54502bSHong Zhang dir->info.shift = shift; 519b54502bSHong Zhang if (shift) dir->info.shift_fraction = 0.0; 529b54502bSHong Zhang PetscFunctionReturn(0); 539b54502bSHong Zhang } 549b54502bSHong Zhang EXTERN_C_END 559b54502bSHong Zhang 569b54502bSHong Zhang EXTERN_C_BEGIN 579b54502bSHong Zhang #undef __FUNCT__ 589b54502bSHong Zhang #define __FUNCT__ "PCICCSetSetZeroPivot_ICC" 599b54502bSHong Zhang PetscErrorCode PCICCSetZeroPivot_ICC(PC pc,PetscReal z) 609b54502bSHong Zhang { 619b54502bSHong Zhang PC_ICC *lu; 629b54502bSHong Zhang 639b54502bSHong Zhang PetscFunctionBegin; 649b54502bSHong Zhang lu = (PC_ICC*)pc->data; 659b54502bSHong Zhang lu->info.zeropivot = z; 669b54502bSHong Zhang PetscFunctionReturn(0); 679b54502bSHong Zhang } 689b54502bSHong Zhang EXTERN_C_END 699b54502bSHong Zhang 709b54502bSHong Zhang EXTERN_C_BEGIN 719b54502bSHong Zhang #undef __FUNCT__ 729b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill_ICC" 739b54502bSHong Zhang PetscErrorCode PCICCSetFill_ICC(PC pc,PetscReal fill) 749b54502bSHong Zhang { 759b54502bSHong Zhang PC_ICC *dir; 769b54502bSHong Zhang 779b54502bSHong Zhang PetscFunctionBegin; 789b54502bSHong Zhang dir = (PC_ICC*)pc->data; 799b54502bSHong Zhang dir->info.fill = fill; 809b54502bSHong Zhang PetscFunctionReturn(0); 819b54502bSHong Zhang } 829b54502bSHong Zhang EXTERN_C_END 839b54502bSHong Zhang 849b54502bSHong Zhang EXTERN_C_BEGIN 859b54502bSHong Zhang #undef __FUNCT__ 869b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels_ICC" 879b54502bSHong Zhang PetscErrorCode PCICCSetLevels_ICC(PC pc,PetscInt levels) 889b54502bSHong Zhang { 899b54502bSHong Zhang PC_ICC *icc; 909b54502bSHong Zhang 919b54502bSHong Zhang PetscFunctionBegin; 929b54502bSHong Zhang icc = (PC_ICC*)pc->data; 939b54502bSHong Zhang icc->info.levels = levels; 949b54502bSHong Zhang PetscFunctionReturn(0); 959b54502bSHong Zhang } 969b54502bSHong Zhang EXTERN_C_END 979b54502bSHong Zhang 989b54502bSHong Zhang #undef __FUNCT__ 999b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering" 1009b54502bSHong Zhang /*@ 1019b54502bSHong Zhang PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to 1029b54502bSHong Zhang be used it the ICC factorization. 1039b54502bSHong Zhang 1049b54502bSHong Zhang Collective on PC 1059b54502bSHong Zhang 1069b54502bSHong Zhang Input Parameters: 1079b54502bSHong Zhang + pc - the preconditioner context 1089b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 1099b54502bSHong Zhang 1109b54502bSHong Zhang Options Database Key: 1119b54502bSHong Zhang . -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine 1129b54502bSHong Zhang 1139b54502bSHong Zhang Level: intermediate 1149b54502bSHong Zhang 1159b54502bSHong Zhang .seealso: PCLUSetMatOrdering() 1169b54502bSHong Zhang 1179b54502bSHong Zhang .keywords: PC, ICC, set, matrix, reordering 1189b54502bSHong Zhang 1199b54502bSHong Zhang @*/ 1209b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering(PC pc,MatOrderingType ordering) 1219b54502bSHong Zhang { 1229b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 1239b54502bSHong Zhang 1249b54502bSHong Zhang PetscFunctionBegin; 1259b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1269b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 1279b54502bSHong Zhang if (f) { 1289b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 1299b54502bSHong Zhang } 1309b54502bSHong Zhang PetscFunctionReturn(0); 1319b54502bSHong Zhang } 1329b54502bSHong Zhang 1339b54502bSHong Zhang #undef __FUNCT__ 1349b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels" 1359b54502bSHong Zhang /*@ 1369b54502bSHong Zhang PCICCSetLevels - Sets the number of levels of fill to use. 1379b54502bSHong Zhang 1389b54502bSHong Zhang Collective on PC 1399b54502bSHong Zhang 1409b54502bSHong Zhang Input Parameters: 1419b54502bSHong Zhang + pc - the preconditioner context 1429b54502bSHong Zhang - levels - number of levels of fill 1439b54502bSHong Zhang 1449b54502bSHong Zhang Options Database Key: 1459b54502bSHong Zhang . -pc_icc_levels <levels> - Sets fill level 1469b54502bSHong Zhang 1479b54502bSHong Zhang Level: intermediate 1489b54502bSHong Zhang 1499b54502bSHong Zhang Concepts: ICC^setting levels of fill 1509b54502bSHong Zhang 1519b54502bSHong Zhang @*/ 1529b54502bSHong Zhang PetscErrorCode PCICCSetLevels(PC pc,PetscInt levels) 1539b54502bSHong Zhang { 1549b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 1559b54502bSHong Zhang 1569b54502bSHong Zhang PetscFunctionBegin; 1579b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1589b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 1599b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 1609b54502bSHong Zhang if (f) { 1619b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang PetscFunctionReturn(0); 1649b54502bSHong Zhang } 1659b54502bSHong Zhang 1669b54502bSHong Zhang #undef __FUNCT__ 1679b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill" 1689b54502bSHong Zhang /*@ 1699b54502bSHong Zhang PCICCSetFill - Indicate the amount of fill you expect in the factored matrix, 1709b54502bSHong Zhang where fill = number nonzeros in factor/number nonzeros in original matrix. 1719b54502bSHong Zhang 1729b54502bSHong Zhang Collective on PC 1739b54502bSHong Zhang 1749b54502bSHong Zhang Input Parameters: 1759b54502bSHong Zhang + pc - the preconditioner context 1769b54502bSHong Zhang - fill - amount of expected fill 1779b54502bSHong Zhang 1789b54502bSHong Zhang Options Database Key: 1799b54502bSHong Zhang $ -pc_icc_fill <fill> 1809b54502bSHong Zhang 1819b54502bSHong Zhang Note: 1829b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 1839b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 1849b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 1859b54502bSHong Zhang future runs. But default PETSc uses a value of 1.0 1869b54502bSHong Zhang 1879b54502bSHong Zhang Level: intermediate 1889b54502bSHong Zhang 1899b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 1909b54502bSHong Zhang 1919b54502bSHong Zhang .seealso: PCLUSetFill() 1929b54502bSHong Zhang @*/ 1939b54502bSHong Zhang PetscErrorCode PCICCSetFill(PC pc,PetscReal fill) 1949b54502bSHong Zhang { 1959b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 1969b54502bSHong Zhang 1979b54502bSHong Zhang PetscFunctionBegin; 1989b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1999b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 2009b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 2019b54502bSHong Zhang if (f) { 2029b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 2039b54502bSHong Zhang } 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang #undef __FUNCT__ 2089b54502bSHong Zhang #define __FUNCT__ "PCICCSetDamping" 2099b54502bSHong Zhang /*@ 2109b54502bSHong Zhang PCICCSetDamping - adds this quantity to the diagonal of the matrix during the 2119b54502bSHong Zhang ICC numerical factorization 2129b54502bSHong Zhang 2139b54502bSHong Zhang Collective on PC 2149b54502bSHong Zhang 2159b54502bSHong Zhang Input Parameters: 2169b54502bSHong Zhang + pc - the preconditioner context 2179b54502bSHong Zhang - damping - amount of damping 2189b54502bSHong Zhang 2199b54502bSHong Zhang Options Database Key: 2209b54502bSHong Zhang . -pc_icc_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 2219b54502bSHong Zhang 2229b54502bSHong Zhang Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 2239b54502bSHong Zhang pivot, then the damping is doubled until this is alleviated. 2249b54502bSHong Zhang 2259b54502bSHong Zhang Level: intermediate 2269b54502bSHong Zhang 2279b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 2289b54502bSHong Zhang 2299b54502bSHong Zhang .seealso: PCICCSetFill(), PCLUSetDamping() 2309b54502bSHong Zhang @*/ 2319b54502bSHong Zhang PetscErrorCode PCICCSetDamping(PC pc,PetscReal damping) 2329b54502bSHong Zhang { 2339b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 2349b54502bSHong Zhang 2359b54502bSHong Zhang PetscFunctionBegin; 2369b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2379b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 2389b54502bSHong Zhang if (f) { 2399b54502bSHong Zhang ierr = (*f)(pc,damping);CHKERRQ(ierr); 2409b54502bSHong Zhang } 2419b54502bSHong Zhang PetscFunctionReturn(0); 2429b54502bSHong Zhang } 2439b54502bSHong Zhang 2449b54502bSHong Zhang #undef __FUNCT__ 2459b54502bSHong Zhang #define __FUNCT__ "PCICCSetShift" 2469b54502bSHong Zhang /*@ 2479b54502bSHong Zhang PCICCSetShift - specify whether to use Manteuffel shifting of ICC. 2489b54502bSHong Zhang If an ICC factorisation breaks down because of nonpositive pivots, 2499b54502bSHong Zhang adding sufficient identity to the diagonal will remedy this. 2509b54502bSHong Zhang 2519b54502bSHong Zhang Manteuffel shifting for ICC uses a different algorithm than the ILU case. 2529b54502bSHong Zhang Here we base the shift on the lack of diagonal dominance when a negative 2539b54502bSHong Zhang pivot occurs. 2549b54502bSHong Zhang 2559b54502bSHong Zhang Input parameters: 2569b54502bSHong Zhang + pc - the preconditioner context 2579b54502bSHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 2589b54502bSHong Zhang 2599b54502bSHong Zhang Options Database Key: 2609b54502bSHong Zhang . -pc_icc_shift - Activate PCICCSetShift() 2619b54502bSHong Zhang 2629b54502bSHong Zhang Level: intermediate 2639b54502bSHong Zhang 2649b54502bSHong Zhang .keywords: PC, indefinite, factorization, incomplete, ICC 2659b54502bSHong Zhang 2669b54502bSHong Zhang .seealso: PCILUSetShift() 2679b54502bSHong Zhang @*/ 2689b54502bSHong Zhang PetscErrorCode PCICCSetShift(PC pc,PetscTruth shift) 2699b54502bSHong Zhang { 2709b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 2719b54502bSHong Zhang 2729b54502bSHong Zhang PetscFunctionBegin; 2739b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2749b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 2759b54502bSHong Zhang if (f) { 2769b54502bSHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 2779b54502bSHong Zhang } 2789b54502bSHong Zhang PetscFunctionReturn(0); 2799b54502bSHong Zhang } 2809b54502bSHong Zhang 2819b54502bSHong Zhang #undef __FUNCT__ 2829b54502bSHong Zhang #define __FUNCT__ "PCICCSetZeroPivot" 2839b54502bSHong Zhang /*@ 2849b54502bSHong Zhang PCICCSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 2859b54502bSHong Zhang 2869b54502bSHong Zhang Collective on PC 2879b54502bSHong Zhang 2889b54502bSHong Zhang Input Parameters: 2899b54502bSHong Zhang + pc - the preconditioner context 2909b54502bSHong Zhang - zero - all pivots smaller than this will be considered zero 2919b54502bSHong Zhang 2929b54502bSHong Zhang Options Database Key: 2939b54502bSHong Zhang . -pc_ilu_zeropivot <zero> - Sets the zero pivot size 2949b54502bSHong Zhang 2959b54502bSHong Zhang Level: intermediate 2969b54502bSHong Zhang 2979b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 2989b54502bSHong Zhang 2999b54502bSHong Zhang .seealso: PCICCSetFill(), PCLUSetDamp(), PCLUSetZeroPivot() 3009b54502bSHong Zhang @*/ 3019b54502bSHong Zhang PetscErrorCode PCICCSetZeroPivot(PC pc,PetscReal zero) 3029b54502bSHong Zhang { 3039b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 3049b54502bSHong Zhang 3059b54502bSHong Zhang PetscFunctionBegin; 3069b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3079b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 3089b54502bSHong Zhang if (f) { 3099b54502bSHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 3109b54502bSHong Zhang } 3119b54502bSHong Zhang PetscFunctionReturn(0); 3129b54502bSHong Zhang } 3139b54502bSHong Zhang 3149b54502bSHong Zhang #undef __FUNCT__ 3159b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 3169b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 3179b54502bSHong Zhang { 3189b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3199b54502bSHong Zhang IS perm,cperm; 3209b54502bSHong Zhang PetscErrorCode ierr; 3219b54502bSHong Zhang 3229b54502bSHong Zhang PetscFunctionBegin; 3239b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 3249b54502bSHong Zhang 3259b54502bSHong Zhang if (!pc->setupcalled) { 3269b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 3279b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 3289b54502bSHong Zhang ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 3299b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 3309b54502bSHong Zhang } 3319b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 3329b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 3339b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 3349b54502bSHong Zhang PetscFunctionReturn(0); 3359b54502bSHong Zhang } 3369b54502bSHong Zhang 3379b54502bSHong Zhang #undef __FUNCT__ 3389b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 3399b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 3409b54502bSHong Zhang { 3419b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3429b54502bSHong Zhang PetscErrorCode ierr; 3439b54502bSHong Zhang 3449b54502bSHong Zhang PetscFunctionBegin; 3459b54502bSHong Zhang if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 3469b54502bSHong Zhang ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 3479b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 3489b54502bSHong Zhang PetscFunctionReturn(0); 3499b54502bSHong Zhang } 3509b54502bSHong Zhang 3519b54502bSHong Zhang #undef __FUNCT__ 3529b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 3539b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 3549b54502bSHong Zhang { 3559b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3569b54502bSHong Zhang PetscErrorCode ierr; 3579b54502bSHong Zhang 3589b54502bSHong Zhang PetscFunctionBegin; 3599b54502bSHong Zhang ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 3609b54502bSHong Zhang PetscFunctionReturn(0); 3619b54502bSHong Zhang } 3629b54502bSHong Zhang 3639b54502bSHong Zhang #undef __FUNCT__ 3649b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 3659b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 3669b54502bSHong Zhang { 3679b54502bSHong Zhang PetscErrorCode ierr; 3689b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3699b54502bSHong Zhang 3709b54502bSHong Zhang PetscFunctionBegin; 3719b54502bSHong Zhang ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 3729b54502bSHong Zhang PetscFunctionReturn(0); 3739b54502bSHong Zhang } 3749b54502bSHong Zhang 3759b54502bSHong Zhang #undef __FUNCT__ 3769b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 3779b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 3789b54502bSHong Zhang { 3799b54502bSHong Zhang PetscErrorCode ierr; 3809b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3819b54502bSHong Zhang 3829b54502bSHong Zhang PetscFunctionBegin; 3839b54502bSHong Zhang ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 3849b54502bSHong Zhang PetscFunctionReturn(0); 3859b54502bSHong Zhang } 3869b54502bSHong Zhang 3879b54502bSHong Zhang #undef __FUNCT__ 3889b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ICC" 3899b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 3909b54502bSHong Zhang { 3919b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 3929b54502bSHong Zhang 3939b54502bSHong Zhang PetscFunctionBegin; 3949b54502bSHong Zhang *mat = icc->fact; 3959b54502bSHong Zhang PetscFunctionReturn(0); 3969b54502bSHong Zhang } 3979b54502bSHong Zhang 3989b54502bSHong Zhang #undef __FUNCT__ 3999b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 4009b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 4019b54502bSHong Zhang { 4029b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 4039b54502bSHong Zhang char tname[256]; 4049b54502bSHong Zhang PetscTruth flg; 4059b54502bSHong Zhang PetscErrorCode ierr; 4069b54502bSHong Zhang PetscFList ordlist; 4079b54502bSHong Zhang 4089b54502bSHong Zhang PetscFunctionBegin; 4099b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 4109b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 4119b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 4129b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 4139b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 4149b54502bSHong Zhang ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 4159b54502bSHong Zhang if (flg) { 4169b54502bSHong Zhang ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr); 4179b54502bSHong Zhang } 4189b54502bSHong Zhang ierr = PetscOptionsName("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",&flg);CHKERRQ(ierr); 4199b54502bSHong Zhang if (flg) { 4209b54502bSHong Zhang ierr = PCICCSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 4219b54502bSHong Zhang } 4229b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",icc->info.damping,&icc->info.damping,0);CHKERRQ(ierr); 4239b54502bSHong Zhang ierr = PetscOptionsName("-pc_icc_shift","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr); 4249b54502bSHong Zhang if (flg) { 4259b54502bSHong Zhang ierr = PCICCSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 4269b54502bSHong Zhang } else { 4279b54502bSHong Zhang ierr = PCICCSetShift(pc,PETSC_FALSE);CHKERRQ(ierr); 4289b54502bSHong Zhang } 4299b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_zeropivot","Pivot is considered zero if less than","PCICCSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 4309b54502bSHong Zhang 4319b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4329b54502bSHong Zhang PetscFunctionReturn(0); 4339b54502bSHong Zhang } 4349b54502bSHong Zhang 4359b54502bSHong Zhang #undef __FUNCT__ 4369b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 4379b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 4389b54502bSHong Zhang { 4399b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 4409b54502bSHong Zhang PetscErrorCode ierr; 4419b54502bSHong Zhang PetscTruth isstring,iascii; 4429b54502bSHong Zhang 4439b54502bSHong Zhang PetscFunctionBegin; 4449b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 4459b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 4469b54502bSHong Zhang if (iascii) { 4479b54502bSHong Zhang if (icc->info.levels == 1) { 4489b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 4499b54502bSHong Zhang } else { 4509b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 4519b54502bSHong Zhang } 4529b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr); 4539b54502bSHong Zhang if (icc->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 4549b54502bSHong Zhang } else if (isstring) { 4559b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 4569b54502bSHong Zhang } else { 4579b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 4589b54502bSHong Zhang } 4599b54502bSHong Zhang PetscFunctionReturn(0); 4609b54502bSHong Zhang } 4619b54502bSHong Zhang 4629b54502bSHong Zhang /*MC 4639b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 4649b54502bSHong Zhang 4659b54502bSHong Zhang Options Database Keys: 4669b54502bSHong Zhang + -pc_icc_levels <k> - number of levels of fill for ICC(k) 4679b54502bSHong Zhang . -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 4689b54502bSHong Zhang its factorization (overwrites original matrix) 4699b54502bSHong Zhang . -pc_icc_damping - add damping to diagonal to prevent zero (or very small) pivots 4709b54502bSHong Zhang . -pc_icc_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner 4719b54502bSHong Zhang . -pc_icc_zeropivot <tol> - set tolerance for what is considered a zero pivot 4729b54502bSHong Zhang . -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 4739b54502bSHong Zhang - -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 4749b54502bSHong Zhang 4759b54502bSHong Zhang Level: beginner 4769b54502bSHong Zhang 4779b54502bSHong Zhang Concepts: incomplete Cholesky factorization 4789b54502bSHong Zhang 4799b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 4809b54502bSHong Zhang 4819b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 4829b54502bSHong Zhang 4839b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 4849b54502bSHong Zhang 4859b54502bSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE); 4869b54502bSHong Zhang to turn off the shift. 4879b54502bSHong Zhang 4889b54502bSHong Zhang 4899b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 4909b54502bSHong Zhang PCICCSetSetZeroPivot(), PCICCSetDamping(), PCICCSetShift(), 4919b54502bSHong Zhang PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(), 4929b54502bSHong Zhang PCICCSetLevels() 4939b54502bSHong Zhang 4949b54502bSHong Zhang M*/ 4959b54502bSHong Zhang 4969b54502bSHong Zhang EXTERN_C_BEGIN 4979b54502bSHong Zhang #undef __FUNCT__ 4989b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 4999b54502bSHong Zhang PetscErrorCode PCCreate_ICC(PC pc) 5009b54502bSHong Zhang { 5019b54502bSHong Zhang PetscErrorCode ierr; 5029b54502bSHong Zhang PC_ICC *icc; 5039b54502bSHong Zhang 5049b54502bSHong Zhang PetscFunctionBegin; 5059b54502bSHong Zhang ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr); 5069b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_ICC)); 5079b54502bSHong Zhang 5089b54502bSHong Zhang icc->fact = 0; 5099b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 5109b54502bSHong Zhang ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 5119b54502bSHong Zhang icc->info.levels = 0; 5129b54502bSHong Zhang icc->info.fill = 1.0; 5139b54502bSHong Zhang icc->implctx = 0; 5149b54502bSHong Zhang 5159b54502bSHong Zhang icc->info.dtcol = PETSC_DEFAULT; 5169b54502bSHong Zhang icc->info.damping = 0.0; 5179b54502bSHong Zhang icc->info.shift = PETSC_TRUE; 5189b54502bSHong Zhang icc->info.shift_fraction = 0.0; 5199b54502bSHong Zhang icc->info.zeropivot = 1.e-12; 5209b54502bSHong Zhang pc->data = (void*)icc; 5219b54502bSHong Zhang 5229b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 5239b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 5249b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 5259b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 5269b54502bSHong Zhang pc->ops->view = PCView_ICC; 5279b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 5289b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 5299b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 5309b54502bSHong Zhang 5319b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC", 5329b54502bSHong Zhang PCICCSetLevels_ICC);CHKERRQ(ierr); 5339b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC", 5349b54502bSHong Zhang PCICCSetFill_ICC);CHKERRQ(ierr); 5359b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetDamping_C","PCICCSetDamping_ICC", 5369b54502bSHong Zhang PCICCSetDamping_ICC);CHKERRQ(ierr); 5379b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetShift_C","PCICCSetShift_ICC", 5389b54502bSHong Zhang PCICCSetShift_ICC);CHKERRQ(ierr); 5399b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC", 5409b54502bSHong Zhang PCICCSetMatOrdering_ICC);CHKERRQ(ierr); 5419b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetZeroPivot_C","PCICCSetZeroPivot_ICC", 5429b54502bSHong Zhang PCICCSetZeroPivot_ICC);CHKERRQ(ierr); 5439b54502bSHong Zhang PetscFunctionReturn(0); 5449b54502bSHong Zhang } 5459b54502bSHong Zhang EXTERN_C_END 5469b54502bSHong Zhang 5479b54502bSHong Zhang 548