xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 57622a8ec272e2944ba90348a8875f087fbb2bc6)
185317021SBarry Smith 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h>     /*I "petscpc.h"  I*/
385317021SBarry Smith 
485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
585317021SBarry Smith 
6f8260c8fSBarry Smith 
7f8260c8fSBarry Smith #undef __FUNCT__
8f8260c8fSBarry Smith #define __FUNCT__ "PCFactorSetUpMatSolverPackage_Factor"
9f8260c8fSBarry Smith PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc)
10f8260c8fSBarry Smith {
11f8260c8fSBarry Smith   PC_Factor      *icc = (PC_Factor*)pc->data;
12f8260c8fSBarry Smith   PetscErrorCode ierr;
13f8260c8fSBarry Smith 
14f8260c8fSBarry Smith   PetscFunctionBegin;
15f8260c8fSBarry Smith   if (!pc->setupcalled && !((PC_Factor*)icc)->fact) {
16f8260c8fSBarry Smith     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
17f8260c8fSBarry Smith   }
18f8260c8fSBarry Smith   PetscFunctionReturn(0);
19f8260c8fSBarry Smith }
20f8260c8fSBarry Smith 
2185317021SBarry Smith #undef __FUNCT__
2285317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor"
237087cfbeSBarry Smith PetscErrorCode  PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
2485317021SBarry Smith {
2585317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
2685317021SBarry Smith 
2785317021SBarry Smith   PetscFunctionBegin;
2885317021SBarry Smith   ilu->info.zeropivot = z;
2985317021SBarry Smith   PetscFunctionReturn(0);
3085317021SBarry Smith }
3185317021SBarry Smith 
3285317021SBarry Smith #undef __FUNCT__
33d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType_Factor"
347087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
35d90ac83dSHong Zhang {
36d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
37d90ac83dSHong Zhang 
38d90ac83dSHong Zhang   PetscFunctionBegin;
392fa5cd67SKarl Rupp   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
402fa5cd67SKarl Rupp   else {
41f4db908eSBarry Smith     dir->info.shifttype = (PetscReal) shifttype;
429dbfc187SHong Zhang     if ((shifttype == MAT_SHIFT_NONZERO || shifttype ==  MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) {
439dbfc187SHong Zhang       dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */
446cc3b3c1SHong Zhang     }
45d90ac83dSHong Zhang   }
46d90ac83dSHong Zhang   PetscFunctionReturn(0);
47d90ac83dSHong Zhang }
48d90ac83dSHong Zhang 
49d90ac83dSHong Zhang #undef __FUNCT__
50d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount_Factor"
517087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
52d90ac83dSHong Zhang {
53d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
54d90ac83dSHong Zhang 
55d90ac83dSHong Zhang   PetscFunctionBegin;
562fa5cd67SKarl Rupp   if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
572fa5cd67SKarl Rupp   else dir->info.shiftamount = shiftamount;
58d90ac83dSHong Zhang   PetscFunctionReturn(0);
59d90ac83dSHong Zhang }
60d90ac83dSHong Zhang 
61d90ac83dSHong Zhang #undef __FUNCT__
62b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_Factor"
637087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
6485317021SBarry Smith {
6585317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
6685317021SBarry Smith 
6785317021SBarry Smith   PetscFunctionBegin;
6885317021SBarry Smith   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
69ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
7085317021SBarry Smith   }
7185317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
7285317021SBarry Smith   ilu->info.dt      = dt;
7385317021SBarry Smith   ilu->info.dtcol   = dtcol;
7485317021SBarry Smith   ilu->info.dtcount = dtcount;
7585317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
7685317021SBarry Smith   PetscFunctionReturn(0);
7785317021SBarry Smith }
7885317021SBarry Smith 
7985317021SBarry Smith #undef __FUNCT__
8085317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor"
817087cfbeSBarry Smith PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
8285317021SBarry Smith {
8385317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
8485317021SBarry Smith 
8585317021SBarry Smith   PetscFunctionBegin;
8685317021SBarry Smith   dir->info.fill = fill;
8785317021SBarry Smith   PetscFunctionReturn(0);
8885317021SBarry Smith }
8985317021SBarry Smith 
9085317021SBarry Smith #undef __FUNCT__
9185317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
9219fd82e9SBarry Smith PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering)
9385317021SBarry Smith {
9485317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
9585317021SBarry Smith   PetscErrorCode ierr;
96ace3abfcSBarry Smith   PetscBool      flg;
9785317021SBarry Smith 
9885317021SBarry Smith   PetscFunctionBegin;
9985317021SBarry Smith   if (!pc->setupcalled) {
100503cfb0cSBarry Smith     ierr = PetscFree(dir->ordering);CHKERRQ(ierr);
10119fd82e9SBarry Smith     ierr = PetscStrallocpy(ordering,(char**)&dir->ordering);CHKERRQ(ierr);
10285317021SBarry Smith   } else {
10385317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
104ce94432eSBarry Smith     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
10585317021SBarry Smith   }
10685317021SBarry Smith   PetscFunctionReturn(0);
10785317021SBarry Smith }
10885317021SBarry Smith 
10985317021SBarry Smith #undef __FUNCT__
1102591b318SToby Isaac #define __FUNCT__ "PCFactorGetLevels_Factor"
1112591b318SToby Isaac PetscErrorCode  PCFactorGetLevels_Factor(PC pc,PetscInt *levels)
1122591b318SToby Isaac {
1132591b318SToby Isaac   PC_Factor      *ilu = (PC_Factor*)pc->data;
1142591b318SToby Isaac 
1152591b318SToby Isaac   PetscFunctionBegin;
1162591b318SToby Isaac   *levels = ilu->info.levels;
1172591b318SToby Isaac   PetscFunctionReturn(0);
1182591b318SToby Isaac }
1192591b318SToby Isaac 
1202591b318SToby Isaac #undef __FUNCT__
12185317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor"
1227087cfbeSBarry Smith PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
12385317021SBarry Smith {
12485317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
1255b9c68c7SBarry Smith   PetscErrorCode ierr;
12685317021SBarry Smith 
12785317021SBarry Smith   PetscFunctionBegin;
1282fa5cd67SKarl Rupp   if (!pc->setupcalled) ilu->info.levels = levels;
1292fa5cd67SKarl Rupp   else if (ilu->info.levels != levels) {
1305b9c68c7SBarry Smith     ierr             = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */
1315b9c68c7SBarry Smith     pc->setupcalled  = 0; /* force a complete rebuild of preconditioner factored matrices */
1325b9c68c7SBarry Smith     ilu->info.levels = levels;
133ce94432eSBarry Smith   } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt");
13485317021SBarry Smith   PetscFunctionReturn(0);
13585317021SBarry Smith }
13685317021SBarry Smith 
13785317021SBarry Smith #undef __FUNCT__
13885317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
1397087cfbeSBarry Smith PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc)
14085317021SBarry Smith {
14185317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
14285317021SBarry Smith 
14385317021SBarry Smith   PetscFunctionBegin;
14485317021SBarry Smith   dir->info.diagonal_fill = 1;
14585317021SBarry Smith   PetscFunctionReturn(0);
14685317021SBarry Smith }
14785317021SBarry Smith 
14885317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
14985317021SBarry Smith 
15085317021SBarry Smith #undef __FUNCT__
15185317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
1527087cfbeSBarry Smith PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
15385317021SBarry Smith {
15485317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
15585317021SBarry Smith 
15685317021SBarry Smith   PetscFunctionBegin;
15785317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
15885317021SBarry Smith   PetscFunctionReturn(0);
15985317021SBarry Smith }
16085317021SBarry Smith 
16185317021SBarry Smith #undef __FUNCT__
16285317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor"
1637087cfbeSBarry Smith PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
16485317021SBarry Smith {
16585317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
16685317021SBarry Smith 
16785317021SBarry Smith   PetscFunctionBegin;
168ce94432eSBarry Smith   if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
16985317021SBarry Smith   *mat = ilu->fact;
17085317021SBarry Smith   PetscFunctionReturn(0);
17185317021SBarry Smith }
17285317021SBarry Smith 
17385317021SBarry Smith #undef __FUNCT__
17485317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
1757087cfbeSBarry Smith PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
17685317021SBarry Smith {
17785317021SBarry Smith   PetscErrorCode ierr;
17885317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
17985317021SBarry Smith 
18085317021SBarry Smith   PetscFunctionBegin;
1817112b564SBarry Smith   if (lu->fact) {
18285317021SBarry Smith     const MatSolverPackage ltype;
183ace3abfcSBarry Smith     PetscBool              flg;
18485317021SBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
18585317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
186ce94432eSBarry Smith     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
187a1f19f5aSHong Zhang   } else {
188503cfb0cSBarry Smith     ierr = PetscFree(lu->solvertype);CHKERRQ(ierr);
18985317021SBarry Smith     ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
190f8260c8fSBarry Smith   }
19185317021SBarry Smith   PetscFunctionReturn(0);
19285317021SBarry Smith }
19385317021SBarry Smith 
19485317021SBarry Smith #undef __FUNCT__
1957112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
1967087cfbeSBarry Smith PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
1977112b564SBarry Smith {
1987112b564SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
1997112b564SBarry Smith 
2007112b564SBarry Smith   PetscFunctionBegin;
2017112b564SBarry Smith   *stype = lu->solvertype;
2027112b564SBarry Smith   PetscFunctionReturn(0);
2037112b564SBarry Smith }
2047112b564SBarry Smith 
2057112b564SBarry Smith #undef __FUNCT__
2068ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot_Factor"
2077087cfbeSBarry Smith PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
20885317021SBarry Smith {
20985317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
21085317021SBarry Smith 
21185317021SBarry Smith   PetscFunctionBegin;
212*57622a8eSBarry Smith   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",(double)dtcol);
21385317021SBarry Smith   dir->info.dtcol = dtcol;
21485317021SBarry Smith   PetscFunctionReturn(0);
21585317021SBarry Smith }
2168ff23777SHong Zhang 
2178ff23777SHong Zhang #undef __FUNCT__
2188ff23777SHong Zhang #define __FUNCT__ "PCSetFromOptions_Factor"
2197087cfbeSBarry Smith PetscErrorCode  PCSetFromOptions_Factor(PC pc)
2208ff23777SHong Zhang {
2218ff23777SHong Zhang   PC_Factor         *factor = (PC_Factor*)pc->data;
2228ff23777SHong Zhang   PetscErrorCode    ierr;
223ace3abfcSBarry Smith   PetscBool         flg = PETSC_FALSE,set;
2248ff23777SHong Zhang   char              tname[256], solvertype[64];
225140e18c1SBarry Smith   PetscFunctionList ordlist;
226018dd85eSSatish Balay   PetscEnum         etmp;
2278ff23777SHong Zhang 
2288ff23777SHong Zhang   PetscFunctionBegin;
229607a6623SBarry Smith   if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll();CHKERRQ(ierr);}
2300298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,NULL);CHKERRQ(ierr);
2318ff23777SHong Zhang   if (flg) {
2328ff23777SHong Zhang     ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
2338ff23777SHong Zhang   }
2348ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);CHKERRQ(ierr);
2358ff23777SHong Zhang 
2361ae3d29cSBarry Smith   ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",
237018dd85eSSatish Balay                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr);
238018dd85eSSatish Balay   if (flg) {
2396cc3b3c1SHong Zhang     ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr);
240018dd85eSSatish Balay   }
241d90ac83dSHong Zhang   ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);CHKERRQ(ierr);
242d90ac83dSHong Zhang 
2438ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);CHKERRQ(ierr);
2448ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_column_pivot","Column pivot tolerance (used only for some factorization)","PCFactorSetColumnPivot",((PC_Factor*)factor)->info.dtcol,&((PC_Factor*)factor)->info.dtcol,&flg);CHKERRQ(ierr);
2458ff23777SHong Zhang 
2468ff23777SHong Zhang   flg  = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
247acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
2488ff23777SHong Zhang   if (set) {
2498ff23777SHong Zhang     ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
2508ff23777SHong Zhang   }
2518ff23777SHong Zhang 
2528ff23777SHong Zhang   flg  = PETSC_FALSE;
2530298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,NULL);CHKERRQ(ierr);
2548ff23777SHong Zhang   if (flg) {
2558ff23777SHong Zhang     ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
2568ff23777SHong Zhang   }
2578ff23777SHong Zhang   flg  = PETSC_FALSE;
2580298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,NULL);CHKERRQ(ierr);
2598ff23777SHong Zhang   if (flg) {
2608ff23777SHong Zhang     ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
2618ff23777SHong Zhang   }
2628ff23777SHong Zhang 
2638ff23777SHong Zhang   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
264a264d7a6SBarry Smith   ierr = PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr);
2658ff23777SHong Zhang   if (flg) {
2668ff23777SHong Zhang     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
2678ff23777SHong Zhang   }
2688ff23777SHong Zhang 
2698ff23777SHong Zhang   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
2708ff23777SHong Zhang   ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
2718ff23777SHong Zhang   if (flg) {
2728ff23777SHong Zhang     ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
2738ff23777SHong Zhang   }
2748ff23777SHong Zhang   PetscFunctionReturn(0);
2758ff23777SHong Zhang }
276914a5d51SHong Zhang 
277914a5d51SHong Zhang #undef __FUNCT__
278914a5d51SHong Zhang #define __FUNCT__ "PCView_Factor"
279914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
280914a5d51SHong Zhang {
281914a5d51SHong Zhang   PC_Factor      *factor = (PC_Factor*)pc->data;
282914a5d51SHong Zhang   PetscErrorCode ierr;
283ace3abfcSBarry Smith   PetscBool      isstring,iascii;
284914a5d51SHong Zhang 
285914a5d51SHong Zhang   PetscFunctionBegin;
286251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
287251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
288914a5d51SHong Zhang   if (iascii) {
289879e8a4dSBarry Smith     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
290914a5d51SHong Zhang       if (factor->info.dt > 0) {
291*57622a8eSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  drop tolerance %g\n",(double)factor->info.dt);CHKERRQ(ierr);
292914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr);
293*57622a8eSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %g\n",(double)factor->info.dtcol);CHKERRQ(ierr);
294914a5d51SHong Zhang       } else if (factor->info.levels == 1) {
295914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
296914a5d51SHong Zhang       } else {
297914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
298914a5d51SHong Zhang       }
299914a5d51SHong Zhang     }
300914a5d51SHong Zhang 
301*57622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %g\n",(double)factor->info.zeropivot);CHKERRQ(ierr);
3025e9742b9SJed Brown     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
3035e9742b9SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);CHKERRQ(ierr);
304914a5d51SHong Zhang     }
305914a5d51SHong Zhang 
306914a5d51SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr);
307914a5d51SHong Zhang 
308914a5d51SHong Zhang     if (factor->fact) {
3096335c336SBarry Smith       MatInfo info;
3106335c336SBarry Smith       ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
311*57622a8eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);CHKERRQ(ierr);
312914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");CHKERRQ(ierr);
313914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
314914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
315914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
316914a5d51SHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
317914a5d51SHong Zhang       ierr = MatView(factor->fact,viewer);CHKERRQ(ierr);
318914a5d51SHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
319914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
320914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
321914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
322914a5d51SHong Zhang     }
323914a5d51SHong Zhang 
324914a5d51SHong Zhang   } else if (isstring) {
3256335c336SBarry Smith     MatFactorType t;
3266335c336SBarry Smith     ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr);
3276335c336SBarry Smith     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
328914a5d51SHong Zhang       ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
329914a5d51SHong Zhang     }
330914a5d51SHong Zhang   }
331914a5d51SHong Zhang   PetscFunctionReturn(0);
332914a5d51SHong Zhang }
333