xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision c3b366b11f3ec31f57fee5376535cdadac340a0d)
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 
77087cfbeSBarry Smith PetscErrorCode  PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
89b54502bSHong Zhang {
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;
159b54502bSHong Zhang   PetscFunctionReturn(0);
169b54502bSHong Zhang }
179b54502bSHong Zhang 
18574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc)
199b54502bSHong Zhang {
209b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
219b54502bSHong Zhang   PetscErrorCode ierr;
229b54502bSHong Zhang 
239b54502bSHong Zhang   PetscFunctionBegin;
243d1c1ea0SBarry Smith   if (!ilu->hdr.inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);}
256bf464f9SBarry Smith   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);}
26fcfd50ebSBarry Smith   ierr = ISDestroy(&ilu->col);CHKERRQ(ierr);
279b54502bSHong Zhang   PetscFunctionReturn(0);
289b54502bSHong Zhang }
299b54502bSHong Zhang 
307087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
319b54502bSHong Zhang {
32075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU*)pc->data;
339b54502bSHong Zhang 
349b54502bSHong Zhang   PetscFunctionBegin;
354c9036c7SBarry Smith   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
36ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
379b54502bSHong Zhang   }
38075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt      = dt;
39075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
40075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount = dtcount;
414c9036c7SBarry Smith   ((PC_Factor*)ilu)->info.usedt   = 1.0;
429b54502bSHong Zhang   PetscFunctionReturn(0);
439b54502bSHong Zhang }
449b54502bSHong Zhang 
454416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc)
469b54502bSHong Zhang {
479b54502bSHong Zhang   PetscErrorCode ierr;
4878fc6b22SHong Zhang   PetscInt       itmp;
498afaa268SBarry Smith   PetscBool      flg,set;
509b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
519b54502bSHong Zhang   PetscReal      tol;
529b54502bSHong Zhang 
539b54502bSHong Zhang   PetscFunctionBegin;
54e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ILU Options");CHKERRQ(ierr);
55e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
568ff23777SHong Zhang 
57075768bcSBarry Smith   ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
58075768bcSBarry Smith   if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
592fa5cd67SKarl Rupp 
608afaa268SBarry Smith   ierr = PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
618afaa268SBarry Smith   if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg;
622401956bSBarry Smith   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
639b54502bSHong Zhang   if (flg) {
649b54502bSHong Zhang     tol  = PETSC_DECIDE;
652401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
662401956bSBarry Smith     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
679b54502bSHong Zhang   }
689b54502bSHong Zhang 
699b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
709b54502bSHong Zhang   PetscFunctionReturn(0);
719b54502bSHong Zhang }
729b54502bSHong Zhang 
739b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc)
749b54502bSHong Zhang {
759b54502bSHong Zhang   PetscErrorCode         ierr;
769b54502bSHong Zhang   PC_ILU                 *ilu = (PC_ILU*)pc->data;
77f3a39becSBarry Smith   MatInfo                info;
78ace3abfcSBarry Smith   PetscBool              flg;
79ea799195SBarry Smith   MatSolverType          stype;
8000e125f8SBarry Smith   MatFactorError         err;
819b54502bSHong Zhang 
829b54502bSHong Zhang   PetscFunctionBegin;
83c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
8492927226SBarry Smith   /* ugly hack to change default, since it is not support by some matrix types */
8592927226SBarry Smith   if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
8622d28d08SBarry Smith     ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr);
8792927226SBarry Smith     if (!flg) {
8822d28d08SBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr);
8992927226SBarry Smith       if (!flg) {
9092927226SBarry Smith         ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
91994fe344SLisandro Dalcin         ierr = PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");CHKERRQ(ierr);
9292927226SBarry Smith       }
9392927226SBarry Smith     }
9492927226SBarry Smith   }
9592927226SBarry Smith 
9684d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
973d1c1ea0SBarry Smith   if (ilu->hdr.inplace) {
989b54502bSHong Zhang     if (!pc->setupcalled) {
999b54502bSHong Zhang 
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 */
102075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
1033bb1ff40SBarry Smith       if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
1043bb1ff40SBarry Smith       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
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 
112*c3b366b1Sprj-     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
11300e125f8SBarry Smith     ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
11400e125f8SBarry Smith     if (err) { /* Factor() fails */
11500e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1166baea169SHong Zhang       PetscFunctionReturn(0);
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 */
121b33856dcSBarry Smith     ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr);
1229b54502bSHong Zhang   } else {
1239b54502bSHong Zhang     if (!pc->setupcalled) {
1249b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
125075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
1263bb1ff40SBarry Smith       if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
1273bb1ff40SBarry Smith       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
1289b54502bSHong Zhang       /*  Remove zeros along diagonal?     */
1299b54502bSHong Zhang       if (ilu->nonzerosalongdiagonal) {
1309b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
1319b54502bSHong Zhang       }
13286daa42cSBarry Smith       if (!((PC_Factor*)ilu)->fact) {
13311bfe483SHong Zhang         ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
134a1f19f5aSHong Zhang       }
135075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
136075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1373d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1382fa5cd67SKarl Rupp 
1393bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
1409b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
1413d1c1ea0SBarry Smith       if (!ilu->hdr.reuseordering) {
1429b54502bSHong Zhang         /* compute a new ordering for the ILU */
143fcfd50ebSBarry Smith         ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);
144fcfd50ebSBarry Smith         ierr = ISDestroy(&ilu->col);CHKERRQ(ierr);
145075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
1463bb1ff40SBarry Smith         if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
1473bb1ff40SBarry Smith         if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
1489b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
1499b54502bSHong Zhang         if (ilu->nonzerosalongdiagonal) {
1509b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
1519b54502bSHong Zhang         }
1529b54502bSHong Zhang       }
1536bf464f9SBarry Smith       ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
15414798fb4SJed Brown       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
155075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
156075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
1573d1c1ea0SBarry Smith       ilu->hdr.actualfill = info.fill_ratio_needed;
1582fa5cd67SKarl Rupp 
1593bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
1609b54502bSHong Zhang     }
16100e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr);
16200e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
16300e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1646baea169SHong Zhang       PetscFunctionReturn(0);
1656baea169SHong Zhang     }
1666baea169SHong Zhang 
167075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
16800e125f8SBarry Smith     ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr);
16900e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
17000e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1716baea169SHong Zhang     }
1729b54502bSHong Zhang   }
17300c67f3bSHong Zhang 
1743ca39a21SBarry Smith   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
17500c67f3bSHong Zhang   if (!stype) {
176ea799195SBarry Smith     MatSolverType solverpackage;
1773ca39a21SBarry Smith     ierr = MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage);CHKERRQ(ierr);
1783ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
17900c67f3bSHong Zhang   }
1809b54502bSHong Zhang   PetscFunctionReturn(0);
1819b54502bSHong Zhang }
1829b54502bSHong Zhang 
1839b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc)
1849b54502bSHong Zhang {
1859b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1869b54502bSHong Zhang   PetscErrorCode ierr;
1879b54502bSHong Zhang 
1889b54502bSHong Zhang   PetscFunctionBegin;
189574deadeSBarry Smith   ierr = PCReset_ILU(pc);CHKERRQ(ierr);
190503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
191503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
192c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1939b54502bSHong Zhang   PetscFunctionReturn(0);
1949b54502bSHong Zhang }
1959b54502bSHong Zhang 
1969b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
1979b54502bSHong Zhang {
1989b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1999b54502bSHong Zhang   PetscErrorCode ierr;
2009b54502bSHong Zhang 
2019b54502bSHong Zhang   PetscFunctionBegin;
202075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2039b54502bSHong Zhang   PetscFunctionReturn(0);
2049b54502bSHong Zhang }
2059b54502bSHong Zhang 
2069b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
2079b54502bSHong Zhang {
2089b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2099b54502bSHong Zhang   PetscErrorCode ierr;
2109b54502bSHong Zhang 
2119b54502bSHong Zhang   PetscFunctionBegin;
212075768bcSBarry Smith   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2139b54502bSHong Zhang   PetscFunctionReturn(0);
2149b54502bSHong Zhang }
2159b54502bSHong Zhang 
216f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
217f0b9ad6cSBarry Smith {
218f0b9ad6cSBarry Smith   PetscErrorCode ierr;
219f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
220f0b9ad6cSBarry Smith 
221f0b9ad6cSBarry Smith   PetscFunctionBegin;
222f0b9ad6cSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
223f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
224f0b9ad6cSBarry Smith }
225f0b9ad6cSBarry Smith 
226f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
227f0b9ad6cSBarry Smith {
228f0b9ad6cSBarry Smith   PetscErrorCode ierr;
229f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
230f0b9ad6cSBarry Smith 
231f0b9ad6cSBarry Smith   PetscFunctionBegin;
232f0b9ad6cSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
233f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
234f0b9ad6cSBarry Smith }
235f0b9ad6cSBarry Smith 
2369b54502bSHong Zhang /*MC
2379b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
2389b54502bSHong Zhang 
2399b54502bSHong Zhang    Options Database Keys:
2402401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
2412401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
2429b54502bSHong Zhang                       its factorization (overwrites original matrix)
2432401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
2442401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
24555ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2462401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
2479b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
2482401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
249967c93d3SBarry Smith -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
2509b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
2519b54502bSHong Zhang                              stability of the ILU factorization
2529b54502bSHong Zhang 
2539b54502bSHong Zhang    Level: beginner
2549b54502bSHong Zhang 
25595452b02SPatrick Sanan    Notes:
25695452b02SPatrick Sanan     Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
2579b54502bSHong Zhang 
2589b54502bSHong Zhang           For BAIJ matrices this implements a point block ILU
2599b54502bSHong Zhang 
260f0b9ad6cSBarry Smith           The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
261f0b9ad6cSBarry Smith           even when the matrix is not symmetric since the U stores the diagonals of the factorization.
262f0b9ad6cSBarry Smith 
263cd0a26f6SPaul Mullowney           If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization
264cd0a26f6SPaul Mullowney           is never done on the GPU).
265ac793be5SBarry Smith 
266c582cd25SBarry Smith    References:
26796a0c994SBarry Smith +  1. - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
26896a0c994SBarry Smith    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
26996a0c994SBarry Smith .  2. -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
27096a0c994SBarry Smith -  3. -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
27196a0c994SBarry Smith       Chapter in Parallel Numerical
272c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
27396a0c994SBarry Smith       Science and Engineering, Kluwer.
274c582cd25SBarry Smith 
275c582cd25SBarry Smith 
2769b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
277145b9266SHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
278b7c853c4SBarry Smith            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
27992e9c092SBarry Smith            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
28092e9c092SBarry Smith            PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace()
2819b54502bSHong Zhang 
2829b54502bSHong Zhang M*/
2839b54502bSHong Zhang 
2848cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
2859b54502bSHong Zhang {
2869b54502bSHong Zhang   PetscErrorCode ierr;
2879b54502bSHong Zhang   PC_ILU         *ilu;
2889b54502bSHong Zhang 
2899b54502bSHong Zhang   PetscFunctionBegin;
290b00a9115SJed Brown   ierr     = PetscNewLog(pc,&ilu);CHKERRQ(ierr);
2913d1c1ea0SBarry Smith   pc->data = (void*)ilu;
2923d1c1ea0SBarry Smith   ierr     = PCFactorInitialize(pc);CHKERRQ(ierr);
2939b54502bSHong Zhang 
294879e8a4dSBarry Smith   ((PC_Factor*)ilu)->factortype         = MAT_FACTOR_ILU;
29575567043SBarry Smith   ((PC_Factor*)ilu)->info.levels        = 0.;
296075768bcSBarry Smith   ((PC_Factor*)ilu)->info.fill          = 1.0;
2979b54502bSHong Zhang   ilu->col                              = 0;
2989b54502bSHong Zhang   ilu->row                              = 0;
29919fd82e9SBarry Smith   ierr                                  = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
300075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt            = PETSC_DEFAULT;
301075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount       = PETSC_DEFAULT;
302075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol         = PETSC_DEFAULT;
3039b54502bSHong Zhang 
304574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
3059b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3069b54502bSHong Zhang   pc->ops->apply               = PCApply_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;
3139b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
314bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
315bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
3169b54502bSHong Zhang   PetscFunctionReturn(0);
3179b54502bSHong Zhang }
318