xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
19b54502bSHong Zhang /*
29b54502bSHong Zhang    Defines a ILU factorization preconditioner for any Mat implementation
39b54502bSHong Zhang */
4c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/ilu/ilu.h> /*I "petscpc.h"  I*/
59b54502bSHong Zhang 
666976f2fSJacob Faibussowitsch static PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc, PetscReal z)
7d71ae5a4SJacob Faibussowitsch {
89b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
99b54502bSHong Zhang 
109b54502bSHong Zhang   PetscFunctionBegin;
119b54502bSHong Zhang   ilu->nonzerosalongdiagonal = PETSC_TRUE;
1213bcc0bdSJacob Faibussowitsch   if (z == (PetscReal)PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
132fa5cd67SKarl Rupp   else ilu->nonzerosalongdiagonaltol = z;
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159b54502bSHong Zhang }
169b54502bSHong Zhang 
1766976f2fSJacob Faibussowitsch static PetscErrorCode PCReset_ILU(PC pc)
18d71ae5a4SJacob Faibussowitsch {
199b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
209b54502bSHong Zhang 
219b54502bSHong Zhang   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   if (!ilu->hdr.inplace) PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
239566063dSJacob Faibussowitsch   if (ilu->row && ilu->col && ilu->row != ilu->col) PetscCall(ISDestroy(&ilu->row));
249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&ilu->col));
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
269b54502bSHong Zhang }
279b54502bSHong Zhang 
28d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount)
29d71ae5a4SJacob Faibussowitsch {
30075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU *)pc->data;
319b54502bSHong Zhang 
329b54502bSHong Zhang   PetscFunctionBegin;
334c9036c7SBarry Smith   if (pc->setupcalled && (((PC_Factor *)ilu)->info.dt != dt || ((PC_Factor *)ilu)->info.dtcol != dtcol || ((PC_Factor *)ilu)->info.dtcount != dtcount)) {
34ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot change drop tolerance after using PC");
359b54502bSHong Zhang   }
36075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = dt;
37075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = dtcol;
38075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = dtcount;
394c9036c7SBarry Smith   ((PC_Factor *)ilu)->info.usedt   = 1.0;
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
419b54502bSHong Zhang }
429b54502bSHong Zhang 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ILU(PC pc, PetscOptionItems *PetscOptionsObject)
44d71ae5a4SJacob Faibussowitsch {
4578fc6b22SHong Zhang   PetscInt  itmp;
468afaa268SBarry Smith   PetscBool flg, set;
479b54502bSHong Zhang   PC_ILU   *ilu = (PC_ILU *)pc->data;
489b54502bSHong Zhang   PetscReal tol;
499b54502bSHong Zhang 
509b54502bSHong Zhang   PetscFunctionBegin;
51d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ILU Options");
52dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
538ff23777SHong Zhang 
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_factor_levels", "levels of fill", "PCFactorSetLevels", (PetscInt)((PC_Factor *)ilu)->info.levels, &itmp, &flg));
55075768bcSBarry Smith   if (flg) ((PC_Factor *)ilu)->info.levels = itmp;
562fa5cd67SKarl Rupp 
579566063dSJacob 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));
588afaa268SBarry Smith   if (set) ((PC_Factor *)ilu)->info.diagonal_fill = (PetscReal)flg;
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
609b54502bSHong Zhang   if (flg) {
619b54502bSHong Zhang     tol = PETSC_DECIDE;
629566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", ilu->nonzerosalongdiagonaltol, &tol, NULL));
639566063dSJacob Faibussowitsch     PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
649b54502bSHong Zhang   }
659b54502bSHong Zhang 
66d0609cedSBarry Smith   PetscOptionsHeadEnd();
673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
689b54502bSHong Zhang }
699b54502bSHong Zhang 
70d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ILU(PC pc)
71d71ae5a4SJacob Faibussowitsch {
729b54502bSHong Zhang   PC_ILU        *ilu = (PC_ILU *)pc->data;
73f3a39becSBarry Smith   MatInfo        info;
74ace3abfcSBarry Smith   PetscBool      flg;
75ea799195SBarry Smith   MatSolverType  stype;
7600e125f8SBarry Smith   MatFactorError err;
77f023e1d5SPierre Jolivet   const char    *prefix;
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 
9326cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
9426cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
9526cc229bSBarry Smith 
969566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
973d1c1ea0SBarry Smith   if (ilu->hdr.inplace) {
989b54502bSHong Zhang     if (!pc->setupcalled) {
999b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
1009b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
1012c7c0729SBarry Smith       /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */
1029566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1039566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
1049566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &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;
1163ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
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;
12603e5aca4SStefano Zampini 
12703e5aca4SStefano Zampini       PetscCall(PCFactorSetUpMatSolverType(pc));
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;
14103e5aca4SStefano Zampini 
1429566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
14303e5aca4SStefano Zampini         PetscCall(PCFactorSetUpMatSolverType(pc));
1449566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
145f73b0415SBarry Smith         if (canuseordering) {
1469b54502bSHong Zhang           /* compute a new ordering for the ILU */
1479566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->row));
1489566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->col));
1499566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1509566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1519b54502bSHong Zhang           /*  Remove zeros along diagonal?     */
1521baa6e33SBarry Smith           if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
1539b54502bSHong Zhang         }
1542c7c0729SBarry Smith       }
1559566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1569566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1573d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1589b54502bSHong Zhang     }
1599566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16000e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
16100e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1623ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1636baea169SHong Zhang     }
1646baea169SHong Zhang 
1659566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)ilu)->fact, pc->pmat, &((PC_Factor *)ilu)->info));
1669566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16700e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
16800e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1696baea169SHong Zhang     }
1709b54502bSHong Zhang   }
17100c67f3bSHong Zhang 
1729566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
17300c67f3bSHong Zhang   if (!stype) {
174ea799195SBarry Smith     MatSolverType solverpackage;
1759566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)ilu)->fact, &solverpackage));
1769566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
17700c67f3bSHong Zhang   }
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1799b54502bSHong Zhang }
1809b54502bSHong Zhang 
181d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ILU(PC pc)
182d71ae5a4SJacob Faibussowitsch {
1839b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1849b54502bSHong Zhang 
1859b54502bSHong Zhang   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(PCReset_ILU(pc));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->solvertype));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->ordering));
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1902e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1929b54502bSHong Zhang }
1939b54502bSHong Zhang 
194d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ILU(PC pc, Vec x, Vec y)
195d71ae5a4SJacob Faibussowitsch {
1969b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1979b54502bSHong Zhang 
1989b54502bSHong Zhang   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)ilu)->fact, x, y));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2019b54502bSHong Zhang }
2029b54502bSHong Zhang 
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ILU(PC pc, Mat X, Mat Y)
204d71ae5a4SJacob Faibussowitsch {
2057b6e2003SPierre Jolivet   PC_ILU *ilu = (PC_ILU *)pc->data;
2067b6e2003SPierre Jolivet 
2077b6e2003SPierre Jolivet   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)ilu)->fact, X, Y));
2093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2107b6e2003SPierre Jolivet }
2117b6e2003SPierre Jolivet 
212d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ILU(PC pc, Vec x, Vec y)
213d71ae5a4SJacob Faibussowitsch {
2149b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
2159b54502bSHong Zhang 
2169b54502bSHong Zhang   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(((PC_Factor *)ilu)->fact, x, y));
2183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2199b54502bSHong Zhang }
2209b54502bSHong Zhang 
221d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc, Vec x, Vec y)
222d71ae5a4SJacob Faibussowitsch {
223f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
224f0b9ad6cSBarry Smith 
225f0b9ad6cSBarry Smith   PetscFunctionBegin;
2269566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228f0b9ad6cSBarry Smith }
229f0b9ad6cSBarry Smith 
230d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ILU(PC pc, Vec x, Vec y)
231d71ae5a4SJacob Faibussowitsch {
232f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
233f0b9ad6cSBarry Smith 
234f0b9ad6cSBarry Smith   PetscFunctionBegin;
2359566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
237f0b9ad6cSBarry Smith }
238f0b9ad6cSBarry Smith 
2399b54502bSHong Zhang /*MC
2409b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
2419b54502bSHong Zhang 
2429b54502bSHong Zhang    Options Database Keys:
2432401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
2442401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
2459b54502bSHong Zhang                       its factorization (overwrites original matrix)
2462401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
2472401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
24855ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2492401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
2509b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
2512401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
252967c93d3SBarry Smith -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
2539b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
2549b54502bSHong Zhang                              stability of the ILU factorization
2559b54502bSHong Zhang 
2569b54502bSHong Zhang    Level: beginner
2579b54502bSHong Zhang 
25895452b02SPatrick Sanan    Notes:
259f1580f4eSBarry Smith    Only implemented for some matrix format and sequential. For parallel see `PCHYPRE` for hypre's ILU
2609b54502bSHong Zhang 
261f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ILU
2629b54502bSHong Zhang 
263f0b9ad6cSBarry Smith    The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
264f0b9ad6cSBarry Smith    even when the matrix is not symmetric since the U stores the diagonals of the factorization.
265f0b9ad6cSBarry Smith 
266f1580f4eSBarry Smith    If you are using `MATSEQAIJCUSPARSE` matrices (or `MATMPIAIJCUSPARSE` matrices with block Jacobi), factorization
267cd0a26f6SPaul Mullowney    is never done on the GPU).
268ac793be5SBarry Smith 
269c582cd25SBarry Smith    References:
270606c0280SSatish Balay +  * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
27196a0c994SBarry Smith    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
272606c0280SSatish Balay .  * -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
273606c0280SSatish Balay -  * -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
27496a0c994SBarry Smith       Chapter in Parallel Numerical
275c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
27696a0c994SBarry Smith       Science and Engineering, Kluwer.
277c582cd25SBarry Smith 
278*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCLU`, `PCICC`, `PCCHOLESKY`,
279db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`,
280c2e3fba1SPatrick Sanan           `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
281db781477SPatrick Sanan           `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`,
282db781477SPatrick Sanan           `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()`
2839b54502bSHong Zhang M*/
2849b54502bSHong Zhang 
285d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
286d71ae5a4SJacob Faibussowitsch {
2879b54502bSHong Zhang   PC_ILU *ilu;
2889b54502bSHong Zhang 
2899b54502bSHong Zhang   PetscFunctionBegin;
2904dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilu));
2913d1c1ea0SBarry Smith   pc->data = (void *)ilu;
2929566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ILU));
2939b54502bSHong Zhang 
29475567043SBarry Smith   ((PC_Factor *)ilu)->info.levels  = 0.;
295075768bcSBarry Smith   ((PC_Factor *)ilu)->info.fill    = 1.0;
2960a545947SLisandro Dalcin   ilu->col                         = NULL;
2970a545947SLisandro Dalcin   ilu->row                         = NULL;
298075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = PETSC_DEFAULT;
299075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = PETSC_DEFAULT;
300075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = PETSC_DEFAULT;
3019b54502bSHong Zhang 
302574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
3039b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3049b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
3057b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ILU;
3069b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
3079b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
3089b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
30992e08861SBarry Smith   pc->ops->view                = PCView_Factor;
310f0b9ad6cSBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
311f0b9ad6cSBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
3120a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", PCFactorSetDropTolerance_ILU));
3149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_ILU));
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3169b54502bSHong Zhang }
317