1 #define PETSCKSP_DLL 2 3 /* 4 Defines a Cholesky factorization preconditioner for any Mat implementation. 5 Presently only provided for MPIRowbs format (i.e. BlockSolve). 6 */ 7 8 #include "src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 9 10 EXTERN_C_BEGIN 11 #undef __FUNCT__ 12 #define __FUNCT__ "PCFactorSetZeroPivot_ICC" 13 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ICC(PC pc,PetscReal z) 14 { 15 PC_ICC *icc; 16 17 PetscFunctionBegin; 18 icc = (PC_ICC*)pc->data; 19 icc->info.zeropivot = z; 20 PetscFunctionReturn(0); 21 } 22 EXTERN_C_END 23 24 EXTERN_C_BEGIN 25 #undef __FUNCT__ 26 #define __FUNCT__ "PCFactorSetShiftNonzero_ICC" 27 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ICC(PC pc,PetscReal shift) 28 { 29 PC_ICC *dir; 30 31 PetscFunctionBegin; 32 dir = (PC_ICC*)pc->data; 33 if (shift == (PetscReal) PETSC_DECIDE) { 34 dir->info.shiftnz = 1.e-12; 35 } else { 36 dir->info.shiftnz = shift; 37 } 38 PetscFunctionReturn(0); 39 } 40 EXTERN_C_END 41 42 EXTERN_C_BEGIN 43 #undef __FUNCT__ 44 #define __FUNCT__ "PCFactorSetShiftPd_ICC" 45 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ICC(PC pc,PetscTruth shift) 46 { 47 PC_ICC *dir; 48 49 PetscFunctionBegin; 50 dir = (PC_ICC*)pc->data; 51 if (shift) { 52 dir->info.shift_fraction = 0.0; 53 dir->info.shiftpd = 1.0; 54 } else { 55 dir->info.shiftpd = 0.0; 56 } 57 PetscFunctionReturn(0); 58 } 59 EXTERN_C_END 60 61 EXTERN_C_BEGIN 62 #undef __FUNCT__ 63 #define __FUNCT__ "PCICCSetMatOrdering_ICC" 64 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering) 65 { 66 PC_ICC *dir = (PC_ICC*)pc->data; 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 71 ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 EXTERN_C_END 75 76 EXTERN_C_BEGIN 77 #undef __FUNCT__ 78 #define __FUNCT__ "PCICCSetFill_ICC" 79 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetFill_ICC(PC pc,PetscReal fill) 80 { 81 PC_ICC *dir; 82 83 PetscFunctionBegin; 84 dir = (PC_ICC*)pc->data; 85 dir->info.fill = fill; 86 PetscFunctionReturn(0); 87 } 88 EXTERN_C_END 89 90 EXTERN_C_BEGIN 91 #undef __FUNCT__ 92 #define __FUNCT__ "PCICCSetLevels_ICC" 93 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetLevels_ICC(PC pc,PetscInt levels) 94 { 95 PC_ICC *icc; 96 97 PetscFunctionBegin; 98 icc = (PC_ICC*)pc->data; 99 icc->info.levels = levels; 100 PetscFunctionReturn(0); 101 } 102 EXTERN_C_END 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "PCICCSetMatOrdering" 106 /*@ 107 PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to 108 be used it the ICC factorization. 109 110 Collective on PC 111 112 Input Parameters: 113 + pc - the preconditioner context 114 - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 115 116 Options Database Key: 117 . -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine 118 119 Level: intermediate 120 121 .seealso: PCLUSetMatOrdering() 122 123 .keywords: PC, ICC, set, matrix, reordering 124 125 @*/ 126 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetMatOrdering(PC pc,MatOrderingType ordering) 127 { 128 PetscErrorCode ierr,(*f)(PC,MatOrderingType); 129 130 PetscFunctionBegin; 131 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 132 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 133 if (f) { 134 ierr = (*f)(pc,ordering);CHKERRQ(ierr); 135 } 136 PetscFunctionReturn(0); 137 } 138 139 #undef __FUNCT__ 140 #define __FUNCT__ "PCICCSetLevels" 141 /*@ 142 PCICCSetLevels - Sets the number of levels of fill to use. 143 144 Collective on PC 145 146 Input Parameters: 147 + pc - the preconditioner context 148 - levels - number of levels of fill 149 150 Options Database Key: 151 . -pc_icc_levels <levels> - Sets fill level 152 153 Level: intermediate 154 155 Concepts: ICC^setting levels of fill 156 157 @*/ 158 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetLevels(PC pc,PetscInt levels) 159 { 160 PetscErrorCode ierr,(*f)(PC,PetscInt); 161 162 PetscFunctionBegin; 163 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 164 if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 165 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 166 if (f) { 167 ierr = (*f)(pc,levels);CHKERRQ(ierr); 168 } 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "PCICCSetFill" 174 /*@ 175 PCICCSetFill - Indicate the amount of fill you expect in the factored matrix, 176 where fill = number nonzeros in factor/number nonzeros in original matrix. 177 178 Collective on PC 179 180 Input Parameters: 181 + pc - the preconditioner context 182 - fill - amount of expected fill 183 184 Options Database Key: 185 $ -pc_icc_fill <fill> 186 187 Note: 188 For sparse matrix factorizations it is difficult to predict how much 189 fill to expect. By running with the option -log_info PETSc will print the 190 actual amount of fill used; allowing you to set the value accurately for 191 future runs. But default PETSc uses a value of 1.0 192 193 Level: intermediate 194 195 .keywords: PC, set, factorization, direct, fill 196 197 .seealso: PCLUSetFill() 198 @*/ 199 PetscErrorCode PETSCKSP_DLLEXPORT PCICCSetFill(PC pc,PetscReal fill) 200 { 201 PetscErrorCode ierr,(*f)(PC,PetscReal); 202 203 PetscFunctionBegin; 204 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 205 if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 206 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 207 if (f) { 208 ierr = (*f)(pc,fill);CHKERRQ(ierr); 209 } 210 PetscFunctionReturn(0); 211 } 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "PCSetup_ICC" 215 static PetscErrorCode PCSetup_ICC(PC pc) 216 { 217 PC_ICC *icc = (PC_ICC*)pc->data; 218 IS perm,cperm; 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 223 224 if (!pc->setupcalled) { 225 ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 226 } else if (pc->flag != SAME_NONZERO_PATTERN) { 227 ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 228 ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 229 } 230 ierr = ISDestroy(cperm);CHKERRQ(ierr); 231 ierr = ISDestroy(perm);CHKERRQ(ierr); 232 ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "PCDestroy_ICC" 238 static PetscErrorCode PCDestroy_ICC(PC pc) 239 { 240 PC_ICC *icc = (PC_ICC*)pc->data; 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 245 ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 246 ierr = PetscFree(icc);CHKERRQ(ierr); 247 PetscFunctionReturn(0); 248 } 249 250 #undef __FUNCT__ 251 #define __FUNCT__ "PCApply_ICC" 252 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 253 { 254 PC_ICC *icc = (PC_ICC*)pc->data; 255 PetscErrorCode ierr; 256 257 PetscFunctionBegin; 258 ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 #undef __FUNCT__ 263 #define __FUNCT__ "PCApplySymmetricLeft_ICC" 264 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 265 { 266 PetscErrorCode ierr; 267 PC_ICC *icc = (PC_ICC*)pc->data; 268 269 PetscFunctionBegin; 270 ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 271 PetscFunctionReturn(0); 272 } 273 274 #undef __FUNCT__ 275 #define __FUNCT__ "PCApplySymmetricRight_ICC" 276 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 277 { 278 PetscErrorCode ierr; 279 PC_ICC *icc = (PC_ICC*)pc->data; 280 281 PetscFunctionBegin; 282 ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 283 PetscFunctionReturn(0); 284 } 285 286 #undef __FUNCT__ 287 #define __FUNCT__ "PCGetFactoredMatrix_ICC" 288 static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 289 { 290 PC_ICC *icc = (PC_ICC*)pc->data; 291 292 PetscFunctionBegin; 293 *mat = icc->fact; 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "PCSetFromOptions_ICC" 299 static PetscErrorCode PCSetFromOptions_ICC(PC pc) 300 { 301 PC_ICC *icc = (PC_ICC*)pc->data; 302 char tname[256]; 303 PetscTruth flg; 304 PetscErrorCode ierr; 305 PetscFList ordlist; 306 307 PetscFunctionBegin; 308 ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 309 ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 310 ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 311 ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 312 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 313 ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 314 if (flg) { 315 ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr); 316 } 317 ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 318 if (flg) { 319 ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 320 } 321 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",icc->info.shiftnz,&icc->info.shiftnz,0);CHKERRQ(ierr); 322 ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr); 323 if (flg) { 324 ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 325 } else { 326 ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 327 } 328 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 329 330 ierr = PetscOptionsTail();CHKERRQ(ierr); 331 PetscFunctionReturn(0); 332 } 333 334 #undef __FUNCT__ 335 #define __FUNCT__ "PCView_ICC" 336 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 337 { 338 PC_ICC *icc = (PC_ICC*)pc->data; 339 PetscErrorCode ierr; 340 PetscTruth isstring,iascii; 341 342 PetscFunctionBegin; 343 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 344 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 345 if (iascii) { 346 if (icc->info.levels == 1) { 347 ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 348 } else { 349 ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 350 } 351 ierr = PetscViewerASCIIPrintf(viewer," ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr); 352 if (icc->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 353 } else if (isstring) { 354 ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 355 } else { 356 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 357 } 358 PetscFunctionReturn(0); 359 } 360 361 /*MC 362 PCICC - Incomplete Cholesky factorization preconditioners. 363 364 Options Database Keys: 365 + -pc_icc_levels <k> - number of levels of fill for ICC(k) 366 . -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 367 its factorization (overwrites original matrix) 368 . -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 369 . -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 370 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 371 - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 372 is optional with PETSC_TRUE being the default 373 374 Level: beginner 375 376 Concepts: incomplete Cholesky factorization 377 378 Notes: Only implemented for some matrix formats. Not implemented in parallel 379 380 For BAIJ matrices this implements a point block ICC. 381 382 The Manteuffel shift is only implemented for matrices with block size 1 383 384 By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE); 385 to turn off the shift. 386 387 388 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 389 PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 390 PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(), 391 PCICCSetLevels(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), 392 393 M*/ 394 395 EXTERN_C_BEGIN 396 #undef __FUNCT__ 397 #define __FUNCT__ "PCCreate_ICC" 398 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 399 { 400 PetscErrorCode ierr; 401 PC_ICC *icc; 402 403 PetscFunctionBegin; 404 ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr); 405 ierr = PetscLogObjectMemory(pc,sizeof(PC_ICC));CHKERRQ(ierr); 406 407 icc->fact = 0; 408 ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 409 ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 410 icc->info.levels = 0; 411 icc->info.fill = 1.0; 412 icc->implctx = 0; 413 414 icc->info.dtcol = PETSC_DEFAULT; 415 icc->info.shiftnz = 0.0; 416 icc->info.shiftpd = 1.0; /* true */ 417 icc->info.shift_fraction = 0.0; 418 icc->info.zeropivot = 1.e-12; 419 pc->data = (void*)icc; 420 421 pc->ops->apply = PCApply_ICC; 422 pc->ops->setup = PCSetup_ICC; 423 pc->ops->destroy = PCDestroy_ICC; 424 pc->ops->setfromoptions = PCSetFromOptions_ICC; 425 pc->ops->view = PCView_ICC; 426 pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 427 pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 428 pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 429 430 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ICC", 431 PCFactorSetZeroPivot_ICC);CHKERRQ(ierr); 432 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ICC", 433 PCFactorSetShiftNonzero_ICC);CHKERRQ(ierr); 434 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ICC", 435 PCFactorSetShiftPd_ICC);CHKERRQ(ierr); 436 437 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC", 438 PCICCSetLevels_ICC);CHKERRQ(ierr); 439 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC", 440 PCICCSetFill_ICC);CHKERRQ(ierr); 441 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC", 442 PCICCSetMatOrdering_ICC);CHKERRQ(ierr); 443 PetscFunctionReturn(0); 444 } 445 EXTERN_C_END 446 447 448