xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision f73b041577f311cd4e522d727240702cd4268ffe)
185317021SBarry Smith 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h>     /*I "petscpc.h"  I*/
385317021SBarry Smith 
485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
585317021SBarry Smith 
6f8260c8fSBarry Smith 
73ca39a21SBarry Smith PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC pc)
8f8260c8fSBarry Smith {
9f8260c8fSBarry Smith   PC_Factor      *icc = (PC_Factor*)pc->data;
10f8260c8fSBarry Smith   PetscErrorCode ierr;
11f8260c8fSBarry Smith 
12f8260c8fSBarry Smith   PetscFunctionBegin;
13ad5697c6SBarry Smith   if (!pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"You can only call this routine after the matrix object has been provided to the solver, for example with KSPSetOperators() or SNESSetJacobian()");
14f8260c8fSBarry Smith   if (!pc->setupcalled && !((PC_Factor*)icc)->fact) {
15f8260c8fSBarry Smith     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
16f8260c8fSBarry Smith   }
17f8260c8fSBarry Smith   PetscFunctionReturn(0);
18f8260c8fSBarry Smith }
19f8260c8fSBarry Smith 
207087cfbeSBarry Smith PetscErrorCode  PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
2185317021SBarry Smith {
2285317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
2385317021SBarry Smith 
2485317021SBarry Smith   PetscFunctionBegin;
2585317021SBarry Smith   ilu->info.zeropivot = z;
2685317021SBarry Smith   PetscFunctionReturn(0);
2785317021SBarry Smith }
2885317021SBarry Smith 
297087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
30d90ac83dSHong Zhang {
31d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
32d90ac83dSHong Zhang 
33d90ac83dSHong Zhang   PetscFunctionBegin;
342fa5cd67SKarl Rupp   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
352fa5cd67SKarl Rupp   else {
36f4db908eSBarry Smith     dir->info.shifttype = (PetscReal) shifttype;
379dbfc187SHong Zhang     if ((shifttype == MAT_SHIFT_NONZERO || shifttype ==  MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) {
389dbfc187SHong Zhang       dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */
396cc3b3c1SHong Zhang     }
40d90ac83dSHong Zhang   }
41d90ac83dSHong Zhang   PetscFunctionReturn(0);
42d90ac83dSHong Zhang }
43d90ac83dSHong Zhang 
447087cfbeSBarry Smith PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
45d90ac83dSHong Zhang {
46d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
47d90ac83dSHong Zhang 
48d90ac83dSHong Zhang   PetscFunctionBegin;
492fa5cd67SKarl Rupp   if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
502fa5cd67SKarl Rupp   else dir->info.shiftamount = shiftamount;
51d90ac83dSHong Zhang   PetscFunctionReturn(0);
52d90ac83dSHong Zhang }
53d90ac83dSHong Zhang 
547087cfbeSBarry Smith PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
5585317021SBarry Smith {
5685317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
5785317021SBarry Smith 
5885317021SBarry Smith   PetscFunctionBegin;
5985317021SBarry 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)) {
60ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
6185317021SBarry Smith   }
6285317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
6385317021SBarry Smith   ilu->info.dt      = dt;
6485317021SBarry Smith   ilu->info.dtcol   = dtcol;
6585317021SBarry Smith   ilu->info.dtcount = dtcount;
6685317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
6785317021SBarry Smith   PetscFunctionReturn(0);
6885317021SBarry Smith }
6985317021SBarry Smith 
707087cfbeSBarry Smith PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
7185317021SBarry Smith {
7285317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
7385317021SBarry Smith 
7485317021SBarry Smith   PetscFunctionBegin;
7585317021SBarry Smith   dir->info.fill = fill;
7685317021SBarry Smith   PetscFunctionReturn(0);
7785317021SBarry Smith }
7885317021SBarry Smith 
7919fd82e9SBarry Smith PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering)
8085317021SBarry Smith {
8185317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
8285317021SBarry Smith   PetscErrorCode ierr;
83ace3abfcSBarry Smith   PetscBool      flg;
8485317021SBarry Smith 
8585317021SBarry Smith   PetscFunctionBegin;
8685317021SBarry Smith   if (!pc->setupcalled) {
87503cfb0cSBarry Smith     ierr = PetscFree(dir->ordering);CHKERRQ(ierr);
8819fd82e9SBarry Smith     ierr = PetscStrallocpy(ordering,(char**)&dir->ordering);CHKERRQ(ierr);
8985317021SBarry Smith   } else {
9085317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
91ce94432eSBarry Smith     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
9285317021SBarry Smith   }
9385317021SBarry Smith   PetscFunctionReturn(0);
9485317021SBarry Smith }
9585317021SBarry Smith 
962591b318SToby Isaac PetscErrorCode  PCFactorGetLevels_Factor(PC pc,PetscInt *levels)
972591b318SToby Isaac {
982591b318SToby Isaac   PC_Factor      *ilu = (PC_Factor*)pc->data;
992591b318SToby Isaac 
1002591b318SToby Isaac   PetscFunctionBegin;
1012591b318SToby Isaac   *levels = ilu->info.levels;
1022591b318SToby Isaac   PetscFunctionReturn(0);
1032591b318SToby Isaac }
1042591b318SToby Isaac 
105c7f610a1SBarry Smith PetscErrorCode  PCFactorGetZeroPivot_Factor(PC pc,PetscReal *pivot)
106c7f610a1SBarry Smith {
107c7f610a1SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
108c7f610a1SBarry Smith 
109c7f610a1SBarry Smith   PetscFunctionBegin;
110c7f610a1SBarry Smith   *pivot = ilu->info.zeropivot;
111c7f610a1SBarry Smith   PetscFunctionReturn(0);
112c7f610a1SBarry Smith }
113c7f610a1SBarry Smith 
114c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftAmount_Factor(PC pc,PetscReal *shift)
115c7f610a1SBarry Smith {
116c7f610a1SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
117c7f610a1SBarry Smith 
118c7f610a1SBarry Smith   PetscFunctionBegin;
119c7f610a1SBarry Smith   *shift = ilu->info.shiftamount;
120c7f610a1SBarry Smith   PetscFunctionReturn(0);
121c7f610a1SBarry Smith }
122c7f610a1SBarry Smith 
123c7f610a1SBarry Smith PetscErrorCode  PCFactorGetShiftType_Factor(PC pc,MatFactorShiftType *type)
124c7f610a1SBarry Smith {
125c7f610a1SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
126c7f610a1SBarry Smith 
127c7f610a1SBarry Smith   PetscFunctionBegin;
128b729e602SBarry Smith   *type = (MatFactorShiftType) (int) ilu->info.shifttype;
129c7f610a1SBarry Smith   PetscFunctionReturn(0);
130c7f610a1SBarry Smith }
131c7f610a1SBarry Smith 
1327087cfbeSBarry Smith PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
13385317021SBarry Smith {
13485317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
1355b9c68c7SBarry Smith   PetscErrorCode ierr;
13685317021SBarry Smith 
13785317021SBarry Smith   PetscFunctionBegin;
1382fa5cd67SKarl Rupp   if (!pc->setupcalled) ilu->info.levels = levels;
1392fa5cd67SKarl Rupp   else if (ilu->info.levels != levels) {
1405b9c68c7SBarry Smith     ierr             = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */
1415b9c68c7SBarry Smith     pc->setupcalled  = 0; /* force a complete rebuild of preconditioner factored matrices */
1425b9c68c7SBarry Smith     ilu->info.levels = levels;
143ce94432eSBarry Smith   } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt");
14485317021SBarry Smith   PetscFunctionReturn(0);
14585317021SBarry Smith }
14685317021SBarry Smith 
14792e9c092SBarry Smith PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg)
14885317021SBarry Smith {
14985317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
15085317021SBarry Smith 
15185317021SBarry Smith   PetscFunctionBegin;
15292e9c092SBarry Smith   dir->info.diagonal_fill = (PetscReal) flg;
15392e9c092SBarry Smith   PetscFunctionReturn(0);
15492e9c092SBarry Smith }
15592e9c092SBarry Smith 
15692e9c092SBarry Smith PetscErrorCode  PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg)
15792e9c092SBarry Smith {
15892e9c092SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
15992e9c092SBarry Smith 
16092e9c092SBarry Smith   PetscFunctionBegin;
16192e9c092SBarry Smith   *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE;
16285317021SBarry Smith   PetscFunctionReturn(0);
16385317021SBarry Smith }
16485317021SBarry Smith 
16585317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
16685317021SBarry Smith 
1677087cfbeSBarry Smith PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
16885317021SBarry Smith {
16985317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
17085317021SBarry Smith 
17185317021SBarry Smith   PetscFunctionBegin;
17285317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
17385317021SBarry Smith   PetscFunctionReturn(0);
17485317021SBarry Smith }
17585317021SBarry Smith 
1767087cfbeSBarry Smith PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
17785317021SBarry Smith {
17885317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
17985317021SBarry Smith 
18085317021SBarry Smith   PetscFunctionBegin;
181ce94432eSBarry Smith   if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
18285317021SBarry Smith   *mat = ilu->fact;
18385317021SBarry Smith   PetscFunctionReturn(0);
18485317021SBarry Smith }
18585317021SBarry Smith 
186ea799195SBarry Smith PetscErrorCode  PCFactorSetMatSolverType_Factor(PC pc,MatSolverType stype)
18785317021SBarry Smith {
18885317021SBarry Smith   PetscErrorCode ierr;
18985317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
19085317021SBarry Smith 
19185317021SBarry Smith   PetscFunctionBegin;
1927112b564SBarry Smith   if (lu->fact) {
193ea799195SBarry Smith     MatSolverType ltype;
194ace3abfcSBarry Smith     PetscBool     flg;
1953ca39a21SBarry Smith     ierr = MatFactorGetSolverType(lu->fact,&ltype);CHKERRQ(ierr);
19685317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
197ce94432eSBarry Smith     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
19800c67f3bSHong Zhang   }
19900c67f3bSHong Zhang 
200503cfb0cSBarry Smith   ierr = PetscFree(lu->solvertype);CHKERRQ(ierr);
20185317021SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
20285317021SBarry Smith   PetscFunctionReturn(0);
20385317021SBarry Smith }
20485317021SBarry Smith 
205ea799195SBarry Smith PetscErrorCode  PCFactorGetMatSolverType_Factor(PC pc,MatSolverType *stype)
2067112b564SBarry Smith {
2077112b564SBarry Smith   PC_Factor *lu = (PC_Factor*)pc->data;
2087112b564SBarry Smith 
2097112b564SBarry Smith   PetscFunctionBegin;
2107112b564SBarry Smith   *stype = lu->solvertype;
2117112b564SBarry Smith   PetscFunctionReturn(0);
2127112b564SBarry Smith }
2137112b564SBarry Smith 
2147087cfbeSBarry Smith PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
21585317021SBarry Smith {
21685317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
21785317021SBarry Smith 
21885317021SBarry Smith   PetscFunctionBegin;
21957622a8eSBarry 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);
22085317021SBarry Smith   dir->info.dtcol = dtcol;
22185317021SBarry Smith   PetscFunctionReturn(0);
22285317021SBarry Smith }
2238ff23777SHong Zhang 
2244416b707SBarry Smith PetscErrorCode  PCSetFromOptions_Factor(PetscOptionItems *PetscOptionsObject,PC pc)
2258ff23777SHong Zhang {
2268ff23777SHong Zhang   PC_Factor         *factor = (PC_Factor*)pc->data;
2278ff23777SHong Zhang   PetscErrorCode    ierr;
2288afaa268SBarry Smith   PetscBool         flg,set;
2298ff23777SHong Zhang   char              tname[256], solvertype[64];
230140e18c1SBarry Smith   PetscFunctionList ordlist;
231018dd85eSSatish Balay   PetscEnum         etmp;
2328e37d05fSBarry Smith   PetscBool         inplace;
2338ff23777SHong Zhang 
2348ff23777SHong Zhang   PetscFunctionBegin;
2358e37d05fSBarry Smith   ierr = PCFactorGetUseInPlace(pc,&inplace);CHKERRQ(ierr);
23692e9c092SBarry Smith   ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set);CHKERRQ(ierr);
2378e37d05fSBarry Smith   if (set) {
2388e37d05fSBarry Smith     ierr = PCFactorSetUseInPlace(pc,flg);CHKERRQ(ierr);
2398ff23777SHong Zhang   }
2408afaa268SBarry Smith   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,NULL);CHKERRQ(ierr);
2418ff23777SHong Zhang 
2428afaa268SBarry Smith   ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr);
243018dd85eSSatish Balay   if (flg) {
2446cc3b3c1SHong Zhang     ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr);
245018dd85eSSatish Balay   }
2460a545947SLisandro Dalcin   ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,NULL);CHKERRQ(ierr);
247d90ac83dSHong Zhang 
2480a545947SLisandro Dalcin   ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,NULL);CHKERRQ(ierr);
2498ff23777SHong 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);
2508ff23777SHong Zhang 
2518afaa268SBarry Smith   ierr = PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
2528ff23777SHong Zhang   if (set) {
2538ff23777SHong Zhang     ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
2548ff23777SHong Zhang   }
2558ff23777SHong Zhang 
2568afaa268SBarry Smith   ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
2578afaa268SBarry Smith   if (set) {
2588afaa268SBarry Smith     ierr = PCFactorSetReuseFill(pc,flg);CHKERRQ(ierr);
2598ff23777SHong Zhang   }
2608afaa268SBarry Smith   ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
2618afaa268SBarry Smith   if (set) {
2628afaa268SBarry Smith     ierr = PCFactorSetReuseOrdering(pc,flg);CHKERRQ(ierr);
2638ff23777SHong Zhang   }
2648ff23777SHong Zhang 
2658ff23777SHong Zhang   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
266589a23caSBarry Smith   ierr = PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,sizeof(tname),&flg);CHKERRQ(ierr);
2678ff23777SHong Zhang   if (flg) {
2688ff23777SHong Zhang     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
2698ff23777SHong Zhang   }
2708ff23777SHong Zhang 
2718ff23777SHong Zhang   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
2729f3a6782SPatrick Sanan   ierr = PetscOptionsDeprecated("-pc_factor_mat_solver_package","-pc_factor_mat_solver_type","3.9",NULL);CHKERRQ(ierr);
273589a23caSBarry Smith   ierr = PetscOptionsString("-pc_factor_mat_solver_type","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,sizeof(solvertype),&flg);CHKERRQ(ierr);
2748ff23777SHong Zhang   if (flg) {
2753ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solvertype);CHKERRQ(ierr);
2768ff23777SHong Zhang   }
2778ff23777SHong Zhang   PetscFunctionReturn(0);
2788ff23777SHong Zhang }
279914a5d51SHong Zhang 
280914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
281914a5d51SHong Zhang {
282914a5d51SHong Zhang   PC_Factor       *factor = (PC_Factor*)pc->data;
283914a5d51SHong Zhang   PetscErrorCode  ierr;
2847a3b9f03SLawrence Mitchell   PetscBool       isstring,iascii,flg;
2859bd791bbSBarry Smith   MatOrderingType ordering;
286914a5d51SHong Zhang 
287914a5d51SHong Zhang   PetscFunctionBegin;
288251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
289251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
290914a5d51SHong Zhang   if (iascii) {
2913d1c1ea0SBarry Smith     if (factor->inplace) {
2923d1c1ea0SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  in-place factorization\n");CHKERRQ(ierr);
2933d1c1ea0SBarry Smith     } else {
2943d1c1ea0SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  out-of-place factorization\n");CHKERRQ(ierr);
2953d1c1ea0SBarry Smith     }
2963d1c1ea0SBarry Smith 
2973d1c1ea0SBarry Smith     if (factor->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"  Reusing fill from past factorization\n");CHKERRQ(ierr);}
2983d1c1ea0SBarry Smith     if (factor->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"  Reusing reordering from past factorization\n");CHKERRQ(ierr);}
299879e8a4dSBarry Smith     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
300914a5d51SHong Zhang       if (factor->info.dt > 0) {
30157622a8eSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  drop tolerance %g\n",(double)factor->info.dt);CHKERRQ(ierr);
302914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr);
30357622a8eSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %g\n",(double)factor->info.dtcol);CHKERRQ(ierr);
304914a5d51SHong Zhang       } else if (factor->info.levels == 1) {
305914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
306914a5d51SHong Zhang       } else {
307914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
308914a5d51SHong Zhang       }
309914a5d51SHong Zhang     }
310914a5d51SHong Zhang 
31157622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %g\n",(double)factor->info.zeropivot);CHKERRQ(ierr);
3125e9742b9SJed Brown     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
3135e9742b9SJed Brown       ierr = PetscViewerASCIIPrintf(viewer,"  using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);CHKERRQ(ierr);
314914a5d51SHong Zhang     }
315914a5d51SHong Zhang 
3169bd791bbSBarry Smith     ierr = PetscStrcmp(factor->ordering,MATORDERINGNATURAL_OR_ND,&flg);CHKERRQ(ierr);
3179bd791bbSBarry Smith     if (flg) {
3189bd791bbSBarry Smith       PetscBool isseqsbaij;
3199bd791bbSBarry Smith       ierr = PetscObjectTypeCompareAny((PetscObject)pc->pmat,&isseqsbaij,MATSEQSBAIJ,MATSEQBAIJ,NULL);CHKERRQ(ierr);
3209bd791bbSBarry Smith       if (isseqsbaij) {
3219bd791bbSBarry Smith         ordering = MATORDERINGNATURAL;
3229bd791bbSBarry Smith       } else {
3239bd791bbSBarry Smith         ordering = MATORDERINGND;
3249bd791bbSBarry Smith       }
3259bd791bbSBarry Smith     } else {
3269bd791bbSBarry Smith       ordering = factor->ordering;
3279bd791bbSBarry Smith     }
3287a3b9f03SLawrence Mitchell     if (!factor->fact) {
3297a3b9f03SLawrence Mitchell       ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s (may be overridden during setup)\n",ordering);CHKERRQ(ierr);
3307a3b9f03SLawrence Mitchell     } else {
331*f73b0415SBarry Smith       PetscBool canuseordering;
3327a3b9f03SLawrence Mitchell       MatInfo   info;
333*f73b0415SBarry Smith       ierr = MatFactorGetCanUseOrdering(factor->fact,&canuseordering);CHKERRQ(ierr);
334*f73b0415SBarry Smith       if (!canuseordering) ordering = "external";
3359bd791bbSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",ordering);CHKERRQ(ierr);
3366335c336SBarry Smith       ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
33757622a8eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);CHKERRQ(ierr);
338914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");CHKERRQ(ierr);
339914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
340914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
341914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
342914a5d51SHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
343914a5d51SHong Zhang       ierr = MatView(factor->fact,viewer);CHKERRQ(ierr);
344914a5d51SHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
345914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
346914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
347914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
348914a5d51SHong Zhang     }
349914a5d51SHong Zhang 
350914a5d51SHong Zhang   } else if (isstring) {
3516335c336SBarry Smith     MatFactorType t;
3526335c336SBarry Smith     ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr);
3536335c336SBarry Smith     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
354c3b366b1Sprj-       ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);
355914a5d51SHong Zhang     }
356914a5d51SHong Zhang   }
357914a5d51SHong Zhang   PetscFunctionReturn(0);
358914a5d51SHong Zhang }
359