1 2 /* 3 Defines a ILU factorization preconditioner for any Mat implementation 4 */ 5 #include <../src/ksp/pc/impls/factor/ilu/ilu.h> /*I "petscpc.h" I*/ 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 9 PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 10 { 11 PC_ILU *ilu = (PC_ILU*)pc->data; 12 13 PetscFunctionBegin; 14 ilu->nonzerosalongdiagonal = PETSC_TRUE; 15 if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10; 16 else ilu->nonzerosalongdiagonaltol = z; 17 PetscFunctionReturn(0); 18 } 19 20 #undef __FUNCT__ 21 #define __FUNCT__ "PCReset_ILU" 22 PetscErrorCode PCReset_ILU(PC pc) 23 { 24 PC_ILU *ilu = (PC_ILU*)pc->data; 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 if (!ilu->hdr.inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 29 if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);} 30 ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 31 PetscFunctionReturn(0); 32 } 33 34 #undef __FUNCT__ 35 #define __FUNCT__ "PCFactorSetDropTolerance_ILU" 36 PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 37 { 38 PC_ILU *ilu = (PC_ILU*)pc->data; 39 40 PetscFunctionBegin; 41 if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 42 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 43 } 44 ((PC_Factor*)ilu)->info.dt = dt; 45 ((PC_Factor*)ilu)->info.dtcol = dtcol; 46 ((PC_Factor*)ilu)->info.dtcount = dtcount; 47 ((PC_Factor*)ilu)->info.usedt = 1.0; 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "PCSetFromOptions_ILU" 53 static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc) 54 { 55 PetscErrorCode ierr; 56 PetscInt itmp; 57 PetscBool flg,set; 58 PC_ILU *ilu = (PC_ILU*)pc->data; 59 PetscReal tol; 60 61 PetscFunctionBegin; 62 ierr = PetscOptionsHead(PetscOptionsObject,"ILU Options");CHKERRQ(ierr); 63 ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 64 65 ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 66 if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 67 68 ierr = PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 69 if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg; 70 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 71 if (flg) { 72 tol = PETSC_DECIDE; 73 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 74 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 75 } 76 77 ierr = PetscOptionsTail();CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "PCView_ILU" 83 static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 84 { 85 PetscErrorCode ierr; 86 87 PetscFunctionBegin; 88 ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 #undef __FUNCT__ 93 #define __FUNCT__ "PCSetUp_ILU" 94 static PetscErrorCode PCSetUp_ILU(PC pc) 95 { 96 PetscErrorCode ierr; 97 PC_ILU *ilu = (PC_ILU*)pc->data; 98 MatInfo info; 99 PetscBool flg; 100 const MatSolverPackage stype; 101 MatFactorError err; 102 103 PetscFunctionBegin; 104 pc->failedreason = PC_NOERROR; 105 /* ugly hack to change default, since it is not support by some matrix types */ 106 if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 107 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 108 if (!flg) { 109 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr); 110 if (!flg) { 111 ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 112 PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");CHKERRQ(ierr); 113 } 114 } 115 } 116 117 ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 118 if (ilu->hdr.inplace) { 119 if (!pc->setupcalled) { 120 121 /* In-place factorization only makes sense with the natural ordering, 122 so we only need to get the ordering once, even if nonzero structure changes */ 123 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 124 if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 125 if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 126 } 127 128 /* In place ILU only makes sense with fill factor of 1.0 because 129 cannot have levels of fill */ 130 ((PC_Factor*)ilu)->info.fill = 1.0; 131 ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 132 133 ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKERRQ(ierr); 134 ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 135 if (err) { /* Factor() fails */ 136 pc->failedreason = (PCFailedReason)err; 137 PetscFunctionReturn(0); 138 } 139 140 ((PC_Factor*)ilu)->fact = pc->pmat; 141 /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 142 ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr); 143 } else { 144 if (!pc->setupcalled) { 145 /* first time in so compute reordering and symbolic factorization */ 146 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 147 if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 148 if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 149 /* Remove zeros along diagonal? */ 150 if (ilu->nonzerosalongdiagonal) { 151 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 152 } 153 if (!((PC_Factor*)ilu)->fact) { 154 ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 155 } 156 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 157 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 158 ilu->hdr.actualfill = info.fill_ratio_needed; 159 160 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 161 } else if (pc->flag != SAME_NONZERO_PATTERN) { 162 if (!ilu->hdr.reuseordering) { 163 /* compute a new ordering for the ILU */ 164 ierr = ISDestroy(&ilu->row);CHKERRQ(ierr); 165 ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 166 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 167 if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 168 if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 169 /* Remove zeros along diagonal? */ 170 if (ilu->nonzerosalongdiagonal) { 171 ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 172 } 173 } 174 ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 175 ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 176 ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 177 ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 178 ilu->hdr.actualfill = info.fill_ratio_needed; 179 180 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 181 } 182 ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 183 if (err) { /* FactorSymbolic() fails */ 184 pc->failedreason = (PCFailedReason)err; 185 PetscFunctionReturn(0); 186 } 187 188 ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 189 ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 190 if (err) { /* FactorNumeric() fails */ 191 pc->failedreason = (PCFailedReason)err; 192 } 193 } 194 195 ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 196 if (!stype) { 197 const MatSolverPackage solverpackage; 198 ierr = MatFactorGetSolverPackage(((PC_Factor*)ilu)->fact,&solverpackage);CHKERRQ(ierr); 199 ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr); 200 } 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "PCDestroy_ILU" 206 static PetscErrorCode PCDestroy_ILU(PC pc) 207 { 208 PC_ILU *ilu = (PC_ILU*)pc->data; 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 ierr = PCReset_ILU(pc);CHKERRQ(ierr); 213 ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 214 ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 215 ierr = PetscFree(pc->data);CHKERRQ(ierr); 216 PetscFunctionReturn(0); 217 } 218 219 #undef __FUNCT__ 220 #define __FUNCT__ "PCApply_ILU" 221 static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 222 { 223 PC_ILU *ilu = (PC_ILU*)pc->data; 224 PetscErrorCode ierr; 225 226 PetscFunctionBegin; 227 ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNCT__ 232 #define __FUNCT__ "PCApplyTranspose_ILU" 233 static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 234 { 235 PC_ILU *ilu = (PC_ILU*)pc->data; 236 PetscErrorCode ierr; 237 238 PetscFunctionBegin; 239 ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 #undef __FUNCT__ 244 #define __FUNCT__ "PCApplySymmetricLeft_ILU" 245 static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 246 { 247 PetscErrorCode ierr; 248 PC_ILU *icc = (PC_ILU*)pc->data; 249 250 PetscFunctionBegin; 251 ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "PCApplySymmetricRight_ILU" 257 static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 258 { 259 PetscErrorCode ierr; 260 PC_ILU *icc = (PC_ILU*)pc->data; 261 262 PetscFunctionBegin; 263 ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 /*MC 268 PCILU - Incomplete factorization preconditioners. 269 270 Options Database Keys: 271 + -pc_factor_levels <k> - number of levels of fill for ILU(k) 272 . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 273 its factorization (overwrites original matrix) 274 . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 275 . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 276 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 277 . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 278 this decreases the chance of getting a zero pivot 279 . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 280 - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 281 than 1 the diagonal blocks are factored with partial pivoting (this increases the 282 stability of the ILU factorization 283 284 Level: beginner 285 286 Concepts: incomplete factorization 287 288 Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 289 290 For BAIJ matrices this implements a point block ILU 291 292 The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 293 even when the matrix is not symmetric since the U stores the diagonals of the factorization. 294 295 If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization 296 is never done on the GPU). 297 298 References: 299 + 1. - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 300 self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 301 . 2. - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 302 - 3. - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 303 Chapter in Parallel Numerical 304 Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 305 Science and Engineering, Kluwer. 306 307 308 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 309 PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 310 PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 311 PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 312 PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace() 313 314 M*/ 315 316 #undef __FUNCT__ 317 #define __FUNCT__ "PCCreate_ILU" 318 PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 319 { 320 PetscErrorCode ierr; 321 PC_ILU *ilu; 322 323 PetscFunctionBegin; 324 ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr); 325 pc->data = (void*)ilu; 326 ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 327 328 ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU; 329 ((PC_Factor*)ilu)->info.levels = 0.; 330 ((PC_Factor*)ilu)->info.fill = 1.0; 331 ilu->col = 0; 332 ilu->row = 0; 333 ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 334 ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 335 ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 336 ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 337 338 pc->ops->reset = PCReset_ILU; 339 pc->ops->destroy = PCDestroy_ILU; 340 pc->ops->apply = PCApply_ILU; 341 pc->ops->applytranspose = PCApplyTranspose_ILU; 342 pc->ops->setup = PCSetUp_ILU; 343 pc->ops->setfromoptions = PCSetFromOptions_ILU; 344 pc->ops->view = PCView_ILU; 345 pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 346 pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 347 pc->ops->applyrichardson = 0; 348 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 349 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352