xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 901f93825bbaa60582af604c6700caf57884a2e1)
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 
PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)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 
PCReset_ILU(PC pc)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 
PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)28523895eeSPierre Jolivet static 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;
33966bd95aSPierre Jolivet   PetscCheck(!pc->setupcalled || !(((PC_Factor *)ilu)->info.dt != dt || ((PC_Factor *)ilu)->info.dtcol != dtcol || ((PC_Factor *)ilu)->info.dtcount != dtcount), PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot change drop tolerance after using PC");
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;
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
399b54502bSHong Zhang }
409b54502bSHong Zhang 
PCSetFromOptions_ILU(PC pc,PetscOptionItems PetscOptionsObject)41ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PC pc, PetscOptionItems PetscOptionsObject)
42d71ae5a4SJacob Faibussowitsch {
4378fc6b22SHong Zhang   PetscInt  itmp;
448afaa268SBarry Smith   PetscBool flg, set;
459b54502bSHong Zhang   PC_ILU   *ilu = (PC_ILU *)pc->data;
469b54502bSHong Zhang   PetscReal tol;
479b54502bSHong Zhang 
489b54502bSHong Zhang   PetscFunctionBegin;
49d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ILU Options");
50dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
518ff23777SHong Zhang 
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_factor_levels", "levels of fill", "PCFactorSetLevels", (PetscInt)((PC_Factor *)ilu)->info.levels, &itmp, &flg));
53075768bcSBarry Smith   if (flg) ((PC_Factor *)ilu)->info.levels = itmp;
542fa5cd67SKarl Rupp 
559566063dSJacob 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));
568afaa268SBarry Smith   if (set) ((PC_Factor *)ilu)->info.diagonal_fill = (PetscReal)flg;
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
589b54502bSHong Zhang   if (flg) {
599b54502bSHong Zhang     tol = PETSC_DECIDE;
609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", ilu->nonzerosalongdiagonaltol, &tol, NULL));
619566063dSJacob Faibussowitsch     PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
629b54502bSHong Zhang   }
639b54502bSHong Zhang 
64d0609cedSBarry Smith   PetscOptionsHeadEnd();
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
669b54502bSHong Zhang }
679b54502bSHong Zhang 
PCSetUp_ILU(PC pc)68d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ILU(PC pc)
69d71ae5a4SJacob Faibussowitsch {
709b54502bSHong Zhang   PC_ILU        *ilu = (PC_ILU *)pc->data;
71f3a39becSBarry Smith   MatInfo        info;
72ace3abfcSBarry Smith   PetscBool      flg;
73ea799195SBarry Smith   MatSolverType  stype;
7400e125f8SBarry Smith   MatFactorError err;
75f023e1d5SPierre Jolivet   const char    *prefix;
769b54502bSHong Zhang 
779b54502bSHong Zhang   PetscFunctionBegin;
78c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
7992927226SBarry Smith   /* ugly hack to change default, since it is not support by some matrix types */
8092927226SBarry Smith   if (((PC_Factor *)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
819566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg));
8292927226SBarry Smith     if (!flg) {
839566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATMPIAIJ, &flg));
8492927226SBarry Smith       if (!flg) {
8592927226SBarry Smith         ((PC_Factor *)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
869566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n"));
8792927226SBarry Smith       }
8892927226SBarry Smith     }
8992927226SBarry Smith   }
9092927226SBarry Smith 
9126cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
9226cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
9326cc229bSBarry Smith 
949566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
953d1c1ea0SBarry Smith   if (ilu->hdr.inplace) {
969b54502bSHong Zhang     if (!pc->setupcalled) {
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));
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;
1143ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
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;
12403e5aca4SStefano Zampini 
12503e5aca4SStefano Zampini       PetscCall(PCFactorSetUpMatSolverType(pc));
1269566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
127f73b0415SBarry Smith       if (canuseordering) {
1289566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1299566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1309b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
1311baa6e33SBarry Smith         if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
132a1f19f5aSHong Zhang       }
1339566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1349566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1353d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1369b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
137*da34b7cdSBarry Smith       PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
138*da34b7cdSBarry Smith       PetscCall(PCFactorSetUpMatSolverType(pc));
1393d1c1ea0SBarry Smith       if (!ilu->hdr.reuseordering) {
140f73b0415SBarry Smith         PetscBool canuseordering;
14103e5aca4SStefano Zampini 
1429566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
143f73b0415SBarry Smith         if (canuseordering) {
1449b54502bSHong Zhang           /* compute a new ordering for the ILU */
1459566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->row));
1469566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->col));
1479566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1489566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1499b54502bSHong Zhang           /*  Remove zeros along diagonal?     */
1501baa6e33SBarry Smith           if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
1519b54502bSHong Zhang         }
1522c7c0729SBarry Smith       }
1539566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1549566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1553d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1569b54502bSHong Zhang     }
1579566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
15800e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
15900e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1603ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1616baea169SHong Zhang     }
1626baea169SHong Zhang 
1639566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)ilu)->fact, pc->pmat, &((PC_Factor *)ilu)->info));
1649566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16500e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
16600e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1676baea169SHong Zhang     }
1689b54502bSHong Zhang   }
16900c67f3bSHong Zhang 
1709566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
17100c67f3bSHong Zhang   if (!stype) {
172ea799195SBarry Smith     MatSolverType solverpackage;
1739566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)ilu)->fact, &solverpackage));
1749566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
17500c67f3bSHong Zhang   }
1763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1779b54502bSHong Zhang }
1789b54502bSHong Zhang 
PCDestroy_ILU(PC pc)179d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ILU(PC pc)
180d71ae5a4SJacob Faibussowitsch {
1819b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1829b54502bSHong Zhang 
1839b54502bSHong Zhang   PetscFunctionBegin;
1849566063dSJacob Faibussowitsch   PetscCall(PCReset_ILU(pc));
1859566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->solvertype));
1869566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->ordering));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1882e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1909b54502bSHong Zhang }
1919b54502bSHong Zhang 
PCApply_ILU(PC pc,Vec x,Vec y)192d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ILU(PC pc, Vec x, Vec y)
193d71ae5a4SJacob Faibussowitsch {
1949b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1959b54502bSHong Zhang 
1969b54502bSHong Zhang   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)ilu)->fact, x, y));
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1999b54502bSHong Zhang }
2009b54502bSHong Zhang 
PCMatApply_ILU(PC pc,Mat X,Mat Y)201d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ILU(PC pc, Mat X, Mat Y)
202d71ae5a4SJacob Faibussowitsch {
2037b6e2003SPierre Jolivet   PC_ILU *ilu = (PC_ILU *)pc->data;
2047b6e2003SPierre Jolivet 
2057b6e2003SPierre Jolivet   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)ilu)->fact, X, Y));
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2087b6e2003SPierre Jolivet }
2097b6e2003SPierre Jolivet 
PCApplyTranspose_ILU(PC pc,Vec x,Vec y)210d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ILU(PC pc, Vec x, Vec y)
211d71ae5a4SJacob Faibussowitsch {
2129b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
2139b54502bSHong Zhang 
2149b54502bSHong Zhang   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(((PC_Factor *)ilu)->fact, x, y));
2163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2179b54502bSHong Zhang }
2189b54502bSHong Zhang 
PCMatApplyTranspose_ILU(PC pc,Mat X,Mat Y)2194dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_ILU(PC pc, Mat X, Mat Y)
2204dbf25a8SPierre Jolivet {
2214dbf25a8SPierre Jolivet   PC_ILU *ilu = (PC_ILU *)pc->data;
2224dbf25a8SPierre Jolivet 
2234dbf25a8SPierre Jolivet   PetscFunctionBegin;
2244dbf25a8SPierre Jolivet   PetscCall(MatMatSolveTranspose(((PC_Factor *)ilu)->fact, X, Y));
2254dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
2264dbf25a8SPierre Jolivet }
2274dbf25a8SPierre Jolivet 
PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)228d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc, Vec x, Vec y)
229d71ae5a4SJacob Faibussowitsch {
230f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
231f0b9ad6cSBarry Smith 
232f0b9ad6cSBarry Smith   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
2343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
235f0b9ad6cSBarry Smith }
236f0b9ad6cSBarry Smith 
PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)237d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ILU(PC pc, Vec x, Vec y)
238d71ae5a4SJacob Faibussowitsch {
239f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
240f0b9ad6cSBarry Smith 
241f0b9ad6cSBarry Smith   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244f0b9ad6cSBarry Smith }
245f0b9ad6cSBarry Smith 
2469b54502bSHong Zhang /*MC
2471d27aa22SBarry Smith      PCILU - Incomplete factorization preconditioners {cite}`dupont1968approximate`, {cite}`oliphant1961implicit`, {cite}`chan1997approximate`
2489b54502bSHong Zhang 
2499b54502bSHong Zhang    Options Database Keys:
2502401956bSBarry Smith +  -pc_factor_levels <k>                                 - number of levels of fill for ILU(k)
2512401956bSBarry Smith .  -pc_factor_in_place                                   - only for ILU(0) with natural ordering, reuses the space of the matrix for
2529b54502bSHong Zhang                                                          its factorization (overwrites original matrix)
2532401956bSBarry Smith .  -pc_factor_diagonal_fill                              - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
2542401956bSBarry Smith .  -pc_factor_reuse_ordering                             - reuse ordering of factorized matrix from previous factorization
25555ba2a51SBarry Smith .  -pc_factor_fill <nfill>                               - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2562401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal                    - reorder the matrix before factorization to remove zeros from the diagonal,
2579b54502bSHong Zhang                                                          this decreases the chance of getting a zero pivot
2582401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
2591d27aa22SBarry Smith -  -pc_factor_pivot_in_blocks                            - for block ILU(k) factorization, i.e. with `MATBAIJ` matrices with block size larger
2609b54502bSHong Zhang                                                          than 1 the diagonal blocks are factored with partial pivoting (this increases the
2619b54502bSHong Zhang                                                          stability of the ILU factorization
2629b54502bSHong Zhang 
2639b54502bSHong Zhang    Level: beginner
2649b54502bSHong Zhang 
26595452b02SPatrick Sanan    Notes:
266f1580f4eSBarry Smith    Only implemented for some matrix format and sequential. For parallel see `PCHYPRE` for hypre's ILU
2679b54502bSHong Zhang 
268f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ILU
2699b54502bSHong Zhang 
270f0b9ad6cSBarry Smith    The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
271f0b9ad6cSBarry Smith    even when the matrix is not symmetric since the U stores the diagonals of the factorization.
272f0b9ad6cSBarry Smith 
273f1580f4eSBarry Smith    If you are using `MATSEQAIJCUSPARSE` matrices (or `MATMPIAIJCUSPARSE` matrices with block Jacobi), factorization
274cd0a26f6SPaul Mullowney    is never done on the GPU).
275ac793be5SBarry Smith 
276562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCLU`, `PCICC`, `PCCHOLESKY`,
27754c05997SPierre Jolivet           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
278c2e3fba1SPatrick Sanan           `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
279db781477SPatrick Sanan           `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`,
280db781477SPatrick Sanan           `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()`
2819b54502bSHong Zhang M*/
2829b54502bSHong Zhang 
PCCreate_ILU(PC pc)283d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
284d71ae5a4SJacob Faibussowitsch {
2859b54502bSHong Zhang   PC_ILU *ilu;
2869b54502bSHong Zhang 
2879b54502bSHong Zhang   PetscFunctionBegin;
2884dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilu));
2893d1c1ea0SBarry Smith   pc->data = (void *)ilu;
2909566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ILU));
2919b54502bSHong Zhang 
29275567043SBarry Smith   ((PC_Factor *)ilu)->info.levels  = 0.;
293075768bcSBarry Smith   ((PC_Factor *)ilu)->info.fill    = 1.0;
2940a545947SLisandro Dalcin   ilu->col                         = NULL;
2950a545947SLisandro Dalcin   ilu->row                         = NULL;
296075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = PETSC_DEFAULT;
297075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = PETSC_DEFAULT;
298075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = PETSC_DEFAULT;
2999b54502bSHong Zhang 
300574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
3019b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3029b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
3037b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ILU;
3049b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
3054dbf25a8SPierre Jolivet   pc->ops->matapplytranspose   = PCMatApplyTranspose_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));
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3159b54502bSHong Zhang }
316