xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision acfcf0e5ba359b58e6a8a7af3f239cd7334278e8)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
39b54502bSHong Zhang /*
49b54502bSHong Zhang    Defines a ILU factorization preconditioner for any Mat implementation
59b54502bSHong Zhang */
67c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/ilu/ilu.h"     /*I "petscpc.h"  I*/
79b54502bSHong Zhang 
89b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/
9afaefe49SHong Zhang EXTERN_C_BEGIN
10afaefe49SHong Zhang #undef __FUNCT__
112401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_ILU"
12ace3abfcSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscBool  flag)
132401956bSBarry Smith {
14075768bcSBarry Smith   PC_ILU *lu = (PC_ILU*)pc->data;
152401956bSBarry Smith 
162401956bSBarry Smith   PetscFunctionBegin;
172401956bSBarry Smith   lu->reusefill = flag;
182401956bSBarry Smith   PetscFunctionReturn(0);
192401956bSBarry Smith }
202401956bSBarry Smith EXTERN_C_END
212401956bSBarry Smith 
222401956bSBarry Smith EXTERN_C_BEGIN
232401956bSBarry Smith #undef __FUNCT__
242401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU"
252401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
269b54502bSHong Zhang {
279b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU*)pc->data;
289b54502bSHong Zhang 
299b54502bSHong Zhang   PetscFunctionBegin;
309b54502bSHong Zhang   ilu->nonzerosalongdiagonal = PETSC_TRUE;
319b54502bSHong Zhang   if (z == PETSC_DECIDE) {
329b54502bSHong Zhang     ilu->nonzerosalongdiagonaltol = 1.e-10;
339b54502bSHong Zhang   } else {
349b54502bSHong Zhang     ilu->nonzerosalongdiagonaltol = z;
359b54502bSHong Zhang   }
369b54502bSHong Zhang   PetscFunctionReturn(0);
379b54502bSHong Zhang }
389b54502bSHong Zhang EXTERN_C_END
399b54502bSHong Zhang 
409b54502bSHong Zhang #undef __FUNCT__
419b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal"
429b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc)
439b54502bSHong Zhang {
449b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
459b54502bSHong Zhang   PetscErrorCode ierr;
469b54502bSHong Zhang 
479b54502bSHong Zhang   PetscFunctionBegin;
48075768bcSBarry Smith   if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);}
499b54502bSHong Zhang   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
509b54502bSHong Zhang   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
519b54502bSHong Zhang   PetscFunctionReturn(0);
529b54502bSHong Zhang }
539b54502bSHong Zhang 
549b54502bSHong Zhang EXTERN_C_BEGIN
559b54502bSHong Zhang #undef __FUNCT__
56b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_ILU"
57b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
589b54502bSHong Zhang {
59075768bcSBarry Smith   PC_ILU         *ilu = (PC_ILU*)pc->data;
609b54502bSHong Zhang 
619b54502bSHong Zhang   PetscFunctionBegin;
624c9036c7SBarry Smith   if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
63e7e72b3dSBarry Smith     SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Cannot change drop tolerance after using PC");
649b54502bSHong Zhang   }
65075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt      = dt;
66075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
67075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount = dtcount;
684c9036c7SBarry Smith   ((PC_Factor*)ilu)->info.usedt   = 1.0;
699b54502bSHong Zhang   PetscFunctionReturn(0);
709b54502bSHong Zhang }
719b54502bSHong Zhang EXTERN_C_END
729b54502bSHong Zhang 
739b54502bSHong Zhang EXTERN_C_BEGIN
749b54502bSHong Zhang #undef __FUNCT__
752401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
76ace3abfcSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscBool  flag)
779b54502bSHong Zhang {
78075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU*)pc->data;
799b54502bSHong Zhang 
809b54502bSHong Zhang   PetscFunctionBegin;
819b54502bSHong Zhang   ilu->reuseordering = flag;
829b54502bSHong Zhang   PetscFunctionReturn(0);
839b54502bSHong Zhang }
849b54502bSHong Zhang EXTERN_C_END
859b54502bSHong Zhang 
869b54502bSHong Zhang EXTERN_C_BEGIN
879b54502bSHong Zhang #undef __FUNCT__
882401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
892401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc)
909b54502bSHong Zhang {
91075768bcSBarry Smith   PC_ILU *dir = (PC_ILU*)pc->data;
929b54502bSHong Zhang 
939b54502bSHong Zhang   PetscFunctionBegin;
949b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
959b54502bSHong Zhang   PetscFunctionReturn(0);
969b54502bSHong Zhang }
979b54502bSHong Zhang EXTERN_C_END
989b54502bSHong Zhang 
999b54502bSHong Zhang #undef __FUNCT__
1009b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU"
1019b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc)
1029b54502bSHong Zhang {
1039b54502bSHong Zhang   PetscErrorCode ierr;
10478fc6b22SHong Zhang   PetscInt       itmp;
105ace3abfcSBarry Smith   PetscBool      flg;
1065a586d82SBarry Smith   /* PetscReal      dt[3]; */
1079b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1089b54502bSHong Zhang   PetscReal      tol;
1099b54502bSHong Zhang 
1109b54502bSHong Zhang   PetscFunctionBegin;
1119b54502bSHong Zhang   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
1128ff23777SHong Zhang     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
1138ff23777SHong Zhang 
114075768bcSBarry Smith     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
115075768bcSBarry Smith     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
11690d69ab7SBarry Smith     flg  = PETSC_FALSE;
117*acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
118075768bcSBarry Smith     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
1195a586d82SBarry Smith     /*
120075768bcSBarry Smith     dt[0] = ((PC_Factor*)ilu)->info.dt;
121075768bcSBarry Smith     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
122075768bcSBarry Smith     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
1235a586d82SBarry Smith 
12478fc6b22SHong Zhang     PetscInt       dtmax = 3;
12578fc6b22SHong Zhang     ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance,","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1269b54502bSHong Zhang     if (flg) {
127b7c853c4SBarry Smith       ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1289b54502bSHong Zhang     }
12978fc6b22SHong Zhang     */
1302401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1319b54502bSHong Zhang     if (flg) {
1329b54502bSHong Zhang       tol = PETSC_DECIDE;
1332401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1342401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1359b54502bSHong Zhang     }
1369b54502bSHong Zhang 
1379b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1389b54502bSHong Zhang   PetscFunctionReturn(0);
1399b54502bSHong Zhang }
1409b54502bSHong Zhang 
1419b54502bSHong Zhang #undef __FUNCT__
1429b54502bSHong Zhang #define __FUNCT__ "PCView_ILU"
1439b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
1449b54502bSHong Zhang {
1459b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1469b54502bSHong Zhang   PetscErrorCode ierr;
147ace3abfcSBarry Smith   PetscBool      iascii;
1489b54502bSHong Zhang 
1499b54502bSHong Zhang   PetscFunctionBegin;
1502692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1519b54502bSHong Zhang   if (iascii) {
152914a5d51SHong Zhang     if (ilu->inplace) {
153914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: in-place factorization\n");CHKERRQ(ierr);
1549b54502bSHong Zhang     } else {
155914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: out-of-place factorization\n");CHKERRQ(ierr);
156145b9266SHong Zhang     }
157145b9266SHong Zhang 
158914a5d51SHong Zhang     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing fill from past factorization\n");CHKERRQ(ierr);}
159914a5d51SHong Zhang     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: Reusing reordering from past factorization\n");CHKERRQ(ierr);}
1609b54502bSHong Zhang   }
161914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1629b54502bSHong Zhang   PetscFunctionReturn(0);
1639b54502bSHong Zhang }
1649b54502bSHong Zhang 
1659b54502bSHong Zhang #undef __FUNCT__
1669b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU"
1679b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc)
1689b54502bSHong Zhang {
1699b54502bSHong Zhang   PetscErrorCode ierr;
1709b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
171f3a39becSBarry Smith   MatInfo        info;
172ace3abfcSBarry Smith   PetscBool      flg;
1739b54502bSHong Zhang 
1749b54502bSHong Zhang   PetscFunctionBegin;
17592927226SBarry Smith   /* ugly hack to change default, since it is not support by some matrix types */
17692927226SBarry Smith   if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
17792927226SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);
17892927226SBarry Smith     if (!flg) {
17992927226SBarry Smith       ierr = PetscTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);
18092927226SBarry Smith       if (!flg) {
18192927226SBarry Smith         ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS;
18292927226SBarry Smith         PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO");
18392927226SBarry Smith       }
18492927226SBarry Smith     }
18592927226SBarry Smith   }
18692927226SBarry Smith 
1879b54502bSHong Zhang   if (ilu->inplace) {
188075768bcSBarry Smith     CHKMEMQ;
1899b54502bSHong Zhang     if (!pc->setupcalled) {
1909b54502bSHong Zhang 
1919b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
1929b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
193075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
194efee365bSSatish Balay       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
195efee365bSSatish Balay       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
1969b54502bSHong Zhang     }
1979b54502bSHong Zhang 
1989b54502bSHong Zhang     /* In place ILU only makes sense with fill factor of 1.0 because
1999b54502bSHong Zhang        cannot have levels of fill */
200075768bcSBarry Smith     ((PC_Factor*)ilu)->info.fill          = 1.0;
20175567043SBarry Smith     ((PC_Factor*)ilu)->info.diagonal_fill = 0.0;
202fe97e370SBarry Smith     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKMEMQ;
203075768bcSBarry Smith     ((PC_Factor*)ilu)->fact = pc->pmat;
2049b54502bSHong Zhang   } else {
2059b54502bSHong Zhang     if (!pc->setupcalled) {
2069b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
207075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
208efee365bSSatish Balay       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
209efee365bSSatish Balay       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
2109b54502bSHong Zhang       /*  Remove zeros along diagonal?     */
2119b54502bSHong Zhang       if (ilu->nonzerosalongdiagonal) {
2129b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
2139b54502bSHong Zhang       }
21486daa42cSBarry Smith       if (!((PC_Factor*)ilu)->fact){
21511bfe483SHong Zhang         ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
216a1f19f5aSHong Zhang       }
217075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
218075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
219f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
220075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2219b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
2229b54502bSHong Zhang       if (!ilu->reuseordering) {
2239b54502bSHong Zhang         /* compute a new ordering for the ILU */
2249b54502bSHong Zhang         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
2259b54502bSHong Zhang         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
226075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
227efee365bSSatish Balay         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
228efee365bSSatish Balay         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
2299b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
2309b54502bSHong Zhang         if (ilu->nonzerosalongdiagonal) {
2319b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
2329b54502bSHong Zhang         }
2339b54502bSHong Zhang       }
234075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2352692d6eeSBarry Smith       ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
236075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
237075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
238f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
239075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2409b54502bSHong Zhang     }
241075768bcSBarry Smith     CHKMEMQ;
242075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
243075768bcSBarry Smith     CHKMEMQ;
2449b54502bSHong Zhang   }
2459b54502bSHong Zhang   PetscFunctionReturn(0);
2469b54502bSHong Zhang }
2479b54502bSHong Zhang 
2489b54502bSHong Zhang #undef __FUNCT__
2499b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU"
2509b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc)
2519b54502bSHong Zhang {
2529b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2539b54502bSHong Zhang   PetscErrorCode ierr;
2549b54502bSHong Zhang 
2559b54502bSHong Zhang   PetscFunctionBegin;
2569b54502bSHong Zhang   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
257503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
258503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
2599b54502bSHong Zhang   ierr = PetscFree(ilu);CHKERRQ(ierr);
2609b54502bSHong Zhang   PetscFunctionReturn(0);
2619b54502bSHong Zhang }
2629b54502bSHong Zhang 
2639b54502bSHong Zhang #undef __FUNCT__
2649b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU"
2659b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
2669b54502bSHong Zhang {
2679b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2689b54502bSHong Zhang   PetscErrorCode ierr;
2699b54502bSHong Zhang 
2709b54502bSHong Zhang   PetscFunctionBegin;
271075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2729b54502bSHong Zhang   PetscFunctionReturn(0);
2739b54502bSHong Zhang }
2749b54502bSHong Zhang 
2759b54502bSHong Zhang #undef __FUNCT__
2769b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU"
2779b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
2789b54502bSHong Zhang {
2799b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
2809b54502bSHong Zhang   PetscErrorCode ierr;
2819b54502bSHong Zhang 
2829b54502bSHong Zhang   PetscFunctionBegin;
283075768bcSBarry Smith   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
2849b54502bSHong Zhang   PetscFunctionReturn(0);
2859b54502bSHong Zhang }
2869b54502bSHong Zhang 
287f0b9ad6cSBarry Smith #undef __FUNCT__
288f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricLeft_ILU"
289f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y)
290f0b9ad6cSBarry Smith {
291f0b9ad6cSBarry Smith   PetscErrorCode ierr;
292f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
293f0b9ad6cSBarry Smith 
294f0b9ad6cSBarry Smith   PetscFunctionBegin;
295f0b9ad6cSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
296f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
297f0b9ad6cSBarry Smith }
298f0b9ad6cSBarry Smith 
299f0b9ad6cSBarry Smith #undef __FUNCT__
300f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricRight_ILU"
301f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y)
302f0b9ad6cSBarry Smith {
303f0b9ad6cSBarry Smith   PetscErrorCode ierr;
304f0b9ad6cSBarry Smith   PC_ILU         *icc = (PC_ILU*)pc->data;
305f0b9ad6cSBarry Smith 
306f0b9ad6cSBarry Smith   PetscFunctionBegin;
307f0b9ad6cSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
308f0b9ad6cSBarry Smith   PetscFunctionReturn(0);
309f0b9ad6cSBarry Smith }
310f0b9ad6cSBarry Smith 
3119b54502bSHong Zhang /*MC
3129b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
3139b54502bSHong Zhang 
3149b54502bSHong Zhang    Options Database Keys:
3152401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
3162401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
3179b54502bSHong Zhang                       its factorization (overwrites original matrix)
3182401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
3192401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
32055ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
3212401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
3229b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
3232401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
3242401956bSBarry Smith .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
3259b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
3269b54502bSHong Zhang                              stability of the ILU factorization
3279b54502bSHong Zhang 
3289b54502bSHong Zhang    Level: beginner
3299b54502bSHong Zhang 
3309b54502bSHong Zhang   Concepts: incomplete factorization
3319b54502bSHong Zhang 
332d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU)
3339b54502bSHong Zhang 
3349b54502bSHong Zhang           For BAIJ matrices this implements a point block ILU
3359b54502bSHong Zhang 
336f0b9ad6cSBarry Smith           The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U)
337f0b9ad6cSBarry Smith           even when the matrix is not symmetric since the U stores the diagonals of the factorization.
338f0b9ad6cSBarry Smith 
339c582cd25SBarry Smith    References:
340c582cd25SBarry Smith    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
341c582cd25SBarry Smith    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
342c582cd25SBarry Smith 
343c582cd25SBarry Smith    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
344c582cd25SBarry Smith    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
345c582cd25SBarry Smith 
346c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
347c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
348c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
349c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
350c582cd25SBarry Smith 
351c582cd25SBarry Smith 
3529b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
353145b9266SHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(),
354b7c853c4SBarry Smith            PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
3558ff23777SHong Zhang            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks()
3569b54502bSHong Zhang 
3579b54502bSHong Zhang M*/
3589b54502bSHong Zhang 
3599b54502bSHong Zhang EXTERN_C_BEGIN
3609b54502bSHong Zhang #undef __FUNCT__
3619b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU"
362dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
3639b54502bSHong Zhang {
3649b54502bSHong Zhang   PetscErrorCode ierr;
3659b54502bSHong Zhang   PC_ILU         *ilu;
3669b54502bSHong Zhang 
3679b54502bSHong Zhang   PetscFunctionBegin;
36838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
3699b54502bSHong Zhang 
370075768bcSBarry Smith   ((PC_Factor*)ilu)->fact                    = 0;
371075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
372879e8a4dSBarry Smith   ((PC_Factor*)ilu)->factortype              = MAT_FACTOR_ILU;
37375567043SBarry Smith   ((PC_Factor*)ilu)->info.levels             = 0.;
374075768bcSBarry Smith   ((PC_Factor*)ilu)->info.fill               = 1.0;
3759b54502bSHong Zhang   ilu->col                     = 0;
3769b54502bSHong Zhang   ilu->row                     = 0;
3779b54502bSHong Zhang   ilu->inplace                 = PETSC_FALSE;
3782692d6eeSBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
3792692d6eeSBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
3809b54502bSHong Zhang   ilu->reuseordering           = PETSC_FALSE;
381075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
382075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
383075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
384357fed3dSBarry Smith   ((PC_Factor*)ilu)->info.shifttype          = (PetscReal)MAT_SHIFT_NONZERO;
385915743fcSHong Zhang   ((PC_Factor*)ilu)->info.shiftamount        = 1.e-12;
386075768bcSBarry Smith   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
387075768bcSBarry Smith   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
3889b54502bSHong Zhang   ilu->reusefill               = PETSC_FALSE;
38975567043SBarry Smith   ((PC_Factor*)ilu)->info.diagonal_fill      = 0.;
3909b54502bSHong Zhang   pc->data                     = (void*)ilu;
3919b54502bSHong Zhang 
3929b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
3939b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
3949b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
3959b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
3969b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
39785317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
3989b54502bSHong Zhang   pc->ops->view                = PCView_ILU;
399f0b9ad6cSBarry Smith   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ILU;
400f0b9ad6cSBarry Smith   pc->ops->applysymmetricright = PCApplySymmetricRight_ILU;
4019b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
4029b54502bSHong Zhang 
40385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
40485317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
405d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
406d90ac83dSHong Zhang                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
407d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
408d90ac83dSHong Zhang                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
4097112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
4107112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
41111bfe483SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
41211bfe483SHong Zhang                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
413b7c853c4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
414b7c853c4SBarry Smith                     PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
41585317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
41685317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
41785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
41885317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
4192401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
4202401956bSBarry Smith                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
4212401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
4222401956bSBarry Smith                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
42385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
42485317021SBarry Smith                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
4252401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
4262401956bSBarry Smith                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
42785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
42885317021SBarry Smith                     PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
42985317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
43085317021SBarry Smith                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
4312401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
4322401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
4339b54502bSHong Zhang   PetscFunctionReturn(0);
4349b54502bSHong Zhang }
4359b54502bSHong Zhang EXTERN_C_END
436