xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 03e5aca47aa9bd9d5a48c628e4311d43ab2119ff)
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;
1313bcc0bdSJacob Faibussowitsch   if (z == (PetscReal)PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
142fa5cd67SKarl Rupp   else ilu->nonzerosalongdiagonaltol = z;
153ba16761SJacob 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));
263ba16761SJacob 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;
413ba16761SJacob 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();
683ba16761SJacob 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;
1173ba16761SJacob 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;
127*03e5aca4SStefano Zampini 
128*03e5aca4SStefano Zampini       PetscCall(PCFactorSetUpMatSolverType(pc));
1299566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
130f73b0415SBarry Smith       if (canuseordering) {
1319566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1329566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1339b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
1341baa6e33SBarry Smith         if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
135a1f19f5aSHong Zhang       }
1369566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1379566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1383d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1399b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1403d1c1ea0SBarry Smith       if (!ilu->hdr.reuseordering) {
141f73b0415SBarry Smith         PetscBool canuseordering;
142*03e5aca4SStefano Zampini 
1439566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)ilu)->fact));
144*03e5aca4SStefano Zampini         PetscCall(PCFactorSetUpMatSolverType(pc));
1459566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)ilu)->fact, &canuseordering));
146f73b0415SBarry Smith         if (canuseordering) {
1479b54502bSHong Zhang           /* compute a new ordering for the ILU */
1489566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->row));
1499566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&ilu->col));
1509566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1519566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)ilu)->ordering, &ilu->row, &ilu->col));
1529b54502bSHong Zhang           /*  Remove zeros along diagonal?     */
1531baa6e33SBarry Smith           if (ilu->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, ilu->nonzerosalongdiagonaltol, ilu->row, ilu->col));
1549b54502bSHong Zhang         }
1552c7c0729SBarry Smith       }
1569566063dSJacob Faibussowitsch       PetscCall(MatILUFactorSymbolic(((PC_Factor *)ilu)->fact, pc->pmat, ilu->row, ilu->col, &((PC_Factor *)ilu)->info));
1579566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)ilu)->fact, MAT_LOCAL, &info));
1583d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1599b54502bSHong Zhang     }
1609566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16100e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
16200e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1633ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1646baea169SHong Zhang     }
1656baea169SHong Zhang 
1669566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)ilu)->fact, pc->pmat, &((PC_Factor *)ilu)->info));
1679566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)ilu)->fact, &err));
16800e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
16900e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1706baea169SHong Zhang     }
1719b54502bSHong Zhang   }
17200c67f3bSHong Zhang 
1739566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
17400c67f3bSHong Zhang   if (!stype) {
175ea799195SBarry Smith     MatSolverType solverpackage;
1769566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)ilu)->fact, &solverpackage));
1779566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
17800c67f3bSHong Zhang   }
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1809b54502bSHong Zhang }
1819b54502bSHong Zhang 
182d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ILU(PC pc)
183d71ae5a4SJacob Faibussowitsch {
1849b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1859b54502bSHong Zhang 
1869b54502bSHong Zhang   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(PCReset_ILU(pc));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->solvertype));
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)ilu)->ordering));
1909566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1912e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1939b54502bSHong Zhang }
1949b54502bSHong Zhang 
195d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ILU(PC pc, Vec x, Vec y)
196d71ae5a4SJacob Faibussowitsch {
1979b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
1989b54502bSHong Zhang 
1999b54502bSHong Zhang   PetscFunctionBegin;
2009566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)ilu)->fact, x, y));
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2029b54502bSHong Zhang }
2039b54502bSHong Zhang 
204d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ILU(PC pc, Mat X, Mat Y)
205d71ae5a4SJacob Faibussowitsch {
2067b6e2003SPierre Jolivet   PC_ILU *ilu = (PC_ILU *)pc->data;
2077b6e2003SPierre Jolivet 
2087b6e2003SPierre Jolivet   PetscFunctionBegin;
2099566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)ilu)->fact, X, Y));
2103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2117b6e2003SPierre Jolivet }
2127b6e2003SPierre Jolivet 
213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ILU(PC pc, Vec x, Vec y)
214d71ae5a4SJacob Faibussowitsch {
2159b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU *)pc->data;
2169b54502bSHong Zhang 
2179b54502bSHong Zhang   PetscFunctionBegin;
2189566063dSJacob Faibussowitsch   PetscCall(MatSolveTranspose(((PC_Factor *)ilu)->fact, x, y));
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2209b54502bSHong Zhang }
2219b54502bSHong Zhang 
222d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc, Vec x, Vec y)
223d71ae5a4SJacob Faibussowitsch {
224f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
225f0b9ad6cSBarry Smith 
226f0b9ad6cSBarry Smith   PetscFunctionBegin;
2279566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229f0b9ad6cSBarry Smith }
230f0b9ad6cSBarry Smith 
231d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ILU(PC pc, Vec x, Vec y)
232d71ae5a4SJacob Faibussowitsch {
233f0b9ad6cSBarry Smith   PC_ILU *icc = (PC_ILU *)pc->data;
234f0b9ad6cSBarry Smith 
235f0b9ad6cSBarry Smith   PetscFunctionBegin;
2369566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
238f0b9ad6cSBarry Smith }
239f0b9ad6cSBarry Smith 
2409b54502bSHong Zhang /*MC
2419b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
2429b54502bSHong Zhang 
2439b54502bSHong Zhang    Options Database Keys:
2442401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
2452401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
2469b54502bSHong Zhang                       its factorization (overwrites original matrix)
2472401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
2482401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
24955ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2502401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
2519b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
2522401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
253967c93d3SBarry Smith -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
2549b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
2559b54502bSHong Zhang                              stability of the ILU factorization
2569b54502bSHong Zhang 
2579b54502bSHong Zhang    Level: beginner
2589b54502bSHong Zhang 
25995452b02SPatrick Sanan    Notes:
260f1580f4eSBarry Smith    Only implemented for some matrix format and sequential. For parallel see `PCHYPRE` for hypre's ILU
2619b54502bSHong Zhang 
262f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ILU
2639b54502bSHong Zhang 
264f0b9ad6cSBarry Smith    The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
265f0b9ad6cSBarry Smith    even when the matrix is not symmetric since the U stores the diagonals of the factorization.
266f0b9ad6cSBarry Smith 
267f1580f4eSBarry Smith    If you are using `MATSEQAIJCUSPARSE` matrices (or `MATMPIAIJCUSPARSE` matrices with block Jacobi), factorization
268cd0a26f6SPaul Mullowney    is never done on the GPU).
269ac793be5SBarry Smith 
270c582cd25SBarry Smith    References:
271606c0280SSatish Balay +  * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
27296a0c994SBarry Smith    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
273606c0280SSatish Balay .  * -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
274606c0280SSatish Balay -  * -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
27596a0c994SBarry Smith       Chapter in Parallel Numerical
276c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
27796a0c994SBarry Smith       Science and Engineering, Kluwer.
278c582cd25SBarry Smith 
279f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCLU`, `PCICC`, `PCCHOLESKY`,
280db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftSetType()`, `PCFactorSetAmount()`,
281c2e3fba1SPatrick Sanan           `PCFactorSetDropTolerance()`, `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
282db781477SPatrick Sanan           `PCFactorSetLevels()`, `PCFactorSetUseInPlace()`, `PCFactorSetAllowDiagonalFill()`, `PCFactorSetPivotInBlocks()`,
283db781477SPatrick Sanan           `PCFactorGetAllowDiagonalFill()`, `PCFactorGetUseInPlace()`
2849b54502bSHong Zhang M*/
2859b54502bSHong Zhang 
286d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
287d71ae5a4SJacob Faibussowitsch {
2889b54502bSHong Zhang   PC_ILU *ilu;
2899b54502bSHong Zhang 
2909b54502bSHong Zhang   PetscFunctionBegin;
2914dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilu));
2923d1c1ea0SBarry Smith   pc->data = (void *)ilu;
2939566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ILU));
2949b54502bSHong Zhang 
29575567043SBarry Smith   ((PC_Factor *)ilu)->info.levels  = 0.;
296075768bcSBarry Smith   ((PC_Factor *)ilu)->info.fill    = 1.0;
2970a545947SLisandro Dalcin   ilu->col                         = NULL;
2980a545947SLisandro Dalcin   ilu->row                         = NULL;
299075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dt      = PETSC_DEFAULT;
300075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcount = PETSC_DEFAULT;
301075768bcSBarry Smith   ((PC_Factor *)ilu)->info.dtcol   = PETSC_DEFAULT;
3029b54502bSHong Zhang 
303574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
3049b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3059b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
3067b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ILU;
3079b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
3089b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
3099b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
31092e08861SBarry Smith   pc->ops->view                = PCView_Factor;
311f0b9ad6cSBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
312f0b9ad6cSBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
3130a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", PCFactorSetDropTolerance_ILU));
3159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_ILU));
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3179b54502bSHong Zhang }
318