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 77087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 89b54502bSHong Zhang { 99b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 109b54502bSHong Zhang 119b54502bSHong Zhang PetscFunctionBegin; 129b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 132fa5cd67SKarl Rupp if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10; 142fa5cd67SKarl Rupp else ilu->nonzerosalongdiagonaltol = z; 159b54502bSHong Zhang PetscFunctionReturn(0); 169b54502bSHong Zhang } 179b54502bSHong Zhang 18574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc) 199b54502bSHong Zhang { 209b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 219b54502bSHong Zhang 229b54502bSHong Zhang PetscFunctionBegin; 239566063dSJacob Faibussowitsch if (!ilu->hdr.inplace) PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact)); 249566063dSJacob Faibussowitsch if (ilu->row && ilu->col && ilu->row != ilu->col) PetscCall(ISDestroy(&ilu->row)); 259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->col)); 269b54502bSHong Zhang PetscFunctionReturn(0); 279b54502bSHong Zhang } 289b54502bSHong Zhang 297087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 309b54502bSHong Zhang { 31075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 329b54502bSHong Zhang 339b54502bSHong Zhang PetscFunctionBegin; 344c9036c7SBarry Smith if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 35ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 369b54502bSHong Zhang } 37075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = dt; 38075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = dtcol; 39075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = dtcount; 404c9036c7SBarry Smith ((PC_Factor*)ilu)->info.usedt = 1.0; 419b54502bSHong Zhang PetscFunctionReturn(0); 429b54502bSHong Zhang } 439b54502bSHong Zhang 444416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc) 459b54502bSHong Zhang { 4678fc6b22SHong Zhang PetscInt itmp; 478afaa268SBarry Smith PetscBool flg,set; 489b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 499b54502bSHong Zhang PetscReal tol; 509b54502bSHong Zhang 519b54502bSHong Zhang PetscFunctionBegin; 52*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"ILU Options"); 539566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc)); 548ff23777SHong Zhang 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg)); 56075768bcSBarry Smith if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 572fa5cd67SKarl Rupp 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set)); 598afaa268SBarry Smith if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg)); 619b54502bSHong Zhang if (flg) { 629b54502bSHong Zhang tol = PETSC_DECIDE; 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,NULL)); 649566063dSJacob Faibussowitsch PetscCall(PCFactorReorderForNonzeroDiagonal(pc,tol)); 659b54502bSHong Zhang } 669b54502bSHong Zhang 67*d0609cedSBarry Smith PetscOptionsHeadEnd(); 689b54502bSHong Zhang PetscFunctionReturn(0); 699b54502bSHong Zhang } 709b54502bSHong Zhang 719b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 729b54502bSHong Zhang { 739b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 74f3a39becSBarry Smith MatInfo info; 75ace3abfcSBarry Smith PetscBool flg; 76ea799195SBarry Smith MatSolverType stype; 7700e125f8SBarry Smith MatFactorError err; 789b54502bSHong Zhang 799b54502bSHong Zhang PetscFunctionBegin; 80c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 8192927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 8292927226SBarry Smith if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg)); 8492927226SBarry Smith if (!flg) { 859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg)); 8692927226SBarry Smith if (!flg) { 8792927226SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 889566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n")); 8992927226SBarry Smith } 9092927226SBarry Smith } 9192927226SBarry Smith } 9292927226SBarry Smith 939566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure)); 943d1c1ea0SBarry Smith if (ilu->hdr.inplace) { 959b54502bSHong Zhang if (!pc->setupcalled) { 969b54502bSHong Zhang 979b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 989b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 992c7c0729SBarry Smith /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */ 1009566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact)); 1029566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1039566063dSJacob Faibussowitsch if (ilu->row) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1049566063dSJacob Faibussowitsch if (ilu->col) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1059b54502bSHong Zhang } 1069b54502bSHong Zhang 1079b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1089b54502bSHong Zhang cannot have levels of fill */ 109075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 11075567043SBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 1112fa5cd67SKarl Rupp 1129566063dSJacob Faibussowitsch PetscCall(MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1139566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat,&err)); 11400e125f8SBarry Smith if (err) { /* Factor() fails */ 11500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1166baea169SHong Zhang PetscFunctionReturn(0); 1176baea169SHong Zhang } 1186baea169SHong Zhang 119075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 120b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 1219566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate)); 1229b54502bSHong Zhang } else { 1239b54502bSHong Zhang if (!pc->setupcalled) { 1249b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 125f73b0415SBarry Smith PetscBool canuseordering; 1262c7c0729SBarry Smith if (!((PC_Factor*)ilu)->fact) { 1279566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact)); 1289566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact)); 1292c7c0729SBarry Smith } 1309566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering)); 131f73b0415SBarry Smith if (canuseordering) { 1329566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1339566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1349566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1359566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1369b54502bSHong Zhang /* Remove zeros along diagonal? */ 1379b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1389566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col)); 1399b54502bSHong Zhang } 140a1f19f5aSHong Zhang } 1419566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1429566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info)); 1433d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1449b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1453d1c1ea0SBarry Smith if (!ilu->hdr.reuseordering) { 146f73b0415SBarry Smith PetscBool canuseordering; 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact)); 1489566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact)); 1499566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact)); 1509566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering)); 151f73b0415SBarry Smith if (canuseordering) { 1529b54502bSHong Zhang /* compute a new ordering for the ILU */ 1539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->row)); 1549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->col)); 1559566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1569566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1579566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1589566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1599b54502bSHong Zhang /* Remove zeros along diagonal? */ 1609b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1619566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col)); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang } 1642c7c0729SBarry Smith } 1659566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1669566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info)); 1673d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1689b54502bSHong Zhang } 1699566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err)); 17000e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 17100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1726baea169SHong Zhang PetscFunctionReturn(0); 1736baea169SHong Zhang } 1746baea169SHong Zhang 1759566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info)); 1769566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err)); 17700e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 17800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1796baea169SHong Zhang } 1809b54502bSHong Zhang } 18100c67f3bSHong Zhang 1829566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 18300c67f3bSHong Zhang if (!stype) { 184ea799195SBarry Smith MatSolverType solverpackage; 1859566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage)); 1869566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 18700c67f3bSHong Zhang } 1889b54502bSHong Zhang PetscFunctionReturn(0); 1899b54502bSHong Zhang } 1909b54502bSHong Zhang 1919b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 1929b54502bSHong Zhang { 1939b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 1949b54502bSHong Zhang 1959b54502bSHong Zhang PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(PCReset_ILU(pc)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)ilu)->solvertype)); 1989566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)ilu)->ordering)); 1999566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2009b54502bSHong Zhang PetscFunctionReturn(0); 2019b54502bSHong Zhang } 2029b54502bSHong Zhang 2039b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 2049b54502bSHong Zhang { 2059b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2069b54502bSHong Zhang 2079b54502bSHong Zhang PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor*)ilu)->fact,x,y)); 2099b54502bSHong Zhang PetscFunctionReturn(0); 2109b54502bSHong Zhang } 2119b54502bSHong Zhang 2127b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ILU(PC pc,Mat X,Mat Y) 2137b6e2003SPierre Jolivet { 2147b6e2003SPierre Jolivet PC_ILU *ilu = (PC_ILU*)pc->data; 2157b6e2003SPierre Jolivet 2167b6e2003SPierre Jolivet PetscFunctionBegin; 2179566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor*)ilu)->fact,X,Y)); 2187b6e2003SPierre Jolivet PetscFunctionReturn(0); 2197b6e2003SPierre Jolivet } 2207b6e2003SPierre Jolivet 2219b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 2229b54502bSHong Zhang { 2239b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2249b54502bSHong Zhang 2259b54502bSHong Zhang PetscFunctionBegin; 2269566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y)); 2279b54502bSHong Zhang PetscFunctionReturn(0); 2289b54502bSHong Zhang } 2299b54502bSHong Zhang 230f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 231f0b9ad6cSBarry Smith { 232f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 233f0b9ad6cSBarry Smith 234f0b9ad6cSBarry Smith PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y)); 236f0b9ad6cSBarry Smith PetscFunctionReturn(0); 237f0b9ad6cSBarry Smith } 238f0b9ad6cSBarry Smith 239f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 240f0b9ad6cSBarry Smith { 241f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 242f0b9ad6cSBarry Smith 243f0b9ad6cSBarry Smith PetscFunctionBegin; 2449566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y)); 245f0b9ad6cSBarry Smith PetscFunctionReturn(0); 246f0b9ad6cSBarry Smith } 247f0b9ad6cSBarry Smith 2489b54502bSHong Zhang /*MC 2499b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2509b54502bSHong Zhang 2519b54502bSHong Zhang Options Database Keys: 2522401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2532401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2549b54502bSHong Zhang its factorization (overwrites original matrix) 2552401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2562401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 25755ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2582401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 2599b54502bSHong Zhang this decreases the chance of getting a zero pivot 2602401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 261967c93d3SBarry Smith - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 2629b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 2639b54502bSHong Zhang stability of the ILU factorization 2649b54502bSHong Zhang 2659b54502bSHong Zhang Level: beginner 2669b54502bSHong Zhang 26795452b02SPatrick Sanan Notes: 26895452b02SPatrick Sanan Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 2699b54502bSHong Zhang 2709b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 2719b54502bSHong Zhang 272f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 273f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 274f0b9ad6cSBarry Smith 2755ea7661aSPierre Jolivet If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARSE matrices with block Jacobi), factorization 276cd0a26f6SPaul Mullowney is never done on the GPU). 277ac793be5SBarry Smith 278c582cd25SBarry Smith References: 279606c0280SSatish Balay + * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 28096a0c994SBarry Smith self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 281606c0280SSatish Balay . * - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 282606c0280SSatish Balay - * - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 28396a0c994SBarry Smith Chapter in Parallel Numerical 284c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 28596a0c994SBarry Smith Science and Engineering, Kluwer. 286c582cd25SBarry Smith 2879b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 288145b9266SHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 289b7c853c4SBarry Smith PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 29092e9c092SBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 29192e9c092SBarry Smith PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace() 2929b54502bSHong Zhang 2939b54502bSHong Zhang M*/ 2949b54502bSHong Zhang 2958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 2969b54502bSHong Zhang { 2979b54502bSHong Zhang PC_ILU *ilu; 2989b54502bSHong Zhang 2999b54502bSHong Zhang PetscFunctionBegin; 3009566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&ilu)); 3013d1c1ea0SBarry Smith pc->data = (void*)ilu; 3029566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc,MAT_FACTOR_ILU)); 3039b54502bSHong Zhang 30475567043SBarry Smith ((PC_Factor*)ilu)->info.levels = 0.; 305075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 3060a545947SLisandro Dalcin ilu->col = NULL; 3070a545947SLisandro Dalcin ilu->row = NULL; 308075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 309075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 310075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 3119b54502bSHong Zhang 312574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3139b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3149b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3157b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ILU; 3169b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3179b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3189b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 31992e08861SBarry Smith pc->ops->view = PCView_Factor; 320f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 321f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3220a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU)); 3249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU)); 3259b54502bSHong Zhang PetscFunctionReturn(0); 3269b54502bSHong Zhang } 327