1*9b54502bSHong Zhang /* 2*9b54502bSHong Zhang Defines a Cholesky factorization preconditioner for any Mat implementation. 3*9b54502bSHong Zhang Presently only provided for MPIRowbs format (i.e. BlockSolve). 4*9b54502bSHong Zhang */ 5*9b54502bSHong Zhang 6*9b54502bSHong Zhang #include "src/ksp/pc/impls/icc/icc.h" /*I "petscpc.h" I*/ 7*9b54502bSHong Zhang 8*9b54502bSHong Zhang EXTERN_C_BEGIN 9*9b54502bSHong Zhang #undef __FUNCT__ 10*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering_ICC" 11*9b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering) 12*9b54502bSHong Zhang { 13*9b54502bSHong Zhang PC_ICC *dir = (PC_ICC*)pc->data; 14*9b54502bSHong Zhang PetscErrorCode ierr; 15*9b54502bSHong Zhang 16*9b54502bSHong Zhang PetscFunctionBegin; 17*9b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 18*9b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 19*9b54502bSHong Zhang PetscFunctionReturn(0); 20*9b54502bSHong Zhang } 21*9b54502bSHong Zhang EXTERN_C_END 22*9b54502bSHong Zhang 23*9b54502bSHong Zhang EXTERN_C_BEGIN 24*9b54502bSHong Zhang #undef __FUNCT__ 25*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetDamping_ICC" 26*9b54502bSHong Zhang PetscErrorCode PCICCSetDamping_ICC(PC pc,PetscReal damping) 27*9b54502bSHong Zhang { 28*9b54502bSHong Zhang PC_ICC *dir; 29*9b54502bSHong Zhang 30*9b54502bSHong Zhang PetscFunctionBegin; 31*9b54502bSHong Zhang dir = (PC_ICC*)pc->data; 32*9b54502bSHong Zhang if (damping == (PetscReal) PETSC_DECIDE) { 33*9b54502bSHong Zhang dir->info.damping = 1.e-12; 34*9b54502bSHong Zhang } else { 35*9b54502bSHong Zhang dir->info.damping = damping; 36*9b54502bSHong Zhang } 37*9b54502bSHong Zhang PetscFunctionReturn(0); 38*9b54502bSHong Zhang } 39*9b54502bSHong Zhang EXTERN_C_END 40*9b54502bSHong Zhang 41*9b54502bSHong Zhang EXTERN_C_BEGIN 42*9b54502bSHong Zhang #undef __FUNCT__ 43*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetShift_ICC" 44*9b54502bSHong Zhang PetscErrorCode PCICCSetShift_ICC(PC pc,PetscTruth shift) 45*9b54502bSHong Zhang { 46*9b54502bSHong Zhang PC_ICC *dir; 47*9b54502bSHong Zhang 48*9b54502bSHong Zhang PetscFunctionBegin; 49*9b54502bSHong Zhang dir = (PC_ICC*)pc->data; 50*9b54502bSHong Zhang dir->info.shift = shift; 51*9b54502bSHong Zhang if (shift) dir->info.shift_fraction = 0.0; 52*9b54502bSHong Zhang PetscFunctionReturn(0); 53*9b54502bSHong Zhang } 54*9b54502bSHong Zhang EXTERN_C_END 55*9b54502bSHong Zhang 56*9b54502bSHong Zhang EXTERN_C_BEGIN 57*9b54502bSHong Zhang #undef __FUNCT__ 58*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetSetZeroPivot_ICC" 59*9b54502bSHong Zhang PetscErrorCode PCICCSetZeroPivot_ICC(PC pc,PetscReal z) 60*9b54502bSHong Zhang { 61*9b54502bSHong Zhang PC_ICC *lu; 62*9b54502bSHong Zhang 63*9b54502bSHong Zhang PetscFunctionBegin; 64*9b54502bSHong Zhang lu = (PC_ICC*)pc->data; 65*9b54502bSHong Zhang lu->info.zeropivot = z; 66*9b54502bSHong Zhang PetscFunctionReturn(0); 67*9b54502bSHong Zhang } 68*9b54502bSHong Zhang EXTERN_C_END 69*9b54502bSHong Zhang 70*9b54502bSHong Zhang EXTERN_C_BEGIN 71*9b54502bSHong Zhang #undef __FUNCT__ 72*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill_ICC" 73*9b54502bSHong Zhang PetscErrorCode PCICCSetFill_ICC(PC pc,PetscReal fill) 74*9b54502bSHong Zhang { 75*9b54502bSHong Zhang PC_ICC *dir; 76*9b54502bSHong Zhang 77*9b54502bSHong Zhang PetscFunctionBegin; 78*9b54502bSHong Zhang dir = (PC_ICC*)pc->data; 79*9b54502bSHong Zhang dir->info.fill = fill; 80*9b54502bSHong Zhang PetscFunctionReturn(0); 81*9b54502bSHong Zhang } 82*9b54502bSHong Zhang EXTERN_C_END 83*9b54502bSHong Zhang 84*9b54502bSHong Zhang EXTERN_C_BEGIN 85*9b54502bSHong Zhang #undef __FUNCT__ 86*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels_ICC" 87*9b54502bSHong Zhang PetscErrorCode PCICCSetLevels_ICC(PC pc,PetscInt levels) 88*9b54502bSHong Zhang { 89*9b54502bSHong Zhang PC_ICC *icc; 90*9b54502bSHong Zhang 91*9b54502bSHong Zhang PetscFunctionBegin; 92*9b54502bSHong Zhang icc = (PC_ICC*)pc->data; 93*9b54502bSHong Zhang icc->info.levels = levels; 94*9b54502bSHong Zhang PetscFunctionReturn(0); 95*9b54502bSHong Zhang } 96*9b54502bSHong Zhang EXTERN_C_END 97*9b54502bSHong Zhang 98*9b54502bSHong Zhang #undef __FUNCT__ 99*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering" 100*9b54502bSHong Zhang /*@ 101*9b54502bSHong Zhang PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to 102*9b54502bSHong Zhang be used it the ICC factorization. 103*9b54502bSHong Zhang 104*9b54502bSHong Zhang Collective on PC 105*9b54502bSHong Zhang 106*9b54502bSHong Zhang Input Parameters: 107*9b54502bSHong Zhang + pc - the preconditioner context 108*9b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 109*9b54502bSHong Zhang 110*9b54502bSHong Zhang Options Database Key: 111*9b54502bSHong Zhang . -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine 112*9b54502bSHong Zhang 113*9b54502bSHong Zhang Level: intermediate 114*9b54502bSHong Zhang 115*9b54502bSHong Zhang .seealso: PCLUSetMatOrdering() 116*9b54502bSHong Zhang 117*9b54502bSHong Zhang .keywords: PC, ICC, set, matrix, reordering 118*9b54502bSHong Zhang 119*9b54502bSHong Zhang @*/ 120*9b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering(PC pc,MatOrderingType ordering) 121*9b54502bSHong Zhang { 122*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 123*9b54502bSHong Zhang 124*9b54502bSHong Zhang PetscFunctionBegin; 125*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 126*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 127*9b54502bSHong Zhang if (f) { 128*9b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 129*9b54502bSHong Zhang } 130*9b54502bSHong Zhang PetscFunctionReturn(0); 131*9b54502bSHong Zhang } 132*9b54502bSHong Zhang 133*9b54502bSHong Zhang #undef __FUNCT__ 134*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels" 135*9b54502bSHong Zhang /*@ 136*9b54502bSHong Zhang PCICCSetLevels - Sets the number of levels of fill to use. 137*9b54502bSHong Zhang 138*9b54502bSHong Zhang Collective on PC 139*9b54502bSHong Zhang 140*9b54502bSHong Zhang Input Parameters: 141*9b54502bSHong Zhang + pc - the preconditioner context 142*9b54502bSHong Zhang - levels - number of levels of fill 143*9b54502bSHong Zhang 144*9b54502bSHong Zhang Options Database Key: 145*9b54502bSHong Zhang . -pc_icc_levels <levels> - Sets fill level 146*9b54502bSHong Zhang 147*9b54502bSHong Zhang Level: intermediate 148*9b54502bSHong Zhang 149*9b54502bSHong Zhang Concepts: ICC^setting levels of fill 150*9b54502bSHong Zhang 151*9b54502bSHong Zhang @*/ 152*9b54502bSHong Zhang PetscErrorCode PCICCSetLevels(PC pc,PetscInt levels) 153*9b54502bSHong Zhang { 154*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 155*9b54502bSHong Zhang 156*9b54502bSHong Zhang PetscFunctionBegin; 157*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 158*9b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 159*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 160*9b54502bSHong Zhang if (f) { 161*9b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 162*9b54502bSHong Zhang } 163*9b54502bSHong Zhang PetscFunctionReturn(0); 164*9b54502bSHong Zhang } 165*9b54502bSHong Zhang 166*9b54502bSHong Zhang #undef __FUNCT__ 167*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill" 168*9b54502bSHong Zhang /*@ 169*9b54502bSHong Zhang PCICCSetFill - Indicate the amount of fill you expect in the factored matrix, 170*9b54502bSHong Zhang where fill = number nonzeros in factor/number nonzeros in original matrix. 171*9b54502bSHong Zhang 172*9b54502bSHong Zhang Collective on PC 173*9b54502bSHong Zhang 174*9b54502bSHong Zhang Input Parameters: 175*9b54502bSHong Zhang + pc - the preconditioner context 176*9b54502bSHong Zhang - fill - amount of expected fill 177*9b54502bSHong Zhang 178*9b54502bSHong Zhang Options Database Key: 179*9b54502bSHong Zhang $ -pc_icc_fill <fill> 180*9b54502bSHong Zhang 181*9b54502bSHong Zhang Note: 182*9b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 183*9b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 184*9b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 185*9b54502bSHong Zhang future runs. But default PETSc uses a value of 1.0 186*9b54502bSHong Zhang 187*9b54502bSHong Zhang Level: intermediate 188*9b54502bSHong Zhang 189*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 190*9b54502bSHong Zhang 191*9b54502bSHong Zhang .seealso: PCLUSetFill() 192*9b54502bSHong Zhang @*/ 193*9b54502bSHong Zhang PetscErrorCode PCICCSetFill(PC pc,PetscReal fill) 194*9b54502bSHong Zhang { 195*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 196*9b54502bSHong Zhang 197*9b54502bSHong Zhang PetscFunctionBegin; 198*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 199*9b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 200*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 201*9b54502bSHong Zhang if (f) { 202*9b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 203*9b54502bSHong Zhang } 204*9b54502bSHong Zhang PetscFunctionReturn(0); 205*9b54502bSHong Zhang } 206*9b54502bSHong Zhang 207*9b54502bSHong Zhang #undef __FUNCT__ 208*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetDamping" 209*9b54502bSHong Zhang /*@ 210*9b54502bSHong Zhang PCICCSetDamping - adds this quantity to the diagonal of the matrix during the 211*9b54502bSHong Zhang ICC numerical factorization 212*9b54502bSHong Zhang 213*9b54502bSHong Zhang Collective on PC 214*9b54502bSHong Zhang 215*9b54502bSHong Zhang Input Parameters: 216*9b54502bSHong Zhang + pc - the preconditioner context 217*9b54502bSHong Zhang - damping - amount of damping 218*9b54502bSHong Zhang 219*9b54502bSHong Zhang Options Database Key: 220*9b54502bSHong Zhang . -pc_icc_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 221*9b54502bSHong Zhang 222*9b54502bSHong Zhang Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 223*9b54502bSHong Zhang pivot, then the damping is doubled until this is alleviated. 224*9b54502bSHong Zhang 225*9b54502bSHong Zhang Level: intermediate 226*9b54502bSHong Zhang 227*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 228*9b54502bSHong Zhang 229*9b54502bSHong Zhang .seealso: PCICCSetFill(), PCLUSetDamping() 230*9b54502bSHong Zhang @*/ 231*9b54502bSHong Zhang PetscErrorCode PCICCSetDamping(PC pc,PetscReal damping) 232*9b54502bSHong Zhang { 233*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 234*9b54502bSHong Zhang 235*9b54502bSHong Zhang PetscFunctionBegin; 236*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 237*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 238*9b54502bSHong Zhang if (f) { 239*9b54502bSHong Zhang ierr = (*f)(pc,damping);CHKERRQ(ierr); 240*9b54502bSHong Zhang } 241*9b54502bSHong Zhang PetscFunctionReturn(0); 242*9b54502bSHong Zhang } 243*9b54502bSHong Zhang 244*9b54502bSHong Zhang #undef __FUNCT__ 245*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetShift" 246*9b54502bSHong Zhang /*@ 247*9b54502bSHong Zhang PCICCSetShift - specify whether to use Manteuffel shifting of ICC. 248*9b54502bSHong Zhang If an ICC factorisation breaks down because of nonpositive pivots, 249*9b54502bSHong Zhang adding sufficient identity to the diagonal will remedy this. 250*9b54502bSHong Zhang 251*9b54502bSHong Zhang Manteuffel shifting for ICC uses a different algorithm than the ILU case. 252*9b54502bSHong Zhang Here we base the shift on the lack of diagonal dominance when a negative 253*9b54502bSHong Zhang pivot occurs. 254*9b54502bSHong Zhang 255*9b54502bSHong Zhang Input parameters: 256*9b54502bSHong Zhang + pc - the preconditioner context 257*9b54502bSHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 258*9b54502bSHong Zhang 259*9b54502bSHong Zhang Options Database Key: 260*9b54502bSHong Zhang . -pc_icc_shift - Activate PCICCSetShift() 261*9b54502bSHong Zhang 262*9b54502bSHong Zhang Level: intermediate 263*9b54502bSHong Zhang 264*9b54502bSHong Zhang .keywords: PC, indefinite, factorization, incomplete, ICC 265*9b54502bSHong Zhang 266*9b54502bSHong Zhang .seealso: PCILUSetShift() 267*9b54502bSHong Zhang @*/ 268*9b54502bSHong Zhang PetscErrorCode PCICCSetShift(PC pc,PetscTruth shift) 269*9b54502bSHong Zhang { 270*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 271*9b54502bSHong Zhang 272*9b54502bSHong Zhang PetscFunctionBegin; 273*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 274*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 275*9b54502bSHong Zhang if (f) { 276*9b54502bSHong Zhang ierr = (*f)(pc,shift);CHKERRQ(ierr); 277*9b54502bSHong Zhang } 278*9b54502bSHong Zhang PetscFunctionReturn(0); 279*9b54502bSHong Zhang } 280*9b54502bSHong Zhang 281*9b54502bSHong Zhang #undef __FUNCT__ 282*9b54502bSHong Zhang #define __FUNCT__ "PCICCSetZeroPivot" 283*9b54502bSHong Zhang /*@ 284*9b54502bSHong Zhang PCICCSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 285*9b54502bSHong Zhang 286*9b54502bSHong Zhang Collective on PC 287*9b54502bSHong Zhang 288*9b54502bSHong Zhang Input Parameters: 289*9b54502bSHong Zhang + pc - the preconditioner context 290*9b54502bSHong Zhang - zero - all pivots smaller than this will be considered zero 291*9b54502bSHong Zhang 292*9b54502bSHong Zhang Options Database Key: 293*9b54502bSHong Zhang . -pc_ilu_zeropivot <zero> - Sets the zero pivot size 294*9b54502bSHong Zhang 295*9b54502bSHong Zhang Level: intermediate 296*9b54502bSHong Zhang 297*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 298*9b54502bSHong Zhang 299*9b54502bSHong Zhang .seealso: PCICCSetFill(), PCLUSetDamp(), PCLUSetZeroPivot() 300*9b54502bSHong Zhang @*/ 301*9b54502bSHong Zhang PetscErrorCode PCICCSetZeroPivot(PC pc,PetscReal zero) 302*9b54502bSHong Zhang { 303*9b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 304*9b54502bSHong Zhang 305*9b54502bSHong Zhang PetscFunctionBegin; 306*9b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 307*9b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 308*9b54502bSHong Zhang if (f) { 309*9b54502bSHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 310*9b54502bSHong Zhang } 311*9b54502bSHong Zhang PetscFunctionReturn(0); 312*9b54502bSHong Zhang } 313*9b54502bSHong Zhang 314*9b54502bSHong Zhang #undef __FUNCT__ 315*9b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 316*9b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 317*9b54502bSHong Zhang { 318*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 319*9b54502bSHong Zhang IS perm,cperm; 320*9b54502bSHong Zhang PetscErrorCode ierr; 321*9b54502bSHong Zhang 322*9b54502bSHong Zhang PetscFunctionBegin; 323*9b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 324*9b54502bSHong Zhang 325*9b54502bSHong Zhang if (!pc->setupcalled) { 326*9b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 327*9b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 328*9b54502bSHong Zhang ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 329*9b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 330*9b54502bSHong Zhang } 331*9b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 332*9b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 333*9b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 334*9b54502bSHong Zhang PetscFunctionReturn(0); 335*9b54502bSHong Zhang } 336*9b54502bSHong Zhang 337*9b54502bSHong Zhang #undef __FUNCT__ 338*9b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 339*9b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 340*9b54502bSHong Zhang { 341*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 342*9b54502bSHong Zhang PetscErrorCode ierr; 343*9b54502bSHong Zhang 344*9b54502bSHong Zhang PetscFunctionBegin; 345*9b54502bSHong Zhang if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 346*9b54502bSHong Zhang ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 347*9b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 348*9b54502bSHong Zhang PetscFunctionReturn(0); 349*9b54502bSHong Zhang } 350*9b54502bSHong Zhang 351*9b54502bSHong Zhang #undef __FUNCT__ 352*9b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 353*9b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 354*9b54502bSHong Zhang { 355*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 356*9b54502bSHong Zhang PetscErrorCode ierr; 357*9b54502bSHong Zhang 358*9b54502bSHong Zhang PetscFunctionBegin; 359*9b54502bSHong Zhang ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 360*9b54502bSHong Zhang PetscFunctionReturn(0); 361*9b54502bSHong Zhang } 362*9b54502bSHong Zhang 363*9b54502bSHong Zhang #undef __FUNCT__ 364*9b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 365*9b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 366*9b54502bSHong Zhang { 367*9b54502bSHong Zhang PetscErrorCode ierr; 368*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 369*9b54502bSHong Zhang 370*9b54502bSHong Zhang PetscFunctionBegin; 371*9b54502bSHong Zhang ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 372*9b54502bSHong Zhang PetscFunctionReturn(0); 373*9b54502bSHong Zhang } 374*9b54502bSHong Zhang 375*9b54502bSHong Zhang #undef __FUNCT__ 376*9b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 377*9b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 378*9b54502bSHong Zhang { 379*9b54502bSHong Zhang PetscErrorCode ierr; 380*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 381*9b54502bSHong Zhang 382*9b54502bSHong Zhang PetscFunctionBegin; 383*9b54502bSHong Zhang ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 384*9b54502bSHong Zhang PetscFunctionReturn(0); 385*9b54502bSHong Zhang } 386*9b54502bSHong Zhang 387*9b54502bSHong Zhang #undef __FUNCT__ 388*9b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ICC" 389*9b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 390*9b54502bSHong Zhang { 391*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 392*9b54502bSHong Zhang 393*9b54502bSHong Zhang PetscFunctionBegin; 394*9b54502bSHong Zhang *mat = icc->fact; 395*9b54502bSHong Zhang PetscFunctionReturn(0); 396*9b54502bSHong Zhang } 397*9b54502bSHong Zhang 398*9b54502bSHong Zhang #undef __FUNCT__ 399*9b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 400*9b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 401*9b54502bSHong Zhang { 402*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 403*9b54502bSHong Zhang char tname[256]; 404*9b54502bSHong Zhang PetscTruth flg; 405*9b54502bSHong Zhang PetscErrorCode ierr; 406*9b54502bSHong Zhang PetscFList ordlist; 407*9b54502bSHong Zhang 408*9b54502bSHong Zhang PetscFunctionBegin; 409*9b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 410*9b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 411*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 412*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 413*9b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 414*9b54502bSHong Zhang ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 415*9b54502bSHong Zhang if (flg) { 416*9b54502bSHong Zhang ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr); 417*9b54502bSHong Zhang } 418*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",&flg);CHKERRQ(ierr); 419*9b54502bSHong Zhang if (flg) { 420*9b54502bSHong Zhang ierr = PCICCSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 421*9b54502bSHong Zhang } 422*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_damping","Damping added to diagonal","PCICCSetDamping",icc->info.damping,&icc->info.damping,0);CHKERRQ(ierr); 423*9b54502bSHong Zhang ierr = PetscOptionsName("-pc_icc_shift","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr); 424*9b54502bSHong Zhang if (flg) { 425*9b54502bSHong Zhang ierr = PCICCSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 426*9b54502bSHong Zhang } else { 427*9b54502bSHong Zhang ierr = PCICCSetShift(pc,PETSC_FALSE);CHKERRQ(ierr); 428*9b54502bSHong Zhang } 429*9b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_zeropivot","Pivot is considered zero if less than","PCICCSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 430*9b54502bSHong Zhang 431*9b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 432*9b54502bSHong Zhang PetscFunctionReturn(0); 433*9b54502bSHong Zhang } 434*9b54502bSHong Zhang 435*9b54502bSHong Zhang #undef __FUNCT__ 436*9b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 437*9b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 438*9b54502bSHong Zhang { 439*9b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 440*9b54502bSHong Zhang PetscErrorCode ierr; 441*9b54502bSHong Zhang PetscTruth isstring,iascii; 442*9b54502bSHong Zhang 443*9b54502bSHong Zhang PetscFunctionBegin; 444*9b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 445*9b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 446*9b54502bSHong Zhang if (iascii) { 447*9b54502bSHong Zhang if (icc->info.levels == 1) { 448*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 449*9b54502bSHong Zhang } else { 450*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 451*9b54502bSHong Zhang } 452*9b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr); 453*9b54502bSHong Zhang if (icc->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 454*9b54502bSHong Zhang } else if (isstring) { 455*9b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 456*9b54502bSHong Zhang } else { 457*9b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 458*9b54502bSHong Zhang } 459*9b54502bSHong Zhang PetscFunctionReturn(0); 460*9b54502bSHong Zhang } 461*9b54502bSHong Zhang 462*9b54502bSHong Zhang /*MC 463*9b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 464*9b54502bSHong Zhang 465*9b54502bSHong Zhang Options Database Keys: 466*9b54502bSHong Zhang + -pc_icc_levels <k> - number of levels of fill for ICC(k) 467*9b54502bSHong Zhang . -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 468*9b54502bSHong Zhang its factorization (overwrites original matrix) 469*9b54502bSHong Zhang . -pc_icc_damping - add damping to diagonal to prevent zero (or very small) pivots 470*9b54502bSHong Zhang . -pc_icc_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner 471*9b54502bSHong Zhang . -pc_icc_zeropivot <tol> - set tolerance for what is considered a zero pivot 472*9b54502bSHong Zhang . -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 473*9b54502bSHong Zhang - -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 474*9b54502bSHong Zhang 475*9b54502bSHong Zhang Level: beginner 476*9b54502bSHong Zhang 477*9b54502bSHong Zhang Concepts: incomplete Cholesky factorization 478*9b54502bSHong Zhang 479*9b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 480*9b54502bSHong Zhang 481*9b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 482*9b54502bSHong Zhang 483*9b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 484*9b54502bSHong Zhang 485*9b54502bSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE); 486*9b54502bSHong Zhang to turn off the shift. 487*9b54502bSHong Zhang 488*9b54502bSHong Zhang 489*9b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 490*9b54502bSHong Zhang PCICCSetSetZeroPivot(), PCICCSetDamping(), PCICCSetShift(), 491*9b54502bSHong Zhang PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(), 492*9b54502bSHong Zhang PCICCSetLevels() 493*9b54502bSHong Zhang 494*9b54502bSHong Zhang M*/ 495*9b54502bSHong Zhang 496*9b54502bSHong Zhang EXTERN_C_BEGIN 497*9b54502bSHong Zhang #undef __FUNCT__ 498*9b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 499*9b54502bSHong Zhang PetscErrorCode PCCreate_ICC(PC pc) 500*9b54502bSHong Zhang { 501*9b54502bSHong Zhang PetscErrorCode ierr; 502*9b54502bSHong Zhang PC_ICC *icc; 503*9b54502bSHong Zhang 504*9b54502bSHong Zhang PetscFunctionBegin; 505*9b54502bSHong Zhang ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr); 506*9b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_ICC)); 507*9b54502bSHong Zhang 508*9b54502bSHong Zhang icc->fact = 0; 509*9b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 510*9b54502bSHong Zhang ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 511*9b54502bSHong Zhang icc->info.levels = 0; 512*9b54502bSHong Zhang icc->info.fill = 1.0; 513*9b54502bSHong Zhang icc->implctx = 0; 514*9b54502bSHong Zhang 515*9b54502bSHong Zhang icc->info.dtcol = PETSC_DEFAULT; 516*9b54502bSHong Zhang icc->info.damping = 0.0; 517*9b54502bSHong Zhang icc->info.shift = PETSC_TRUE; 518*9b54502bSHong Zhang icc->info.shift_fraction = 0.0; 519*9b54502bSHong Zhang icc->info.zeropivot = 1.e-12; 520*9b54502bSHong Zhang pc->data = (void*)icc; 521*9b54502bSHong Zhang 522*9b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 523*9b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 524*9b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 525*9b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 526*9b54502bSHong Zhang pc->ops->view = PCView_ICC; 527*9b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 528*9b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 529*9b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 530*9b54502bSHong Zhang 531*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC", 532*9b54502bSHong Zhang PCICCSetLevels_ICC);CHKERRQ(ierr); 533*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC", 534*9b54502bSHong Zhang PCICCSetFill_ICC);CHKERRQ(ierr); 535*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetDamping_C","PCICCSetDamping_ICC", 536*9b54502bSHong Zhang PCICCSetDamping_ICC);CHKERRQ(ierr); 537*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetShift_C","PCICCSetShift_ICC", 538*9b54502bSHong Zhang PCICCSetShift_ICC);CHKERRQ(ierr); 539*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC", 540*9b54502bSHong Zhang PCICCSetMatOrdering_ICC);CHKERRQ(ierr); 541*9b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetZeroPivot_C","PCICCSetZeroPivot_ICC", 542*9b54502bSHong Zhang PCICCSetZeroPivot_ICC);CHKERRQ(ierr); 543*9b54502bSHong Zhang PetscFunctionReturn(0); 544*9b54502bSHong Zhang } 545*9b54502bSHong Zhang EXTERN_C_END 546*9b54502bSHong Zhang 547*9b54502bSHong Zhang 548