xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 11bfe4838c662b8f8aa3d085a86728ee2cd97e24)
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"
122401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth 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__
562401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU"
572401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
589b54502bSHong Zhang {
59075768bcSBarry Smith   PC_ILU         *ilu = (PC_ILU*)pc->data;
609b54502bSHong Zhang   PetscErrorCode ierr;
619b54502bSHong Zhang 
629b54502bSHong Zhang   PetscFunctionBegin;
63075768bcSBarry Smith   if (pc->setupcalled && (!ilu->usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
649b54502bSHong Zhang     pc->setupcalled   = 0;
659b54502bSHong Zhang     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
669b54502bSHong Zhang   }
679b54502bSHong Zhang   ilu->usedt                      = PETSC_TRUE;
68075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt      = dt;
69075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol   = dtcol;
70075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount = dtcount;
71075768bcSBarry Smith   ((PC_Factor*)ilu)->info.fill    = PETSC_DEFAULT;
729b54502bSHong Zhang   PetscFunctionReturn(0);
739b54502bSHong Zhang }
749b54502bSHong Zhang EXTERN_C_END
759b54502bSHong Zhang 
769b54502bSHong Zhang EXTERN_C_BEGIN
779b54502bSHong Zhang #undef __FUNCT__
782401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU"
792401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag)
809b54502bSHong Zhang {
81075768bcSBarry Smith   PC_ILU *ilu = (PC_ILU*)pc->data;
829b54502bSHong Zhang 
839b54502bSHong Zhang   PetscFunctionBegin;
849b54502bSHong Zhang   ilu->reuseordering = flag;
859b54502bSHong Zhang   PetscFunctionReturn(0);
869b54502bSHong Zhang }
879b54502bSHong Zhang EXTERN_C_END
889b54502bSHong Zhang 
899b54502bSHong Zhang EXTERN_C_BEGIN
909b54502bSHong Zhang #undef __FUNCT__
912401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU"
922401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc)
939b54502bSHong Zhang {
94075768bcSBarry Smith   PC_ILU *dir = (PC_ILU*)pc->data;
959b54502bSHong Zhang 
969b54502bSHong Zhang   PetscFunctionBegin;
979b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
989b54502bSHong Zhang   PetscFunctionReturn(0);
999b54502bSHong Zhang }
1009b54502bSHong Zhang EXTERN_C_END
1019b54502bSHong Zhang 
1029b54502bSHong Zhang #undef __FUNCT__
1039b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU"
1049b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc)
1059b54502bSHong Zhang {
1069b54502bSHong Zhang   PetscErrorCode ierr;
1079b54502bSHong Zhang   PetscInt       dtmax = 3,itmp;
1089b54502bSHong Zhang   PetscTruth     flg,set;
1099b54502bSHong Zhang   PetscReal      dt[3];
110*11bfe483SHong Zhang   char           tname[256], solvertype[64];
1119b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1129b54502bSHong Zhang   PetscFList     ordlist;
1139b54502bSHong Zhang   PetscReal      tol;
1149b54502bSHong Zhang 
1159b54502bSHong Zhang   PetscFunctionBegin;
116cce0b1b2SLisandro Dalcin   if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1179b54502bSHong Zhang   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
118075768bcSBarry Smith     ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr);
119075768bcSBarry Smith     if (flg) ((PC_Factor*)ilu)->info.levels = itmp;
12090d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",ilu->inplace,&ilu->inplace,PETSC_NULL);CHKERRQ(ierr);
12190d69ab7SBarry Smith     flg  = PETSC_FALSE;
12290d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
123075768bcSBarry Smith     ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg;
12490d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",ilu->reusefill,&ilu->reusefill,PETSC_NULL);CHKERRQ(ierr);
12590d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",ilu->reuseordering,&ilu->reuseordering,PETSC_NULL);CHKERRQ(ierr);
12690d69ab7SBarry Smith     flg  = PETSC_FALSE;
12790d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1289b54502bSHong Zhang     if (flg) {
129afaefe49SHong Zhang       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
1309b54502bSHong Zhang     }
131075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)ilu)->info.shiftnz,&((PC_Factor*)ilu)->info.shiftnz,0);CHKERRQ(ierr);
132075768bcSBarry Smith     flg = (((PC_Factor*)ilu)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE;
13362bba022SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
13462bba022SBarry Smith     ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr);
135075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)ilu)->info.zeropivot,&((PC_Factor*)ilu)->info.zeropivot,0);CHKERRQ(ierr);
1369b54502bSHong Zhang 
137075768bcSBarry Smith     dt[0] = ((PC_Factor*)ilu)->info.dt;
138075768bcSBarry Smith     dt[1] = ((PC_Factor*)ilu)->info.dtcol;
139075768bcSBarry Smith     dt[2] = ((PC_Factor*)ilu)->info.dtcount;
1402401956bSBarry Smith     ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1419b54502bSHong Zhang     if (flg) {
1422401956bSBarry Smith       ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1439b54502bSHong Zhang     }
144075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)ilu)->info.fill,&((PC_Factor*)ilu)->info.fill,&flg);CHKERRQ(ierr);
1452401956bSBarry Smith     ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
1469b54502bSHong Zhang     if (flg) {
1479b54502bSHong Zhang       tol = PETSC_DECIDE;
1482401956bSBarry Smith       ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
1492401956bSBarry Smith       ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
1509b54502bSHong Zhang     }
1519b54502bSHong Zhang 
1529b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
153075768bcSBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)ilu)->ordering,tname,256,&flg);CHKERRQ(ierr);
1549b54502bSHong Zhang     if (flg) {
155e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
1569b54502bSHong Zhang     }
157*11bfe483SHong Zhang 
158*11bfe483SHong Zhang     /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
159*11bfe483SHong Zhang     ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific ILU solver to use","MatGetFactor",((PC_Factor*)ilu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
160*11bfe483SHong Zhang     if (flg) {
161*11bfe483SHong Zhang       ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
162*11bfe483SHong Zhang     }
163*11bfe483SHong Zhang 
164075768bcSBarry Smith     flg = ((PC_Factor*)ilu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
1652401956bSBarry Smith     ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
1669b54502bSHong Zhang     if (set) {
1672401956bSBarry Smith       ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
1689b54502bSHong Zhang     }
16990d69ab7SBarry Smith     flg  = PETSC_FALSE;
17090d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
17162bba022SBarry Smith     if (flg) {
17262bba022SBarry Smith       ierr = PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
17362bba022SBarry Smith     }
174075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",((PC_Factor*)ilu)->info.shiftinblocks,&((PC_Factor*)ilu)->info.shiftinblocks,0);CHKERRQ(ierr);
17562bba022SBarry Smith 
1769b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1779b54502bSHong Zhang   PetscFunctionReturn(0);
1789b54502bSHong Zhang }
1799b54502bSHong Zhang 
1809b54502bSHong Zhang #undef __FUNCT__
1819b54502bSHong Zhang #define __FUNCT__ "PCView_ILU"
1829b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
1839b54502bSHong Zhang {
1849b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
1859b54502bSHong Zhang   PetscErrorCode ierr;
1869b54502bSHong Zhang   PetscTruth     isstring,iascii;
1879b54502bSHong Zhang 
1889b54502bSHong Zhang   PetscFunctionBegin;
1899b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1909b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1919b54502bSHong Zhang   if (iascii) {
1929b54502bSHong Zhang     if (ilu->usedt) {
193075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %G\n",((PC_Factor*)ilu)->info.dt);CHKERRQ(ierr);
194075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)((PC_Factor*)ilu)->info.dtcount);CHKERRQ(ierr);
195075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %G\n",((PC_Factor*)ilu)->info.dtcol);CHKERRQ(ierr);
196075768bcSBarry Smith     } else if (((PC_Factor*)ilu)->info.levels == 1) {
197075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr);
1989b54502bSHong Zhang     } else {
199075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr);
2009b54502bSHong Zhang     }
201075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio allocated %G\n",((PC_Factor*)ilu)->info.fill);CHKERRQ(ierr);
202075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %G\n",((PC_Factor*)ilu)->info.zeropivot);CHKERRQ(ierr);
203075768bcSBarry Smith     if (((PC_Factor*)ilu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
204075768bcSBarry Smith     if (((PC_Factor*)ilu)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
205075768bcSBarry Smith     if (((PC_Factor*)ilu)->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);}
2069b54502bSHong Zhang     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
2079b54502bSHong Zhang     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
208075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
2099b54502bSHong Zhang     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
2109b54502bSHong Zhang     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
211075768bcSBarry Smith     if (((PC_Factor*)ilu)->fact) {
212f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr);
2139b54502bSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
2149b54502bSHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
215f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
216f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2179b54502bSHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
218075768bcSBarry Smith       ierr = MatView(((PC_Factor*)ilu)->fact,viewer);CHKERRQ(ierr);
2199b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
2209b54502bSHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
221f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
222f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2239b54502bSHong Zhang     }
2249b54502bSHong Zhang   } else if (isstring) {
225075768bcSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
2269b54502bSHong Zhang   } else {
2279b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
2289b54502bSHong Zhang   }
2299b54502bSHong Zhang   PetscFunctionReturn(0);
2309b54502bSHong Zhang }
2319b54502bSHong Zhang 
2329b54502bSHong Zhang #undef __FUNCT__
2339b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU"
2349b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc)
2359b54502bSHong Zhang {
2369b54502bSHong Zhang   PetscErrorCode ierr;
2379b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
238f3a39becSBarry Smith   MatInfo        info;
2399b54502bSHong Zhang 
2409b54502bSHong Zhang   PetscFunctionBegin;
2419b54502bSHong Zhang   if (ilu->inplace) {
242075768bcSBarry Smith     CHKMEMQ;
2439b54502bSHong Zhang     if (!pc->setupcalled) {
2449b54502bSHong Zhang 
2459b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
2469b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
247075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
248efee365bSSatish Balay       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
249efee365bSSatish Balay       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
2509b54502bSHong Zhang     }
2519b54502bSHong Zhang 
2529b54502bSHong Zhang     /* In place ILU only makes sense with fill factor of 1.0 because
2539b54502bSHong Zhang        cannot have levels of fill */
254075768bcSBarry Smith     ((PC_Factor*)ilu)->info.fill          = 1.0;
255075768bcSBarry Smith     ((PC_Factor*)ilu)->info.diagonal_fill = 0;
256075768bcSBarry Smith     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
257075768bcSBarry Smith     CHKMEMQ;
258075768bcSBarry Smith     ((PC_Factor*)ilu)->fact = pc->pmat;
2599b54502bSHong Zhang   } else if (ilu->usedt) {
2609b54502bSHong Zhang     if (!pc->setupcalled) {
261075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
262075768bcSBarry Smith     CHKMEMQ;
2631da05e5aSHong Zhang       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
264efee365bSSatish Balay       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
265efee365bSSatish Balay       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
266075768bcSBarry Smith       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
267075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2689b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
269075768bcSBarry Smith     CHKMEMQ;
270075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
271075768bcSBarry Smith     CHKMEMQ;
2729b54502bSHong Zhang       if (!ilu->reuseordering) {
2739b54502bSHong Zhang         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
2749b54502bSHong Zhang         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
275075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
276efee365bSSatish Balay         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
277efee365bSSatish Balay         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
2789b54502bSHong Zhang       }
279b78a477dSHong Zhang       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
280075768bcSBarry Smith       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
281075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2829b54502bSHong Zhang     } else if (!ilu->reusefill) {
283075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
284b78a477dSHong Zhang       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILUDT,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
285075768bcSBarry Smith       ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
286075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
2879b54502bSHong Zhang     } else {
288075768bcSBarry Smith       ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
2899b54502bSHong Zhang     }
2909b54502bSHong Zhang   } else {
2919b54502bSHong Zhang     if (!pc->setupcalled) {
2929b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
293075768bcSBarry Smith       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
294efee365bSSatish Balay       if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
295efee365bSSatish Balay       if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
2969b54502bSHong Zhang       /*  Remove zeros along diagonal?     */
2979b54502bSHong Zhang       if (ilu->nonzerosalongdiagonal) {
2989b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
2999b54502bSHong Zhang       }
300075768bcSBarry Smith     CHKMEMQ;
301*11bfe483SHong Zhang       ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
302075768bcSBarry Smith     CHKMEMQ;
303075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
304075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
305f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
306075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
3079b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
3089b54502bSHong Zhang       if (!ilu->reuseordering) {
3099b54502bSHong Zhang         /* compute a new ordering for the ILU */
3109b54502bSHong Zhang         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
3119b54502bSHong Zhang         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
312075768bcSBarry Smith         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
313efee365bSSatish Balay         if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);}
314efee365bSSatish Balay         if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);}
3159b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
3169b54502bSHong Zhang         if (ilu->nonzerosalongdiagonal) {
3179b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
3189b54502bSHong Zhang         }
3199b54502bSHong Zhang       }
320075768bcSBarry Smith       ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
321075768bcSBarry Smith       ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
322075768bcSBarry Smith       ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
323075768bcSBarry Smith       ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
324f3a39becSBarry Smith       ilu->actualfill = info.fill_ratio_needed;
325075768bcSBarry Smith       ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr);
3269b54502bSHong Zhang     }
327075768bcSBarry Smith     CHKMEMQ;
328075768bcSBarry Smith     ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
329075768bcSBarry Smith     CHKMEMQ;
3309b54502bSHong Zhang   }
3319b54502bSHong Zhang   PetscFunctionReturn(0);
3329b54502bSHong Zhang }
3339b54502bSHong Zhang 
3349b54502bSHong Zhang #undef __FUNCT__
3359b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU"
3369b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc)
3379b54502bSHong Zhang {
3389b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
3399b54502bSHong Zhang   PetscErrorCode ierr;
3409b54502bSHong Zhang 
3419b54502bSHong Zhang   PetscFunctionBegin;
3429b54502bSHong Zhang   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
343075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
3449b54502bSHong Zhang   ierr = PetscFree(ilu);CHKERRQ(ierr);
3459b54502bSHong Zhang   PetscFunctionReturn(0);
3469b54502bSHong Zhang }
3479b54502bSHong Zhang 
3489b54502bSHong Zhang #undef __FUNCT__
3499b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU"
3509b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
3519b54502bSHong Zhang {
3529b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
3539b54502bSHong Zhang   PetscErrorCode ierr;
3549b54502bSHong Zhang 
3559b54502bSHong Zhang   PetscFunctionBegin;
356075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
3579b54502bSHong Zhang   PetscFunctionReturn(0);
3589b54502bSHong Zhang }
3599b54502bSHong Zhang 
3609b54502bSHong Zhang #undef __FUNCT__
3619b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU"
3629b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
3639b54502bSHong Zhang {
3649b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
3659b54502bSHong Zhang   PetscErrorCode ierr;
3669b54502bSHong Zhang 
3679b54502bSHong Zhang   PetscFunctionBegin;
368075768bcSBarry Smith   ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr);
3699b54502bSHong Zhang   PetscFunctionReturn(0);
3709b54502bSHong Zhang }
3719b54502bSHong Zhang 
3729b54502bSHong Zhang /*MC
3739b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
3749b54502bSHong Zhang 
3759b54502bSHong Zhang    Options Database Keys:
3762401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ILU(k)
3772401956bSBarry Smith .  -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
3789b54502bSHong Zhang                       its factorization (overwrites original matrix)
3792401956bSBarry Smith .  -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
3802401956bSBarry Smith .  -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization
3812401956bSBarry Smith .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
38255ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
3832401956bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
3849b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
3852401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
3862401956bSBarry Smith .  -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
3879b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
3889b54502bSHong Zhang                              stability of the ILU factorization
38962bba022SBarry Smith .  -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization
390f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
39153ae6812SBarry Smith -  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
39253ae6812SBarry Smith    is optional with true being the default
3939b54502bSHong Zhang 
3949b54502bSHong Zhang    Level: beginner
3959b54502bSHong Zhang 
3969b54502bSHong Zhang   Concepts: incomplete factorization
3979b54502bSHong Zhang 
398c582cd25SBarry Smith    Notes: Only implemented for some matrix formats. (for parallel use you
399c582cd25SBarry Smith              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended
400c582cd25SBarry Smith              unless you really want a parallel ILU).
4019b54502bSHong Zhang 
4029b54502bSHong Zhang           For BAIJ matrices this implements a point block ILU
4039b54502bSHong Zhang 
404c582cd25SBarry Smith    References:
405c582cd25SBarry Smith    T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving
406c582cd25SBarry Smith    self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968.
407c582cd25SBarry Smith 
408c582cd25SBarry Smith    T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif-
409c582cd25SBarry Smith    fusion problems. Quart. Appl. Math., 19:221--229, 1961.
410c582cd25SBarry Smith 
411c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
412c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
413c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
414c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
415c582cd25SBarry Smith 
416c582cd25SBarry Smith 
4179b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
4182401956bSBarry Smith            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(),
419e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
4202401956bSBarry Smith            PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(),
421f251bdbdSHong Zhang            PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
4229b54502bSHong Zhang 
4239b54502bSHong Zhang M*/
4249b54502bSHong Zhang 
4259b54502bSHong Zhang EXTERN_C_BEGIN
4269b54502bSHong Zhang #undef __FUNCT__
4279b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU"
428dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc)
4299b54502bSHong Zhang {
4309b54502bSHong Zhang   PetscErrorCode ierr;
4319b54502bSHong Zhang   PC_ILU         *ilu;
4329b54502bSHong Zhang 
4339b54502bSHong Zhang   PetscFunctionBegin;
43438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr);
4359b54502bSHong Zhang 
436075768bcSBarry Smith   ((PC_Factor*)ilu)->fact                    = 0;
437075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr);
438075768bcSBarry Smith   ((PC_Factor*)ilu)->info.levels             = 0;
439075768bcSBarry Smith   ((PC_Factor*)ilu)->info.fill               = 1.0;
4409b54502bSHong Zhang   ilu->col                     = 0;
4419b54502bSHong Zhang   ilu->row                     = 0;
4429b54502bSHong Zhang   ilu->inplace                 = PETSC_FALSE;
443*11bfe483SHong Zhang   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr);
444075768bcSBarry Smith   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);
4459b54502bSHong Zhang   ilu->reuseordering           = PETSC_FALSE;
4469b54502bSHong Zhang   ilu->usedt                   = PETSC_FALSE;
447075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dt                 = PETSC_DEFAULT;
448075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcount            = PETSC_DEFAULT;
449075768bcSBarry Smith   ((PC_Factor*)ilu)->info.dtcol              = PETSC_DEFAULT;
450075768bcSBarry Smith   ((PC_Factor*)ilu)->info.shiftnz            = 1.e-12;
451075768bcSBarry Smith   ((PC_Factor*)ilu)->info.shiftpd            = 0.0; /* false */
452075768bcSBarry Smith   ((PC_Factor*)ilu)->info.zeropivot          = 1.e-12;
453075768bcSBarry Smith   ((PC_Factor*)ilu)->info.pivotinblocks      = 1.0;
454075768bcSBarry Smith   ((PC_Factor*)ilu)->info.shiftinblocks      = 1.e-12;
4559b54502bSHong Zhang   ilu->reusefill               = PETSC_FALSE;
456075768bcSBarry Smith   ((PC_Factor*)ilu)->info.diagonal_fill      = 0;
4579b54502bSHong Zhang   pc->data                     = (void*)ilu;
4589b54502bSHong Zhang 
4599b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
4609b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
4619b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
4629b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
4639b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
46485317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
4659b54502bSHong Zhang   pc->ops->view                = PCView_ILU;
4669b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
4679b54502bSHong Zhang 
46885317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
46985317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
47085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
47185317021SBarry Smith                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
47285317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
47385317021SBarry Smith                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
474afaefe49SHong Zhang 
4757112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
4767112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
477*11bfe483SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
478*11bfe483SHong Zhang                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
4792401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
4802401956bSBarry Smith                     PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr);
48185317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
48285317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
48385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
48485317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
4852401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU",
4862401956bSBarry Smith                     PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr);
4872401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU",
4882401956bSBarry Smith                     PCFactorSetReuseFill_ILU);CHKERRQ(ierr);
48985317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
49085317021SBarry Smith                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
4912401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU",
4922401956bSBarry Smith                     PCFactorSetUseInPlace_ILU);CHKERRQ(ierr);
49385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_Factor",
49485317021SBarry Smith                     PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr);
49585317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor",
49685317021SBarry Smith                     PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
49785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
49885317021SBarry Smith                     PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr);
4992401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU",
5002401956bSBarry Smith                     PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
5019b54502bSHong Zhang   PetscFunctionReturn(0);
5029b54502bSHong Zhang }
5039b54502bSHong Zhang EXTERN_C_END
504