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; 52d0609cedSBarry 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 67d0609cedSBarry 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; 78f023e1d5SPierre Jolivet const char *prefix; 799b54502bSHong Zhang 809b54502bSHong Zhang PetscFunctionBegin; 81f023e1d5SPierre Jolivet PetscCall(PCGetOptionsPrefix(pc,&prefix)); 82f023e1d5SPierre Jolivet PetscCall(MatSetOptionsPrefix(pc->pmat,prefix)); 83c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 8492927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 8592927226SBarry Smith if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 869566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg)); 8792927226SBarry Smith if (!flg) { 889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg)); 8992927226SBarry Smith if (!flg) { 9092927226SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 919566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n")); 9292927226SBarry Smith } 9392927226SBarry Smith } 9492927226SBarry Smith } 9592927226SBarry Smith 969566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure)); 973d1c1ea0SBarry Smith if (ilu->hdr.inplace) { 989b54502bSHong Zhang if (!pc->setupcalled) { 999b54502bSHong Zhang 1009b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 1019b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 1022c7c0729SBarry Smith /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */ 1039566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact)); 1059566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1069566063dSJacob Faibussowitsch if (ilu->row) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1079566063dSJacob Faibussowitsch if (ilu->col) PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1089b54502bSHong Zhang } 1099b54502bSHong Zhang 1109b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1119b54502bSHong Zhang cannot have levels of fill */ 112075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 11375567043SBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 1142fa5cd67SKarl Rupp 1159566063dSJacob Faibussowitsch PetscCall(MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1169566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat,&err)); 11700e125f8SBarry Smith if (err) { /* Factor() fails */ 11800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1196baea169SHong Zhang PetscFunctionReturn(0); 1206baea169SHong Zhang } 1216baea169SHong Zhang 122075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 123b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 1249566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate)); 1259b54502bSHong Zhang } else { 1269b54502bSHong Zhang if (!pc->setupcalled) { 1279b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 128f73b0415SBarry Smith PetscBool canuseordering; 1292c7c0729SBarry Smith if (!((PC_Factor*)ilu)->fact) { 1309566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact)); 1319566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact)); 1322c7c0729SBarry Smith } 1339566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering)); 134f73b0415SBarry Smith if (canuseordering) { 1359566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1369566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1379566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1389566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1399b54502bSHong Zhang /* Remove zeros along diagonal? */ 1409b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1419566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col)); 1429b54502bSHong Zhang } 143a1f19f5aSHong Zhang } 1449566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1459566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info)); 1463d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1479b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1483d1c1ea0SBarry Smith if (!ilu->hdr.reuseordering) { 149f73b0415SBarry Smith PetscBool canuseordering; 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)ilu)->fact)); 1519566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact)); 1529566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact)); 1539566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering)); 154f73b0415SBarry Smith if (canuseordering) { 1559b54502bSHong Zhang /* compute a new ordering for the ILU */ 1569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->row)); 1579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->col)); 1589566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1599566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col)); 1609566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row)); 1619566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col)); 1629b54502bSHong Zhang /* Remove zeros along diagonal? */ 1639b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1649566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col)); 1659b54502bSHong Zhang } 1669b54502bSHong Zhang } 1672c7c0729SBarry Smith } 1689566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info)); 1699566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info)); 1703d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1719b54502bSHong Zhang } 1729566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err)); 17300e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 17400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1756baea169SHong Zhang PetscFunctionReturn(0); 1766baea169SHong Zhang } 1776baea169SHong Zhang 1789566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info)); 1799566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)ilu)->fact,&err)); 18000e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 18100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1826baea169SHong Zhang } 1839b54502bSHong Zhang } 18400c67f3bSHong Zhang 1859566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 18600c67f3bSHong Zhang if (!stype) { 187ea799195SBarry Smith MatSolverType solverpackage; 1889566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage)); 1899566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 19000c67f3bSHong Zhang } 1919b54502bSHong Zhang PetscFunctionReturn(0); 1929b54502bSHong Zhang } 1939b54502bSHong Zhang 1949b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 1959b54502bSHong Zhang { 1969b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 1979b54502bSHong Zhang 1989b54502bSHong Zhang PetscFunctionBegin; 1999566063dSJacob Faibussowitsch PetscCall(PCReset_ILU(pc)); 2009566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)ilu)->solvertype)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)ilu)->ordering)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2039b54502bSHong Zhang PetscFunctionReturn(0); 2049b54502bSHong Zhang } 2059b54502bSHong Zhang 2069b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 2079b54502bSHong Zhang { 2089b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2099b54502bSHong Zhang 2109b54502bSHong Zhang PetscFunctionBegin; 2119566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor*)ilu)->fact,x,y)); 2129b54502bSHong Zhang PetscFunctionReturn(0); 2139b54502bSHong Zhang } 2149b54502bSHong Zhang 2157b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ILU(PC pc,Mat X,Mat Y) 2167b6e2003SPierre Jolivet { 2177b6e2003SPierre Jolivet PC_ILU *ilu = (PC_ILU*)pc->data; 2187b6e2003SPierre Jolivet 2197b6e2003SPierre Jolivet PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor*)ilu)->fact,X,Y)); 2217b6e2003SPierre Jolivet PetscFunctionReturn(0); 2227b6e2003SPierre Jolivet } 2237b6e2003SPierre Jolivet 2249b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 2259b54502bSHong Zhang { 2269b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2279b54502bSHong Zhang 2289b54502bSHong Zhang PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y)); 2309b54502bSHong Zhang PetscFunctionReturn(0); 2319b54502bSHong Zhang } 2329b54502bSHong Zhang 233f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 234f0b9ad6cSBarry Smith { 235f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 236f0b9ad6cSBarry Smith 237f0b9ad6cSBarry Smith PetscFunctionBegin; 2389566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y)); 239f0b9ad6cSBarry Smith PetscFunctionReturn(0); 240f0b9ad6cSBarry Smith } 241f0b9ad6cSBarry Smith 242f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 243f0b9ad6cSBarry Smith { 244f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 245f0b9ad6cSBarry Smith 246f0b9ad6cSBarry Smith PetscFunctionBegin; 2479566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y)); 248f0b9ad6cSBarry Smith PetscFunctionReturn(0); 249f0b9ad6cSBarry Smith } 250f0b9ad6cSBarry Smith 2519b54502bSHong Zhang /*MC 2529b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2539b54502bSHong Zhang 2549b54502bSHong Zhang Options Database Keys: 2552401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2562401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2579b54502bSHong Zhang its factorization (overwrites original matrix) 2582401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2592401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 26055ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2612401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 2629b54502bSHong Zhang this decreases the chance of getting a zero pivot 2632401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 264967c93d3SBarry Smith - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 2659b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 2669b54502bSHong Zhang stability of the ILU factorization 2679b54502bSHong Zhang 2689b54502bSHong Zhang Level: beginner 2699b54502bSHong Zhang 27095452b02SPatrick Sanan Notes: 27195452b02SPatrick Sanan Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 2729b54502bSHong Zhang 2739b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 2749b54502bSHong Zhang 275f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 276f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 277f0b9ad6cSBarry Smith 2785ea7661aSPierre Jolivet If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARSE matrices with block Jacobi), factorization 279cd0a26f6SPaul Mullowney is never done on the GPU). 280ac793be5SBarry Smith 281c582cd25SBarry Smith References: 282606c0280SSatish Balay + * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 28396a0c994SBarry Smith self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 284606c0280SSatish Balay . * - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 285606c0280SSatish Balay - * - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 28696a0c994SBarry Smith Chapter in Parallel Numerical 287c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 28896a0c994SBarry Smith Science and Engineering, Kluwer. 289c582cd25SBarry Smith 290db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, 291db781477SPatrick Sanan `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`, 292*c2e3fba1SPatrick Sanan `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`, 293db781477SPatrick Sanan `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`, 294db781477SPatrick Sanan `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()` 2959b54502bSHong Zhang 2969b54502bSHong Zhang M*/ 2979b54502bSHong Zhang 2988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 2999b54502bSHong Zhang { 3009b54502bSHong Zhang PC_ILU *ilu; 3019b54502bSHong Zhang 3029b54502bSHong Zhang PetscFunctionBegin; 3039566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&ilu)); 3043d1c1ea0SBarry Smith pc->data = (void*)ilu; 3059566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc,MAT_FACTOR_ILU)); 3069b54502bSHong Zhang 30775567043SBarry Smith ((PC_Factor*)ilu)->info.levels = 0.; 308075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 3090a545947SLisandro Dalcin ilu->col = NULL; 3100a545947SLisandro Dalcin ilu->row = NULL; 311075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 312075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 313075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 3149b54502bSHong Zhang 315574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3169b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3179b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3187b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ILU; 3199b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3209b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3219b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 32292e08861SBarry Smith pc->ops->view = PCView_Factor; 323f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 324f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3250a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU)); 3279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU)); 3289b54502bSHong Zhang PetscFunctionReturn(0); 3299b54502bSHong Zhang } 330