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 7d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc, PetscReal z) 8d71ae5a4SJacob Faibussowitsch { 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; 15*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 169b54502bSHong Zhang } 179b54502bSHong Zhang 18d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_ILU(PC pc) 19d71ae5a4SJacob Faibussowitsch { 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)); 26*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 279b54502bSHong Zhang } 289b54502bSHong Zhang 29d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount) 30d71ae5a4SJacob Faibussowitsch { 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; 41*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429b54502bSHong Zhang } 439b54502bSHong Zhang 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ILU(PC pc, PetscOptionItems *PetscOptionsObject) 45d71ae5a4SJacob Faibussowitsch { 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"); 53dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject)); 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(); 68*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 699b54502bSHong Zhang } 709b54502bSHong Zhang 71d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ILU(PC pc) 72d71ae5a4SJacob Faibussowitsch { 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; 81c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 8292927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 8392927226SBarry Smith if (((PC_Factor *)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 849566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg)); 8592927226SBarry Smith if (!flg) { 869566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATMPIAIJ, &flg)); 8792927226SBarry Smith if (!flg) { 8892927226SBarry Smith ((PC_Factor *)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 899566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n")); 9092927226SBarry Smith } 9192927226SBarry Smith } 9292927226SBarry Smith } 9392927226SBarry Smith 9426cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 9526cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 9626cc229bSBarry Smith 979566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 983d1c1ea0SBarry Smith if (ilu->hdr.inplace) { 999b54502bSHong Zhang if (!pc->setupcalled) { 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)); 1069b54502bSHong Zhang } 1079b54502bSHong Zhang 1089b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1099b54502bSHong Zhang cannot have levels of fill */ 110075768bcSBarry Smith ((PC_Factor *)ilu)->info.fill = 1.0; 11175567043SBarry Smith ((PC_Factor *)ilu)->info.diagonal_fill = 0.0; 1122fa5cd67SKarl Rupp 1139566063dSJacob Faibussowitsch PetscCall(MatILUFactor(pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info)); 1149566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err)); 11500e125f8SBarry Smith if (err) { /* Factor() fails */ 11600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 117*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1186baea169SHong Zhang } 1196baea169SHong Zhang 120075768bcSBarry Smith ((PC_Factor *)ilu)->fact = pc->pmat; 121b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 1229566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &pc->matstate)); 1239b54502bSHong Zhang } else { 1249b54502bSHong Zhang if (!pc->setupcalled) { 1259b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 126f73b0415SBarry Smith PetscBool canuseordering; 1274dfa11a4SJacob Faibussowitsch if (!((PC_Factor *)ilu)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)ilu)->solvertype, MAT_FACTOR_ILU, &((PC_Factor *)ilu)->fact)); } 1289566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering)); 129f73b0415SBarry Smith if (canuseordering) { 1309566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1319566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col)); 1329b54502bSHong Zhang /* Remove zeros along diagonal? */ 1331baa6e33SBarry Smith if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col)); 134a1f19f5aSHong Zhang } 1359566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info)); 1369566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info)); 1373d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1389b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1393d1c1ea0SBarry Smith if (!ilu->hdr.reuseordering) { 140f73b0415SBarry Smith PetscBool canuseordering; 1419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact)); 1429566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)ilu)->solvertype, MAT_FACTOR_ILU, &((PC_Factor *)ilu)->fact)); 1439566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering)); 144f73b0415SBarry Smith if (canuseordering) { 1459b54502bSHong Zhang /* compute a new ordering for the ILU */ 1469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->row)); 1479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ilu->col)); 1489566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 1499566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col)); 1509b54502bSHong Zhang /* Remove zeros along diagonal? */ 1511baa6e33SBarry Smith if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col)); 1529b54502bSHong Zhang } 1532c7c0729SBarry Smith } 1549566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info)); 1559566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info)); 1563d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1579b54502bSHong Zhang } 1589566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err)); 15900e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 16000e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 161*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1626baea169SHong Zhang } 1636baea169SHong Zhang 1649566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor *)ilu)->fact, pc->pmat, &((PC_Factor *)ilu)->info)); 1659566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err)); 16600e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 16700e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1686baea169SHong Zhang } 1699b54502bSHong Zhang } 17000c67f3bSHong Zhang 1719566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 17200c67f3bSHong Zhang if (!stype) { 173ea799195SBarry Smith MatSolverType solverpackage; 1749566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)ilu)->fact, &solverpackage)); 1759566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 17600c67f3bSHong Zhang } 177*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1789b54502bSHong Zhang } 1799b54502bSHong Zhang 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ILU(PC pc) 181d71ae5a4SJacob Faibussowitsch { 1829b54502bSHong Zhang PC_ILU *ilu = (PC_ILU *)pc->data; 1839b54502bSHong Zhang 1849b54502bSHong Zhang PetscFunctionBegin; 1859566063dSJacob Faibussowitsch PetscCall(PCReset_ILU(pc)); 1869566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)ilu)->solvertype)); 1879566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)ilu)->ordering)); 1889566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1892e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 190*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1919b54502bSHong Zhang } 1929b54502bSHong Zhang 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ILU(PC pc, Vec x, Vec y) 194d71ae5a4SJacob Faibussowitsch { 1959b54502bSHong Zhang PC_ILU *ilu = (PC_ILU *)pc->data; 1969b54502bSHong Zhang 1979b54502bSHong Zhang PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)ilu)->fact, x, y)); 199*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2009b54502bSHong Zhang } 2019b54502bSHong Zhang 202d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ILU(PC pc, Mat X, Mat Y) 203d71ae5a4SJacob Faibussowitsch { 2047b6e2003SPierre Jolivet PC_ILU *ilu = (PC_ILU *)pc->data; 2057b6e2003SPierre Jolivet 2067b6e2003SPierre Jolivet PetscFunctionBegin; 2079566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)ilu)->fact, X, Y)); 208*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2097b6e2003SPierre Jolivet } 2107b6e2003SPierre Jolivet 211d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ILU(PC pc, Vec x, Vec y) 212d71ae5a4SJacob Faibussowitsch { 2139b54502bSHong Zhang PC_ILU *ilu = (PC_ILU *)pc->data; 2149b54502bSHong Zhang 2159b54502bSHong Zhang PetscFunctionBegin; 2169566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)ilu)->fact, x, y)); 217*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2189b54502bSHong Zhang } 2199b54502bSHong Zhang 220d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc, Vec x, Vec y) 221d71ae5a4SJacob Faibussowitsch { 222f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU *)pc->data; 223f0b9ad6cSBarry Smith 224f0b9ad6cSBarry Smith PetscFunctionBegin; 2259566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y)); 226*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227f0b9ad6cSBarry Smith } 228f0b9ad6cSBarry Smith 229d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ILU(PC pc, Vec x, Vec y) 230d71ae5a4SJacob Faibussowitsch { 231f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU *)pc->data; 232f0b9ad6cSBarry Smith 233f0b9ad6cSBarry Smith PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y)); 235*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236f0b9ad6cSBarry Smith } 237f0b9ad6cSBarry Smith 2389b54502bSHong Zhang /*MC 2399b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2409b54502bSHong Zhang 2419b54502bSHong Zhang Options Database Keys: 2422401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2432401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2449b54502bSHong Zhang its factorization (overwrites original matrix) 2452401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2462401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 24755ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2482401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 2499b54502bSHong Zhang this decreases the chance of getting a zero pivot 2502401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 251967c93d3SBarry Smith - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 2529b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 2539b54502bSHong Zhang stability of the ILU factorization 2549b54502bSHong Zhang 2559b54502bSHong Zhang Level: beginner 2569b54502bSHong Zhang 25795452b02SPatrick Sanan Notes: 258f1580f4eSBarry Smith Only implemented for some matrix format and sequential. For parallel see `PCHYPRE` for hypre's ILU 2599b54502bSHong Zhang 260f1580f4eSBarry Smith For `MATSEQBAIJ` matrices this implements a point block ILU 2619b54502bSHong Zhang 262f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 263f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 264f0b9ad6cSBarry Smith 265f1580f4eSBarry Smith If you are using `MATSEQAIJCUSPARSE` matrices (or `MATMPIAIJCUSPARSE` matrices with block Jacobi), factorization 266cd0a26f6SPaul Mullowney is never done on the GPU). 267ac793be5SBarry Smith 268c582cd25SBarry Smith References: 269606c0280SSatish Balay + * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 27096a0c994SBarry Smith self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 271606c0280SSatish Balay . * - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 272606c0280SSatish Balay - * - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 27396a0c994SBarry Smith Chapter in Parallel Numerical 274c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 27596a0c994SBarry Smith Science and Engineering, Kluwer. 276c582cd25SBarry Smith 277f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCLU`, `PCICC`, `PCCHOLESKY`, 278db781477SPatrick Sanan `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`, 279c2e3fba1SPatrick Sanan `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`, 280db781477SPatrick Sanan `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`, 281db781477SPatrick Sanan `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()` 2829b54502bSHong Zhang M*/ 2839b54502bSHong Zhang 284d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 285d71ae5a4SJacob Faibussowitsch { 2869b54502bSHong Zhang PC_ILU *ilu; 2879b54502bSHong Zhang 2889b54502bSHong Zhang PetscFunctionBegin; 2894dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ilu)); 2903d1c1ea0SBarry Smith pc->data = (void *)ilu; 2919566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ILU)); 2929b54502bSHong Zhang 29375567043SBarry Smith ((PC_Factor *)ilu)->info.levels = 0.; 294075768bcSBarry Smith ((PC_Factor *)ilu)->info.fill = 1.0; 2950a545947SLisandro Dalcin ilu->col = NULL; 2960a545947SLisandro Dalcin ilu->row = NULL; 297075768bcSBarry Smith ((PC_Factor *)ilu)->info.dt = PETSC_DEFAULT; 298075768bcSBarry Smith ((PC_Factor *)ilu)->info.dtcount = PETSC_DEFAULT; 299075768bcSBarry Smith ((PC_Factor *)ilu)->info.dtcol = PETSC_DEFAULT; 3009b54502bSHong Zhang 301574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3029b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3039b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3047b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ILU; 3059b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3069b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3079b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 30892e08861SBarry Smith pc->ops->view = PCView_Factor; 309f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 310f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3110a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", PCFactorSetDropTolerance_ILU)); 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_ILU)); 314*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3159b54502bSHong Zhang } 316