xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
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 
79371c9d4SSatish Balay PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc, PetscReal z) {
89b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
99b54502bSHong Zhang 
109b54502bSHong Zhang   PetscFunctionBegin;
119b54502bSHong Zhang   ilu->nonzerosalongdiagonal = PETSC_TRUE;
122fa5cd67SKarl Rupp   if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
132fa5cd67SKarl Rupp   else ilu->nonzerosalongdiagonaltol = z;
149b54502bSHong Zhang   PetscFunctionReturn(0);
159b54502bSHong Zhang }
169b54502bSHong Zhang 
179371c9d4SSatish Balay PetscErrorCode PCReset_ILU(PC pc) {
189b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
199b54502bSHong Zhang 
209b54502bSHong Zhang   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   if (!ilu->hdr.inplace) PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
229566063dSJacob Faibussowitsch   if (ilu->row && ilu->col && ilu->row != ilu->col) PetscCall(ISDestroy(&ilu->row));
239566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ilu->col));
249b54502bSHong Zhang   PetscFunctionReturn(0);
259b54502bSHong Zhang }
269b54502bSHong Zhang 
279371c9d4SSatish Balay PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount) {
28075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU *)pc->data;
299b54502bSHong Zhang 
309b54502bSHong Zhang   PetscFunctionBegin;
314c9036c7SBarry Smith   if (pc->setupcalled && (((PC_Factor *)ilu)->info.dt != dt || ((PC_Factor *)ilu)->info.dtcol != dtcol || ((PC_Factor *)ilu)->info.dtcount != dtcount)) {
32ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot change drop tolerance after using PC");
339b54502bSHong Zhang   }
34075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = dt;
35075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = dtcol;
36075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = dtcount;
374c9036c7SBarry Smith   ((PC_Factor *)ilu)->info.usedt   = 1.0;
389b54502bSHong Zhang   PetscFunctionReturn(0);
399b54502bSHong Zhang }
409b54502bSHong Zhang 
419371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_ILU(PC pc, PetscOptionItems *PetscOptionsObject) {
4278fc6b22SHong Zhang   PetscInt  itmp;
438afaa268SBarry Smith   PetscBool flg, set;
449b54502bSHong Zhang   PC_ILU   *ilu = (PC_ILU *)pc->data;
459b54502bSHong Zhang   PetscReal tol;
469b54502bSHong Zhang 
479b54502bSHong Zhang   PetscFunctionBegin;
48d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ILU Options");
49dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
508ff23777SHong Zhang 
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_factor_levels", "levels of fill", "PCFactorSetLevels", (PetscInt)((PC_Factor *)ilu)->info.levels, &itmp, &flg));
52075768bcSBarry Smith   if (flg) ((PC_Factor *)ilu)->info.levels = itmp;
532fa5cd67SKarl Rupp 
549566063dSJacob 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));
558afaa268SBarry Smith   if (set) ((PC_Factor *)ilu)->info.diagonal_fill = (PetscReal)flg;
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
579b54502bSHong Zhang   if (flg) {
589b54502bSHong Zhang     tol = PETSC_DECIDE;
599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", ilu->nonzerosalongdiagonaltol, &tol, NULL));
609566063dSJacob Faibussowitsch     PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
619b54502bSHong Zhang   }
629b54502bSHong Zhang 
63d0609cedSBarry Smith   PetscOptionsHeadEnd();
649b54502bSHong Zhang   PetscFunctionReturn(0);
659b54502bSHong Zhang }
669b54502bSHong Zhang 
679371c9d4SSatish Balay static PetscErrorCode PCSetUp_ILU(PC pc) {
689b54502bSHong Zhang   PC_ILU        *ilu = (PC_ILU *)pc->data;
69f3a39becSBarry Smith   MatInfo        info;
70ace3abfcSBarry Smith   PetscBool      flg;
71ea799195SBarry Smith   MatSolverType  stype;
7200e125f8SBarry Smith   MatFactorError err;
73f023e1d5SPierre Jolivet   const char    *prefix;
749b54502bSHong Zhang 
759b54502bSHong Zhang   PetscFunctionBegin;
76c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
7792927226SBarry Smith   /* ugly hack to change default, since it is not support by some matrix types */
7892927226SBarry Smith   if (((PC_Factor *)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
799566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg));
8092927226SBarry Smith     if (!flg) {
819566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATMPIAIJ, &flg));
8292927226SBarry Smith       if (!flg) {
8392927226SBarry Smith         ((PC_Factor *)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
849566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n"));
8592927226SBarry Smith       }
8692927226SBarry Smith     }
8792927226SBarry Smith   }
8892927226SBarry Smith 
8926cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
9026cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
9126cc229bSBarry Smith 
929566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
933d1c1ea0SBarry Smith   if (ilu->hdr.inplace) {
949b54502bSHong Zhang     if (!pc->setupcalled) {
959b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
969b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
972c7c0729SBarry Smith       /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */
989566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
999566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
1009566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1019566063dSJacob Faibussowitsch       if (ilu->row) PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->row));
1029566063dSJacob Faibussowitsch       if (ilu->col) PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->col));
1039b54502bSHong Zhang     }
1049b54502bSHong Zhang 
1059b54502bSHong Zhang     /* In place ILU only makes sense with fill factor of 1.0 because
1069b54502bSHong Zhang        cannot have levels of fill */
107075768bcSBarry Smith     ((PC_Factor *)ilu)->info.fill          = 1.0;
10875567043SBarry Smith     ((PC_Factor *)ilu)->info.diagonal_fill = 0.0;
1092fa5cd67SKarl Rupp 
1109566063dSJacob Faibussowitsch     PetscCall(MatILUFactor(pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1119566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(pc->pmat, &err));
11200e125f8SBarry Smith     if (err) { /* Factor() fails */
11300e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1146baea169SHong Zhang       PetscFunctionReturn(0);
1156baea169SHong Zhang     }
1166baea169SHong Zhang 
117075768bcSBarry Smith     ((PC_Factor *)ilu)->fact = pc->pmat;
118b33856dcSBarry Smith     /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
1199566063dSJacob Faibussowitsch     PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &pc->matstate));
1209b54502bSHong Zhang   } else {
1219b54502bSHong Zhang     if (!pc->setupcalled) {
1229b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
123f73b0415SBarry Smith       PetscBool canuseordering;
1242c7c0729SBarry Smith       if (!((PC_Factor *)ilu)->fact) {
1259566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)ilu)->solvertype, MAT_FACTOR_ILU, &((PC_Factor *)ilu)->fact));
1269566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)ilu)->fact));
1272c7c0729SBarry Smith       }
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));
1329566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->row));
1339566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->col));
1349b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
1351baa6e33SBarry Smith         if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
136a1f19f5aSHong Zhang       }
1379566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1389566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1393d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1409b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1413d1c1ea0SBarry Smith       if (!ilu->hdr.reuseordering) {
142f73b0415SBarry Smith         PetscBool canuseordering;
1439566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
1449566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)ilu)->solvertype, MAT_FACTOR_ILU, &((PC_Factor *)ilu)->fact));
1459566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)ilu)->fact));
1469566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
147f73b0415SBarry Smith         if (canuseordering) {
1489b54502bSHong Zhang           /* compute a new ordering for the ILU */
1499566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->row));
1509566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->col));
1519566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1529566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1539566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->row));
1549566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ilu->col));
1559b54502bSHong Zhang           /*  Remove zeros along diagonal?     */
1561baa6e33SBarry Smith           if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
1579b54502bSHong Zhang         }
1582c7c0729SBarry Smith       }
1599566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1609566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1613d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1629b54502bSHong Zhang     }
1639566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16400e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
16500e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1666baea169SHong Zhang       PetscFunctionReturn(0);
1676baea169SHong Zhang     }
1686baea169SHong Zhang 
1699566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)ilu)->fact, pc->pmat, &((PC_Factor *)ilu)->info));
1709566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
17100e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
17200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1736baea169SHong Zhang     }
1749b54502bSHong Zhang   }
17500c67f3bSHong Zhang 
1769566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
17700c67f3bSHong Zhang   if (!stype) {
178ea799195SBarry Smith     MatSolverType solverpackage;
1799566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)ilu)->fact, &solverpackage));
1809566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
18100c67f3bSHong Zhang   }
1829b54502bSHong Zhang   PetscFunctionReturn(0);
1839b54502bSHong Zhang }
1849b54502bSHong Zhang 
1859371c9d4SSatish Balay static PetscErrorCode PCDestroy_ILU(PC pc) {
1869b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1879b54502bSHong Zhang 
1889b54502bSHong Zhang   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(PCReset_ILU(pc));
1909566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->solvertype));
1919566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->ordering));
1929566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1932e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1949b54502bSHong Zhang   PetscFunctionReturn(0);
1959b54502bSHong Zhang }
1969b54502bSHong Zhang 
1979371c9d4SSatish Balay static PetscErrorCode PCApply_ILU(PC pc, Vec x, Vec y) {
1989b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1999b54502bSHong Zhang 
2009b54502bSHong Zhang   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)ilu)->fact, x, y));
2029b54502bSHong Zhang   PetscFunctionReturn(0);
2039b54502bSHong Zhang }
2049b54502bSHong Zhang 
2059371c9d4SSatish Balay static PetscErrorCode PCMatApply_ILU(PC pc, Mat X, Mat Y) {
2067b6e2003SPierre Jolivet   PC_ILU *ilu = (PC_ILU *)pc->data;
2077b6e2003SPierre Jolivet 
2087b6e2003SPierre Jolivet   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)ilu)->fact, X, Y));
2107b6e2003SPierre Jolivet   PetscFunctionReturn(0);
2117b6e2003SPierre Jolivet }
2127b6e2003SPierre Jolivet 
2139371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_ILU(PC pc, Vec x, Vec y) {
2149b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
2159b54502bSHong Zhang 
2169b54502bSHong Zhang   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(((PC_Factor *)ilu)->fact, x, y));
2189b54502bSHong Zhang   PetscFunctionReturn(0);
2199b54502bSHong Zhang }
2209b54502bSHong Zhang 
2219371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc, Vec x, Vec y) {
222f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
223f0b9ad6cSBarry Smith 
224f0b9ad6cSBarry Smith   PetscFunctionBegin;
2259566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
226f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
227f0b9ad6cSBarry Smith }
228f0b9ad6cSBarry Smith 
2299371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricRight_ILU(PC pc, Vec x, Vec y) {
230f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
231f0b9ad6cSBarry Smith 
232f0b9ad6cSBarry Smith   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
234f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
235f0b9ad6cSBarry Smith }
236f0b9ad6cSBarry Smith 
2379b54502bSHong Zhang /*MC
2389b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
2399b54502bSHong Zhang 
2409b54502bSHong Zhang    Options Database Keys:
2412401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
2422401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
2439b54502bSHong Zhang                       its factorization (overwrites original matrix)
2442401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
2452401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
24655ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2472401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
2489b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
2492401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
250967c93d3SBarry Smith -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
2519b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
2529b54502bSHong Zhang                              stability of the ILU factorization
2539b54502bSHong Zhang 
2549b54502bSHong Zhang    Level: beginner
2559b54502bSHong Zhang 
25695452b02SPatrick Sanan    Notes:
257*f1580f4eSBarry Smith    Only implemented for some matrix format and sequential. For parallel see `PCHYPRE` for hypre's ILU
2589b54502bSHong Zhang 
259*f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ILU
2609b54502bSHong Zhang 
261f0b9ad6cSBarry Smith    The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
262f0b9ad6cSBarry Smith    even when the matrix is not symmetric since the U stores the diagonals of the factorization.
263f0b9ad6cSBarry Smith 
264*f1580f4eSBarry Smith    If you are using `MATSEQAIJCUSPARSE` matrices (or `MATMPIAIJCUSPARSE` matrices with block Jacobi), factorization
265cd0a26f6SPaul Mullowney    is never done on the GPU).
266ac793be5SBarry Smith 
267c582cd25SBarry Smith    References:
268606c0280SSatish Balay +  * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
26996a0c994SBarry Smith    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
270606c0280SSatish Balay .  * -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
271606c0280SSatish Balay -  * -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
27296a0c994SBarry Smith       Chapter in Parallel Numerical
273c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
27496a0c994SBarry Smith       Science and Engineering, Kluwer.
275c582cd25SBarry Smith 
276*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCLU`, `PCICC`, `PCCHOLESKY`,
277db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`,
278c2e3fba1SPatrick Sanan           `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
279db781477SPatrick Sanan           `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`,
280db781477SPatrick Sanan           `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()`
2819b54502bSHong Zhang M*/
2829b54502bSHong Zhang 
2839371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) {
2849b54502bSHong Zhang   PC_ILU *ilu;
2859b54502bSHong Zhang 
2869b54502bSHong Zhang   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &ilu));
2883d1c1ea0SBarry Smith   pc->data = (void *)ilu;
2899566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ILU));
2909b54502bSHong Zhang 
29175567043SBarry Smith   ((PC_Factor *)ilu)->info.levels  = 0.;
292075768bcSBarry Smith   ((PC_Factor *)ilu)->info.fill    = 1.0;
2930a545947SLisandro Dalcin   ilu->col                         = NULL;
2940a545947SLisandro Dalcin   ilu->row                         = NULL;
295075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = PETSC_DEFAULT;
296075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = PETSC_DEFAULT;
297075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = PETSC_DEFAULT;
2989b54502bSHong Zhang 
299574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
3009b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3019b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
3027b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ILU;
3039b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
3049b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
3059b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
30692e08861SBarry Smith   pc->ops->view                = PCView_Factor;
307f0b9ad6cSBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
308f0b9ad6cSBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
3090a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", PCFactorSetDropTolerance_ILU));
3119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_ILU));
3129b54502bSHong Zhang   PetscFunctionReturn(0);
3139b54502bSHong Zhang }
314