1dba47a55SKris Buschelman 29b54502bSHong Zhang /* 39b54502bSHong Zhang Defines a ILU factorization preconditioner for any Mat implementation 49b54502bSHong Zhang */ 5c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/ilu/ilu.h> /*I "petscpc.h" I*/ 69b54502bSHong Zhang 72401956bSBarry Smith #undef __FUNCT__ 82401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 97087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 109b54502bSHong Zhang { 119b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 129b54502bSHong Zhang 139b54502bSHong Zhang PetscFunctionBegin; 149b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 152fa5cd67SKarl Rupp if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10; 162fa5cd67SKarl Rupp else ilu->nonzerosalongdiagonaltol = z; 179b54502bSHong Zhang PetscFunctionReturn(0); 189b54502bSHong Zhang } 199b54502bSHong Zhang 209b54502bSHong Zhang #undef __FUNCT__ 21574deadeSBarry Smith #define __FUNCT__ "PCReset_ILU" 22574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc) 239b54502bSHong Zhang { 249b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 259b54502bSHong Zhang PetscErrorCode ierr; 269b54502bSHong Zhang 279b54502bSHong Zhang PetscFunctionBegin; 28*3d1c1ea0SBarry Smith if (!ilu->hdr.inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 296bf464f9SBarry Smith if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);} 30fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 319b54502bSHong Zhang PetscFunctionReturn(0); 329b54502bSHong Zhang } 339b54502bSHong Zhang 349b54502bSHong Zhang #undef __FUNCT__ 35b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_ILU" 367087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 379b54502bSHong Zhang { 38075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 399b54502bSHong Zhang 409b54502bSHong Zhang PetscFunctionBegin; 414c9036c7SBarry Smith if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 42ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 439b54502bSHong Zhang } 44075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = dt; 45075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = dtcol; 46075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = dtcount; 474c9036c7SBarry Smith ((PC_Factor*)ilu)->info.usedt = 1.0; 489b54502bSHong Zhang PetscFunctionReturn(0); 499b54502bSHong Zhang } 509b54502bSHong Zhang 519b54502bSHong Zhang #undef __FUNCT__ 529b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 534416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc) 549b54502bSHong Zhang { 559b54502bSHong Zhang PetscErrorCode ierr; 5678fc6b22SHong Zhang PetscInt itmp; 578afaa268SBarry Smith PetscBool flg,set; 589b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 599b54502bSHong Zhang PetscReal tol; 609b54502bSHong Zhang 619b54502bSHong Zhang PetscFunctionBegin; 62e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ILU Options");CHKERRQ(ierr); 63e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 648ff23777SHong Zhang 65075768bcSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 66075768bcSBarry Smith if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 672fa5cd67SKarl Rupp 688afaa268SBarry Smith 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); 698afaa268SBarry Smith if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg; 702401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 719b54502bSHong Zhang if (flg) { 729b54502bSHong Zhang tol = PETSC_DECIDE; 732401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 742401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 759b54502bSHong Zhang } 769b54502bSHong Zhang 779b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 789b54502bSHong Zhang PetscFunctionReturn(0); 799b54502bSHong Zhang } 809b54502bSHong Zhang 819b54502bSHong Zhang #undef __FUNCT__ 829b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 839b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 849b54502bSHong Zhang { 859b54502bSHong Zhang PetscErrorCode ierr; 869b54502bSHong Zhang 879b54502bSHong Zhang PetscFunctionBegin; 88914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 899b54502bSHong Zhang PetscFunctionReturn(0); 909b54502bSHong Zhang } 919b54502bSHong Zhang 929b54502bSHong Zhang #undef __FUNCT__ 939b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 949b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 959b54502bSHong Zhang { 969b54502bSHong Zhang PetscErrorCode ierr; 979b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 98f3a39becSBarry Smith MatInfo info; 99ace3abfcSBarry Smith PetscBool flg; 10000c67f3bSHong Zhang const MatSolverPackage stype; 10100e125f8SBarry Smith MatFactorError err; 1029b54502bSHong Zhang 1039b54502bSHong Zhang PetscFunctionBegin; 10492927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 10592927226SBarry Smith if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 10622d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 10792927226SBarry Smith if (!flg) { 10822d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr); 10992927226SBarry Smith if (!flg) { 11092927226SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 111955c1f14SBarry Smith PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");CHKERRQ(ierr); 11292927226SBarry Smith } 11392927226SBarry Smith } 11492927226SBarry Smith } 11592927226SBarry Smith 11684d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 117*3d1c1ea0SBarry Smith if (ilu->hdr.inplace) { 1189b54502bSHong Zhang if (!pc->setupcalled) { 1199b54502bSHong Zhang 1209b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 1219b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 122075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1233bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1243bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1259b54502bSHong Zhang } 1269b54502bSHong Zhang 1279b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1289b54502bSHong Zhang cannot have levels of fill */ 129075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 13075567043SBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 1312fa5cd67SKarl Rupp 132365a8a9eSBarry Smith ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKERRQ(ierr); 13300e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 13400e125f8SBarry Smith if (err) { /* Factor() fails */ 13500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1366baea169SHong Zhang PetscFunctionReturn(0); 1376baea169SHong Zhang } 1386baea169SHong Zhang 139075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 140b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 141b33856dcSBarry Smith ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr); 1429b54502bSHong Zhang } else { 1439b54502bSHong Zhang if (!pc->setupcalled) { 1449b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 145075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1463bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1473bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1489b54502bSHong Zhang /* Remove zeros along diagonal? */ 1499b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1509b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 1519b54502bSHong Zhang } 15286daa42cSBarry Smith if (!((PC_Factor*)ilu)->fact) { 15311bfe483SHong Zhang ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 154a1f19f5aSHong Zhang } 155075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 156075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 157*3d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1582fa5cd67SKarl Rupp 1593bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1609b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 161*3d1c1ea0SBarry Smith if (!ilu->hdr.reuseordering) { 1629b54502bSHong Zhang /* compute a new ordering for the ILU */ 163fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->row);CHKERRQ(ierr); 164fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 165075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1663bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1673bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1689b54502bSHong Zhang /* Remove zeros along diagonal? */ 1699b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1709b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 1719b54502bSHong Zhang } 1729b54502bSHong Zhang } 1736bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 17414798fb4SJed Brown ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 175075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 176075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 177*3d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1782fa5cd67SKarl Rupp 1793bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1809b54502bSHong Zhang } 18100e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 18200e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 18300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1846baea169SHong Zhang PetscFunctionReturn(0); 1856baea169SHong Zhang } 1866baea169SHong Zhang 187075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 18800e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 18900e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 19000e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1916baea169SHong Zhang } 1929b54502bSHong Zhang } 19300c67f3bSHong Zhang 19400c67f3bSHong Zhang ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 19500c67f3bSHong Zhang if (!stype) { 19600e125f8SBarry Smith const MatSolverPackage solverpackage; 19700e125f8SBarry Smith ierr = MatFactorGetSolverPackage(((PC_Factor*)ilu)->fact,&solverpackage);CHKERRQ(ierr); 19800e125f8SBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr); 19900c67f3bSHong Zhang } 2009b54502bSHong Zhang PetscFunctionReturn(0); 2019b54502bSHong Zhang } 2029b54502bSHong Zhang 2039b54502bSHong Zhang #undef __FUNCT__ 2049b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 2059b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 2069b54502bSHong Zhang { 2079b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2089b54502bSHong Zhang PetscErrorCode ierr; 2099b54502bSHong Zhang 2109b54502bSHong Zhang PetscFunctionBegin; 211574deadeSBarry Smith ierr = PCReset_ILU(pc);CHKERRQ(ierr); 212503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 213503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 214c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2159b54502bSHong Zhang PetscFunctionReturn(0); 2169b54502bSHong Zhang } 2179b54502bSHong Zhang 2189b54502bSHong Zhang #undef __FUNCT__ 2199b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 2209b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 2219b54502bSHong Zhang { 2229b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2239b54502bSHong Zhang PetscErrorCode ierr; 2249b54502bSHong Zhang 2259b54502bSHong Zhang PetscFunctionBegin; 226075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2279b54502bSHong Zhang PetscFunctionReturn(0); 2289b54502bSHong Zhang } 2299b54502bSHong Zhang 2309b54502bSHong Zhang #undef __FUNCT__ 2319b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 2329b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 2339b54502bSHong Zhang { 2349b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2359b54502bSHong Zhang PetscErrorCode ierr; 2369b54502bSHong Zhang 2379b54502bSHong Zhang PetscFunctionBegin; 238075768bcSBarry Smith ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2399b54502bSHong Zhang PetscFunctionReturn(0); 2409b54502bSHong Zhang } 2419b54502bSHong Zhang 242f0b9ad6cSBarry Smith #undef __FUNCT__ 243f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricLeft_ILU" 244f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 245f0b9ad6cSBarry Smith { 246f0b9ad6cSBarry Smith PetscErrorCode ierr; 247f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 248f0b9ad6cSBarry Smith 249f0b9ad6cSBarry Smith PetscFunctionBegin; 250f0b9ad6cSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 251f0b9ad6cSBarry Smith PetscFunctionReturn(0); 252f0b9ad6cSBarry Smith } 253f0b9ad6cSBarry Smith 254f0b9ad6cSBarry Smith #undef __FUNCT__ 255f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricRight_ILU" 256f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 257f0b9ad6cSBarry Smith { 258f0b9ad6cSBarry Smith PetscErrorCode ierr; 259f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 260f0b9ad6cSBarry Smith 261f0b9ad6cSBarry Smith PetscFunctionBegin; 262f0b9ad6cSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 263f0b9ad6cSBarry Smith PetscFunctionReturn(0); 264f0b9ad6cSBarry Smith } 265f0b9ad6cSBarry Smith 2669b54502bSHong Zhang /*MC 2679b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2689b54502bSHong Zhang 2699b54502bSHong Zhang Options Database Keys: 2702401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2712401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2729b54502bSHong Zhang its factorization (overwrites original matrix) 2732401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2742401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 27555ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2762401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 2779b54502bSHong Zhang this decreases the chance of getting a zero pivot 2782401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 279967c93d3SBarry Smith - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 2809b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 2819b54502bSHong Zhang stability of the ILU factorization 2829b54502bSHong Zhang 2839b54502bSHong Zhang Level: beginner 2849b54502bSHong Zhang 2859b54502bSHong Zhang Concepts: incomplete factorization 2869b54502bSHong Zhang 287d9ba8547SBarry Smith Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 2889b54502bSHong Zhang 2899b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 2909b54502bSHong Zhang 291f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 292f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 293f0b9ad6cSBarry Smith 294cd0a26f6SPaul Mullowney If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization 295cd0a26f6SPaul Mullowney is never done on the GPU). 296ac793be5SBarry Smith 297c582cd25SBarry Smith References: 29896a0c994SBarry Smith + 1. - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 29996a0c994SBarry Smith self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 30096a0c994SBarry Smith . 2. - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 30196a0c994SBarry Smith - 3. - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 30296a0c994SBarry Smith Chapter in Parallel Numerical 303c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 30496a0c994SBarry Smith Science and Engineering, Kluwer. 305c582cd25SBarry Smith 306c582cd25SBarry Smith 3079b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 308145b9266SHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 309b7c853c4SBarry Smith PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 31092e9c092SBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 31192e9c092SBarry Smith PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace() 3129b54502bSHong Zhang 3139b54502bSHong Zhang M*/ 3149b54502bSHong Zhang 3159b54502bSHong Zhang #undef __FUNCT__ 3169b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 3178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 3189b54502bSHong Zhang { 3199b54502bSHong Zhang PetscErrorCode ierr; 3209b54502bSHong Zhang PC_ILU *ilu; 3219b54502bSHong Zhang 3229b54502bSHong Zhang PetscFunctionBegin; 323b00a9115SJed Brown ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr); 324*3d1c1ea0SBarry Smith pc->data = (void*)ilu; 325*3d1c1ea0SBarry Smith ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 3269b54502bSHong Zhang 327879e8a4dSBarry Smith ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU; 32875567043SBarry Smith ((PC_Factor*)ilu)->info.levels = 0.; 329075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 3309b54502bSHong Zhang ilu->col = 0; 3319b54502bSHong Zhang ilu->row = 0; 33219fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 333075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 334075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 335075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 3369b54502bSHong Zhang 337574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3389b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3399b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3409b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3419b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3429b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 3439b54502bSHong Zhang pc->ops->view = PCView_ILU; 344f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 345f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3469b54502bSHong Zhang pc->ops->applyrichardson = 0; 347bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 348bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 3499b54502bSHong Zhang PetscFunctionReturn(0); 3509b54502bSHong Zhang } 351