xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 00c67f3b09c0bcda06af5ed306d845d9138e5003)
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 
79b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/
8afaefe49SHong Zhang #undef __FUNCT__
92401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_ILU"
107087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseFill_ILU(PC pc,PetscBool flag)
112401956bSBarry Smith {
12075768bcSBarry Smith   PC_ILU *lu = (PC_ILU*)pc->data;
132401956bSBarry Smith 
142401956bSBarry Smith   PetscFunctionBegin;
152401956bSBarry Smith   lu->reusefill = flag;
162401956bSBarry Smith   PetscFunctionReturn(0);
172401956bSBarry Smith }
182401956bSBarry Smith 
192401956bSBarry Smith #undef __FUNCT__
202401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU"
217087cfbeSBarry Smith PetscErrorCode  PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
229b54502bSHong Zhang {
239b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU*)pc->data;
249b54502bSHong Zhang 
259b54502bSHong Zhang   PetscFunctionBegin;
269b54502bSHong Zhang   ilu->nonzerosalongdiagonal = PETSC_TRUE;
272fa5cd67SKarl Rupp   if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10;
282fa5cd67SKarl Rupp   else ilu->nonzerosalongdiagonaltol = z;
299b54502bSHong Zhang   PetscFunctionReturn(0);
309b54502bSHong Zhang }
319b54502bSHong Zhang 
329b54502bSHong Zhang #undef __FUNCT__
33574deadeSBarry Smith #define __FUNCT__ "PCReset_ILU"
34574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc)
359b54502bSHong Zhang {
369b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
379b54502bSHong Zhang   PetscErrorCode ierr;
389b54502bSHong Zhang 
399b54502bSHong Zhang   PetscFunctionBegin;
40298cc208SBarry Smith   if (!ilu->inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);}
416bf464f9SBarry Smith   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);}
42fcfd50ebSBarry Smith   ierr = ISDestroy(&ilu->col);CHKERRQ(ierr);
439b54502bSHong Zhang   PetscFunctionReturn(0);
449b54502bSHong Zhang }
459b54502bSHong Zhang 
469b54502bSHong Zhang #undef __FUNCT__
47b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_ILU"
487087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
499b54502bSHong Zhang {
50075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU*)pc->data;
519b54502bSHong Zhang 
529b54502bSHong Zhang   PetscFunctionBegin;
534c9036c7SBarry Smith   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
54ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
559b54502bSHong Zhang   }
56075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt      = dt;
57075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
58075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount = dtcount;
594c9036c7SBarry Smith   ((PC_Factor*)ilu)->info.usedt   = 1.0;
609b54502bSHong Zhang   PetscFunctionReturn(0);
619b54502bSHong Zhang }
629b54502bSHong Zhang 
639b54502bSHong Zhang #undef __FUNCT__
642401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
657087cfbeSBarry Smith PetscErrorCode  PCFactorSetReuseOrdering_ILU(PC pc,PetscBool flag)
669b54502bSHong Zhang {
67075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU*)pc->data;
689b54502bSHong Zhang 
699b54502bSHong Zhang   PetscFunctionBegin;
709b54502bSHong Zhang   ilu->reuseordering = flag;
719b54502bSHong Zhang   PetscFunctionReturn(0);
729b54502bSHong Zhang }
739b54502bSHong Zhang 
749b54502bSHong Zhang #undef __FUNCT__
752401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
768e37d05fSBarry Smith PetscErrorCode  PCFactorSetUseInPlace_ILU(PC pc,PetscBool flg)
779b54502bSHong Zhang {
78075768bcSBarry Smith   PC_ILU *dir = (PC_ILU*)pc->data;
799b54502bSHong Zhang 
809b54502bSHong Zhang   PetscFunctionBegin;
818e37d05fSBarry Smith   dir->inplace = flg;
828e37d05fSBarry Smith   PetscFunctionReturn(0);
838e37d05fSBarry Smith }
848e37d05fSBarry Smith 
858e37d05fSBarry Smith #undef __FUNCT__
868e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_ILU"
878e37d05fSBarry Smith PetscErrorCode  PCFactorGetUseInPlace_ILU(PC pc,PetscBool *flg)
888e37d05fSBarry Smith {
898e37d05fSBarry Smith   PC_ILU *dir = (PC_ILU*)pc->data;
908e37d05fSBarry Smith 
918e37d05fSBarry Smith   PetscFunctionBegin;
928e37d05fSBarry Smith   *flg = dir->inplace;
939b54502bSHong Zhang   PetscFunctionReturn(0);
949b54502bSHong Zhang }
959b54502bSHong Zhang 
969b54502bSHong Zhang #undef __FUNCT__
979b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU"
984416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc)
999b54502bSHong Zhang {
1009b54502bSHong Zhang   PetscErrorCode ierr;
10178fc6b22SHong Zhang   PetscInt       itmp;
1028afaa268SBarry Smith   PetscBool      flg,set;
1039b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1049b54502bSHong Zhang   PetscReal      tol;
1059b54502bSHong Zhang 
1069b54502bSHong Zhang   PetscFunctionBegin;
107e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ILU Options");CHKERRQ(ierr);
108e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1098ff23777SHong Zhang 
110075768bcSBarry Smith   ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
111075768bcSBarry Smith   if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
1122fa5cd67SKarl Rupp 
1138afaa268SBarry 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);
1148afaa268SBarry Smith   if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg;
1152401956bSBarry Smith   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1169b54502bSHong Zhang   if (flg) {
1179b54502bSHong Zhang     tol  = PETSC_DECIDE;
1182401956bSBarry Smith     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1192401956bSBarry Smith     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1209b54502bSHong Zhang   }
1219b54502bSHong Zhang 
1229b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1239b54502bSHong Zhang   PetscFunctionReturn(0);
1249b54502bSHong Zhang }
1259b54502bSHong Zhang 
1269b54502bSHong Zhang #undef __FUNCT__
1279b54502bSHong Zhang #define __FUNCT__ "PCView_ILU"
1289b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
1299b54502bSHong Zhang {
1309b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1319b54502bSHong Zhang   PetscErrorCode ierr;
132ace3abfcSBarry Smith   PetscBool      iascii;
1339b54502bSHong Zhang 
1349b54502bSHong Zhang   PetscFunctionBegin;
135251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1369b54502bSHong Zhang   if (iascii) {
137914a5d51SHong Zhang     if (ilu->inplace) {
138914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: in-place factorization\n");CHKERRQ(ierr);
1399b54502bSHong Zhang     } else {
140914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: out-of-place factorization\n");CHKERRQ(ierr);
141145b9266SHong Zhang     }
142145b9266SHong Zhang 
143914a5d51SHong Zhang     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing fill from past factorization\n");CHKERRQ(ierr);}
144914a5d51SHong Zhang     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1459b54502bSHong Zhang   }
146914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1479b54502bSHong Zhang   PetscFunctionReturn(0);
1489b54502bSHong Zhang }
1499b54502bSHong Zhang 
1509b54502bSHong Zhang #undef __FUNCT__
1519b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU"
1529b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc)
1539b54502bSHong Zhang {
1549b54502bSHong Zhang   PetscErrorCode ierr;
1559b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
156f3a39becSBarry Smith   MatInfo        info;
157ace3abfcSBarry Smith   PetscBool      flg;
158*00c67f3bSHong Zhang   const MatSolverPackage stype;
1599b54502bSHong Zhang 
1609b54502bSHong Zhang   PetscFunctionBegin;
16192927226SBarry Smith   /* ugly hack to change default, since it is not support by some matrix types */
16292927226SBarry Smith   if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
16322d28d08SBarry Smith     ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr);
16492927226SBarry Smith     if (!flg) {
16522d28d08SBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr);
16692927226SBarry Smith       if (!flg) {
16792927226SBarry Smith         ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
168955c1f14SBarry Smith         PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");CHKERRQ(ierr);
16992927226SBarry Smith       }
17092927226SBarry Smith     }
17192927226SBarry Smith   }
17292927226SBarry Smith 
17384d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
1749b54502bSHong Zhang   if (ilu->inplace) {
1759b54502bSHong Zhang     if (!pc->setupcalled) {
1769b54502bSHong Zhang 
1779b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
1789b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
179075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
1803bb1ff40SBarry Smith       if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
1813bb1ff40SBarry Smith       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
1829b54502bSHong Zhang     }
1839b54502bSHong Zhang 
1849b54502bSHong Zhang     /* In place ILU only makes sense with fill factor of 1.0 because
1859b54502bSHong Zhang        cannot have levels of fill */
186075768bcSBarry Smith     ((PC_Factor*)ilu)->info.fill          = 1.0;
18775567043SBarry Smith     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
1882fa5cd67SKarl Rupp 
189365a8a9eSBarry Smith     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKERRQ(ierr);
1904cccfbddSHong Zhang     if (pc->pmat->errortype) { /* Factor() fails */
1914cccfbddSHong Zhang       pc->failedreason = (PCFailedReason)pc->pmat->errortype;
1926baea169SHong Zhang       PetscFunctionReturn(0);
1936baea169SHong Zhang     }
1946baea169SHong Zhang 
195075768bcSBarry Smith     ((PC_Factor*)ilu)->fact = pc->pmat;
196b33856dcSBarry Smith     /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */
197b33856dcSBarry Smith     ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr);
1989b54502bSHong Zhang   } else {
1994cccfbddSHong Zhang     Mat F;
2009b54502bSHong Zhang     if (!pc->setupcalled) {
2019b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
202075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
2033bb1ff40SBarry Smith       if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
2043bb1ff40SBarry Smith       if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
2059b54502bSHong Zhang       /*  Remove zeros along diagonal?     */
2069b54502bSHong Zhang       if (ilu->nonzerosalongdiagonal) {
2079b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
2089b54502bSHong Zhang       }
20986daa42cSBarry Smith       if (!((PC_Factor*)ilu)->fact) {
21011bfe483SHong Zhang         ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
211a1f19f5aSHong Zhang       }
212075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
213075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2142fa5cd67SKarl Rupp 
215f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
2162fa5cd67SKarl Rupp 
2173bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2189b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2199b54502bSHong Zhang       if (!ilu->reuseordering) {
2209b54502bSHong Zhang         /* compute a new ordering for the ILU */
221fcfd50ebSBarry Smith         ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);
222fcfd50ebSBarry Smith         ierr = ISDestroy(&ilu->col);CHKERRQ(ierr);
223075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
2243bb1ff40SBarry Smith         if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);}
2253bb1ff40SBarry Smith         if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);}
2269b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
2279b54502bSHong Zhang         if (ilu->nonzerosalongdiagonal) {
2289b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
2299b54502bSHong Zhang         }
2309b54502bSHong Zhang       }
2316bf464f9SBarry Smith       ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
23214798fb4SJed Brown       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
233075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
234075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
2352fa5cd67SKarl Rupp 
236f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
2372fa5cd67SKarl Rupp 
2383bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2399b54502bSHong Zhang     }
2404cccfbddSHong Zhang     F = ((PC_Factor*)ilu)->fact;
2414cccfbddSHong Zhang     if (F->errortype) { /* FactorSymbolic() fails */
2424cccfbddSHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
2436baea169SHong Zhang       PetscFunctionReturn(0);
2446baea169SHong Zhang     }
2456baea169SHong Zhang 
246075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
2474cccfbddSHong Zhang     if (F->errortype) { /* FactorNumeric() fails */
2484cccfbddSHong Zhang       pc->failedreason = (PCFailedReason)F->errortype;
2496baea169SHong Zhang     }
2509b54502bSHong Zhang   }
251*00c67f3bSHong Zhang 
252*00c67f3bSHong Zhang   ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
253*00c67f3bSHong Zhang   if (!stype) {
254*00c67f3bSHong Zhang     ierr = PCFactorSetMatSolverPackage(pc,((PC_Factor*)ilu)->fact->solvertype);CHKERRQ(ierr);
255*00c67f3bSHong Zhang   }
2569b54502bSHong Zhang   PetscFunctionReturn(0);
2579b54502bSHong Zhang }
2589b54502bSHong Zhang 
2599b54502bSHong Zhang #undef __FUNCT__
2609b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU"
2619b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc)
2629b54502bSHong Zhang {
2639b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2649b54502bSHong Zhang   PetscErrorCode ierr;
2659b54502bSHong Zhang 
2669b54502bSHong Zhang   PetscFunctionBegin;
267574deadeSBarry Smith   ierr = PCReset_ILU(pc);CHKERRQ(ierr);
268503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
269503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
270c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2719b54502bSHong Zhang   PetscFunctionReturn(0);
2729b54502bSHong Zhang }
2739b54502bSHong Zhang 
2749b54502bSHong Zhang #undef __FUNCT__
2759b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU"
2769b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
2779b54502bSHong Zhang {
2789b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2799b54502bSHong Zhang   PetscErrorCode ierr;
2809b54502bSHong Zhang 
2819b54502bSHong Zhang   PetscFunctionBegin;
282075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2839b54502bSHong Zhang   PetscFunctionReturn(0);
2849b54502bSHong Zhang }
2859b54502bSHong Zhang 
2869b54502bSHong Zhang #undef __FUNCT__
2879b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU"
2889b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
2899b54502bSHong Zhang {
2909b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2919b54502bSHong Zhang   PetscErrorCode ierr;
2929b54502bSHong Zhang 
2939b54502bSHong Zhang   PetscFunctionBegin;
294075768bcSBarry Smith   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2959b54502bSHong Zhang   PetscFunctionReturn(0);
2969b54502bSHong Zhang }
2979b54502bSHong Zhang 
298f0b9ad6cSBarry Smith #undef __FUNCT__
299f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricLeft_ILU"
300f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
301f0b9ad6cSBarry Smith {
302f0b9ad6cSBarry Smith   PetscErrorCode ierr;
303f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
304f0b9ad6cSBarry Smith 
305f0b9ad6cSBarry Smith   PetscFunctionBegin;
306f0b9ad6cSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
307f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
308f0b9ad6cSBarry Smith }
309f0b9ad6cSBarry Smith 
310f0b9ad6cSBarry Smith #undef __FUNCT__
311f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricRight_ILU"
312f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
313f0b9ad6cSBarry Smith {
314f0b9ad6cSBarry Smith   PetscErrorCode ierr;
315f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
316f0b9ad6cSBarry Smith 
317f0b9ad6cSBarry Smith   PetscFunctionBegin;
318f0b9ad6cSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
319f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
320f0b9ad6cSBarry Smith }
321f0b9ad6cSBarry Smith 
3229b54502bSHong Zhang /*MC
3239b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
3249b54502bSHong Zhang 
3259b54502bSHong Zhang    Options Database Keys:
3262401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
3272401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
3289b54502bSHong Zhang                       its factorization (overwrites original matrix)
3292401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
3302401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
33155ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
3322401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
3339b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
3342401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
335967c93d3SBarry Smith -  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
3369b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
3379b54502bSHong Zhang                              stability of the ILU factorization
3389b54502bSHong Zhang 
3399b54502bSHong Zhang    Level: beginner
3409b54502bSHong Zhang 
3419b54502bSHong Zhang   Concepts: incomplete factorization
3429b54502bSHong Zhang 
343d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
3449b54502bSHong Zhang 
3459b54502bSHong Zhang           For BAIJ matrices this implements a point block ILU
3469b54502bSHong Zhang 
347f0b9ad6cSBarry Smith           The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
348f0b9ad6cSBarry Smith           even when the matrix is not symmetric since the U stores the diagonals of the factorization.
349f0b9ad6cSBarry Smith 
350cd0a26f6SPaul Mullowney           If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization
351cd0a26f6SPaul Mullowney           is never done on the GPU).
352ac793be5SBarry Smith 
353c582cd25SBarry Smith    References:
35496a0c994SBarry Smith +  1. - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
35596a0c994SBarry Smith    self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968.
35696a0c994SBarry Smith .  2. -  T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961.
35796a0c994SBarry Smith -  3. -  TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
35896a0c994SBarry Smith       Chapter in Parallel Numerical
359c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
36096a0c994SBarry Smith       Science and Engineering, Kluwer.
361c582cd25SBarry Smith 
362c582cd25SBarry Smith 
3639b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
364145b9266SHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
365b7c853c4SBarry Smith            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
36692e9c092SBarry Smith            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
36792e9c092SBarry Smith            PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace()
3689b54502bSHong Zhang 
3699b54502bSHong Zhang M*/
3709b54502bSHong Zhang 
3719b54502bSHong Zhang #undef __FUNCT__
3729b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU"
3738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc)
3749b54502bSHong Zhang {
3759b54502bSHong Zhang   PetscErrorCode ierr;
3769b54502bSHong Zhang   PC_ILU         *ilu;
3779b54502bSHong Zhang 
3789b54502bSHong Zhang   PetscFunctionBegin;
379b00a9115SJed Brown   ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr);
3809b54502bSHong Zhang 
381075768bcSBarry Smith   ((PC_Factor*)ilu)->fact               = 0;
382075768bcSBarry Smith   ierr                                  = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
383879e8a4dSBarry Smith   ((PC_Factor*)ilu)->factortype         = MAT_FACTOR_ILU;
38475567043SBarry Smith   ((PC_Factor*)ilu)->info.levels        = 0.;
385075768bcSBarry Smith   ((PC_Factor*)ilu)->info.fill          = 1.0;
3869b54502bSHong Zhang   ilu->col                              = 0;
3879b54502bSHong Zhang   ilu->row                              = 0;
3889b54502bSHong Zhang   ilu->inplace                          = PETSC_FALSE;
38919fd82e9SBarry Smith   ierr                                  = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
3909b54502bSHong Zhang   ilu->reuseordering                    = PETSC_FALSE;
391075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt            = PETSC_DEFAULT;
392075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount       = PETSC_DEFAULT;
393075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol         = PETSC_DEFAULT;
394c378d667SBarry Smith   ((PC_Factor*)ilu)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
3950ed735ceSHong Zhang   ((PC_Factor*)ilu)->info.shiftamount   = 100.0*PETSC_MACHINE_EPSILON;
3960ed735ceSHong Zhang   ((PC_Factor*)ilu)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
397075768bcSBarry Smith   ((PC_Factor*)ilu)->info.pivotinblocks = 1.0;
3989b54502bSHong Zhang   ilu->reusefill                        = PETSC_FALSE;
3990ed735ceSHong Zhang   ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
4009b54502bSHong Zhang   pc->data                              = (void*)ilu;
4019b54502bSHong Zhang 
402574deadeSBarry Smith   pc->ops->reset               = PCReset_ILU;
4039b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
4049b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
4059b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
4069b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
4079b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
40885317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
4099b54502bSHong Zhang   pc->ops->view                = PCView_ILU;
410f0b9ad6cSBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
411f0b9ad6cSBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
4129b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
4139b54502bSHong Zhang 
414bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
415bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
416bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
417bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
418bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
419bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
420bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
421bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
422bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
423bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
424bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
425bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
4262591b318SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
427bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
4288e37d05fSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ILU);CHKERRQ(ierr);
429bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
43092e9c092SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);CHKERRQ(ierr);
431bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
432bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
4339b54502bSHong Zhang   PetscFunctionReturn(0);
4349b54502bSHong Zhang }
435