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