1 #define PETSCKSP_DLL 2 3 /* 4 Defines a ILU factorization preconditioner for any Mat implementation 5 */ 6 #include "../src/ksp/pc/impls/factor/ilu/ilu.h" /*I "petscpc.h" I*/ 7 8 /* ------------------------------------------------------------------------------------------*/ 9 EXTERN_C_BEGIN 10 #undef __FUNCT__ 11 #define __FUNCT__ "PCFactorSetReuseFill_ILU" 12 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag) 13 { 14 PC_ILU *lu = (PC_ILU*)pc->data; 15 16 PetscFunctionBegin; 17 lu->reusefill = flag; 18 PetscFunctionReturn(0); 19 } 20 EXTERN_C_END 21 22 EXTERN_C_BEGIN 23 #undef __FUNCT__ 24 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 25 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 26 { 27 PC_ILU *ilu = (PC_ILU*)pc->data; 28 29 PetscFunctionBegin; 30 ilu->nonzerosalongdiagonal = PETSC_TRUE; 31 if (z == PETSC_DECIDE) { 32 ilu->nonzerosalongdiagonaltol = 1.e-10; 33 } else { 34 ilu->nonzerosalongdiagonaltol = z; 35 } 36 PetscFunctionReturn(0); 37 } 38 EXTERN_C_END 39 40 #undef __FUNCT__ 41 #define __FUNCT__ "PCDestroy_ILU_Internal" 42 PetscErrorCode PCDestroy_ILU_Internal(PC pc) 43 { 44 PC_ILU *ilu = (PC_ILU*)pc->data; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 49 if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 50 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 51 PetscFunctionReturn(0); 52 } 53 54 EXTERN_C_BEGIN 55 #undef __FUNCT__ 56 #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU" 57 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 58 { 59 PC_ILU *ilu = (PC_ILU*)pc->data; 60 PetscErrorCode ierr; 61 62 PetscFunctionBegin; 63 if (pc->setupcalled && (!ilu->usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 64 pc->setupcalled = 0; 65 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 66 } 67 ilu->usedt = PETSC_TRUE; 68 ((PC_Factor*)ilu)->info.dt = dt; 69 ((PC_Factor*)ilu)->info.dtcol = dtcol; 70 ((PC_Factor*)ilu)->info.dtcount = dtcount; 71 ((PC_Factor*)ilu)->info.fill = PETSC_DEFAULT; 72 PetscFunctionReturn(0); 73 } 74 EXTERN_C_END 75 76 EXTERN_C_BEGIN 77 #undef __FUNCT__ 78 #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 79 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag) 80 { 81 PC_ILU *ilu = (PC_ILU*)pc->data; 82 83 PetscFunctionBegin; 84 ilu->reuseordering = flag; 85 PetscFunctionReturn(0); 86 } 87 EXTERN_C_END 88 89 EXTERN_C_BEGIN 90 #undef __FUNCT__ 91 #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 92 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc) 93 { 94 PC_ILU *dir = (PC_ILU*)pc->data; 95 96 PetscFunctionBegin; 97 dir->inplace = PETSC_TRUE; 98 PetscFunctionReturn(0); 99 } 100 EXTERN_C_END 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "PCSetFromOptions_ILU" 104 static PetscErrorCode PCSetFromOptions_ILU(PC pc) 105 { 106 PetscErrorCode ierr; 107 PetscInt dtmax = 3,itmp; 108 PetscTruth flg,set; 109 PetscReal dt[3]; 110 char tname[256]; 111 PC_ILU *ilu = (PC_ILU*)pc->data; 112 PetscFList ordlist; 113 PetscReal tol; 114 115 PetscFunctionBegin; 116 if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 117 ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 118 ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 119 if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 120 ierr = PetscOptionsTruth("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",ilu->inplace,&ilu->inplace,PETSC_NULL);CHKERRQ(ierr); 121 flg = PETSC_FALSE; 122 ierr = PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 123 ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg; 124 ierr = PetscOptionsTruth("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",ilu->reusefill,&ilu->reusefill,PETSC_NULL);CHKERRQ(ierr); 125 ierr = PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",ilu->reuseordering,&ilu->reuseordering,PETSC_NULL);CHKERRQ(ierr); 126 flg = PETSC_FALSE; 127 ierr = PetscOptionsTruth("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 128 if (flg) { 129 ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 130 } 131 ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)ilu)->info.shiftnz,&((PC_Factor*)ilu)->info.shiftnz,0);CHKERRQ(ierr); 132 flg = (((PC_Factor*)ilu)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE; 133 ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 134 ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr); 135 ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)ilu)->info.zeropivot,&((PC_Factor*)ilu)->info.zeropivot,0);CHKERRQ(ierr); 136 137 dt[0] = ((PC_Factor*)ilu)->info.dt; 138 dt[1] = ((PC_Factor*)ilu)->info.dtcol; 139 dt[2] = ((PC_Factor*)ilu)->info.dtcount; 140 ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 141 if (flg) { 142 ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 143 } 144 ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)ilu)->info.fill,&((PC_Factor*)ilu)->info.fill,&flg);CHKERRQ(ierr); 145 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 146 if (flg) { 147 tol = PETSC_DECIDE; 148 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 149 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 150 } 151 152 ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 153 ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)ilu)->ordering,tname,256,&flg);CHKERRQ(ierr); 154 if (flg) { 155 ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 156 } 157 flg = ((PC_Factor*)ilu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 158 ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 159 if (set) { 160 ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 161 } 162 flg = PETSC_FALSE; 163 ierr = PetscOptionsTruth("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 164 if (flg) { 165 ierr = PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 166 } 167 ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",((PC_Factor*)ilu)->info.shiftinblocks,&((PC_Factor*)ilu)->info.shiftinblocks,0);CHKERRQ(ierr); 168 169 ierr = PetscOptionsTail();CHKERRQ(ierr); 170 PetscFunctionReturn(0); 171 } 172 173 #undef __FUNCT__ 174 #define __FUNCT__ "PCView_ILU" 175 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 176 { 177 PC_ILU *ilu = (PC_ILU*)pc->data; 178 PetscErrorCode ierr; 179 PetscTruth isstring,iascii; 180 181 PetscFunctionBegin; 182 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 183 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 184 if (iascii) { 185 if (ilu->usedt) { 186 ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",((PC_Factor*)ilu)->info.dt);CHKERRQ(ierr); 187 ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)((PC_Factor*)ilu)->info.dtcount);CHKERRQ(ierr); 188 ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",((PC_Factor*)ilu)->info.dtcol);CHKERRQ(ierr); 189 } else if (((PC_Factor*)ilu)->info.levels == 1) { 190 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr); 191 } else { 192 ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr); 193 } 194 ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio allocated %G\n",((PC_Factor*)ilu)->info.fill);CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",((PC_Factor*)ilu)->info.zeropivot);CHKERRQ(ierr); 196 if (((PC_Factor*)ilu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 197 if (((PC_Factor*)ilu)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);} 198 if (((PC_Factor*)ilu)->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);} 199 if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 200 else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 201 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 202 if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 203 if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 204 if (((PC_Factor*)ilu)->fact) { 205 ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr); 206 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 207 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 208 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 209 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 210 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 211 ierr = MatView(((PC_Factor*)ilu)->fact,viewer);CHKERRQ(ierr); 212 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 213 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 214 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 215 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 216 } 217 } else if (isstring) { 218 ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 219 } else { 220 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 221 } 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "PCSetUp_ILU" 227 static PetscErrorCode PCSetUp_ILU(PC pc) 228 { 229 PetscErrorCode ierr; 230 PC_ILU *ilu = (PC_ILU*)pc->data; 231 MatInfo info; 232 233 PetscFunctionBegin; 234 if (ilu->inplace) { 235 CHKMEMQ; 236 if (!pc->setupcalled) { 237 238 /* In-place factorization only makes sense with the natural ordering, 239 so we only need to get the ordering once, even if nonzero structure changes */ 240 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 241 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 242 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 243 } 244 245 /* In place ILU only makes sense with fill factor of 1.0 because 246 cannot have levels of fill */ 247 ((PC_Factor*)ilu)->info.fill = 1.0; 248 ((PC_Factor*)ilu)->info.diagonal_fill = 0; 249 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 250 CHKMEMQ; 251 ((PC_Factor*)ilu)->fact = pc->pmat; 252 } else if (ilu->usedt) { 253 if (!pc->setupcalled) { 254 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 255 CHKMEMQ; 256 ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 257 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 258 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 259 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 260 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 261 } else if (pc->flag != SAME_NONZERO_PATTERN) { 262 CHKMEMQ; 263 ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 264 CHKMEMQ; 265 if (!ilu->reuseordering) { 266 if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 267 if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 268 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 269 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 270 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 271 } 272 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 273 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 274 } else if (!ilu->reusefill) { 275 ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 276 ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 277 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 278 } else { 279 ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 280 } 281 } else { 282 if (!pc->setupcalled) { 283 /* first time in so compute reordering and symbolic factorization */ 284 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 285 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 286 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 287 /* Remove zeros along diagonal? */ 288 if (ilu->nonzerosalongdiagonal) { 289 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 290 } 291 CHKMEMQ; 292 ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 293 CHKMEMQ; 294 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 295 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 296 ilu->actualfill = info.fill_ratio_needed; 297 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 298 } else if (pc->flag != SAME_NONZERO_PATTERN) { 299 if (!ilu->reuseordering) { 300 /* compute a new ordering for the ILU */ 301 ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 302 ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 303 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 304 if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 305 if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 306 /* Remove zeros along diagonal? */ 307 if (ilu->nonzerosalongdiagonal) { 308 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 309 } 310 } 311 ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 312 ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 313 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 314 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 315 ilu->actualfill = info.fill_ratio_needed; 316 ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 317 } 318 CHKMEMQ; 319 ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 320 CHKMEMQ; 321 } 322 PetscFunctionReturn(0); 323 } 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "PCDestroy_ILU" 327 static PetscErrorCode PCDestroy_ILU(PC pc) 328 { 329 PC_ILU *ilu = (PC_ILU*)pc->data; 330 PetscErrorCode ierr; 331 332 PetscFunctionBegin; 333 ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 334 ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 335 ierr = PetscFree(ilu);CHKERRQ(ierr); 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "PCApply_ILU" 341 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 342 { 343 PC_ILU *ilu = (PC_ILU*)pc->data; 344 PetscErrorCode ierr; 345 346 PetscFunctionBegin; 347 ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 348 PetscFunctionReturn(0); 349 } 350 351 #undef __FUNCT__ 352 #define __FUNCT__ "PCApplyTranspose_ILU" 353 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 354 { 355 PC_ILU *ilu = (PC_ILU*)pc->data; 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 /*MC 364 PCILU - Incomplete factorization preconditioners. 365 366 Options Database Keys: 367 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 368 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 369 its factorization (overwrites original matrix) 370 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 371 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 372 . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 373 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 374 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 375 this decreases the chance of getting a zero pivot 376 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 377 . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 378 than 1 the diagonal blocks are factored with partial pivoting (this increases the 379 stability of the ILU factorization 380 . -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization 381 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 382 - -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value 383 is optional with true being the default 384 385 Level: beginner 386 387 Concepts: incomplete factorization 388 389 Notes: Only implemented for some matrix formats. (for parallel use you 390 must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended 391 unless you really want a parallel ILU). 392 393 For BAIJ matrices this implements a point block ILU 394 395 References: 396 T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 397 self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 398 399 T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 400 fusion problems. Quart. Appl. Math., 19:221--229, 1961. 401 402 Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 403 http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 404 Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 405 Science and Engineering, Kluwer, pp. 167--202. 406 407 408 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 409 PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 410 PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 411 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 412 PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 413 414 M*/ 415 416 EXTERN_C_BEGIN 417 #undef __FUNCT__ 418 #define __FUNCT__ "PCCreate_ILU" 419 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 420 { 421 PetscErrorCode ierr; 422 PC_ILU *ilu; 423 424 PetscFunctionBegin; 425 ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 426 427 ((PC_Factor*)ilu)->fact = 0; 428 ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 429 ((PC_Factor*)ilu)->info.levels = 0; 430 ((PC_Factor*)ilu)->info.fill = 1.0; 431 ilu->col = 0; 432 ilu->row = 0; 433 ilu->inplace = PETSC_FALSE; 434 ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 435 ilu->reuseordering = PETSC_FALSE; 436 ilu->usedt = PETSC_FALSE; 437 ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 438 ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 439 ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 440 ((PC_Factor*)ilu)->info.shiftnz = 1.e-12; 441 ((PC_Factor*)ilu)->info.shiftpd = 0.0; /* false */ 442 ((PC_Factor*)ilu)->info.zeropivot = 1.e-12; 443 ((PC_Factor*)ilu)->info.pivotinblocks = 1.0; 444 ((PC_Factor*)ilu)->info.shiftinblocks = 1.e-12; 445 ilu->reusefill = PETSC_FALSE; 446 ((PC_Factor*)ilu)->info.diagonal_fill = 0; 447 pc->data = (void*)ilu; 448 449 pc->ops->destroy = PCDestroy_ILU; 450 pc->ops->apply = PCApply_ILU; 451 pc->ops->applytranspose = PCApplyTranspose_ILU; 452 pc->ops->setup = PCSetUp_ILU; 453 pc->ops->setfromoptions = PCSetFromOptions_ILU; 454 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 455 pc->ops->view = PCView_ILU; 456 pc->ops->applyrichardson = 0; 457 458 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 459 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 460 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor", 461 PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr); 462 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor", 463 PCFactorSetShiftPd_Factor);CHKERRQ(ierr); 464 465 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 466 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 467 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 468 PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 469 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 470 PCFactorSetFill_Factor);CHKERRQ(ierr); 471 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 472 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 473 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 474 PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 475 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 476 PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 477 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 478 PCFactorSetLevels_Factor);CHKERRQ(ierr); 479 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 480 PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 481 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor", 482 PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr); 483 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 484 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 485 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor", 486 PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr); 487 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 488 PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 EXTERN_C_END 492