xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 018dd85ee6e4adaaa90244c554663d0d584c3b2f)
185317021SBarry Smith #define PETSCKSP_DLL
285317021SBarry Smith 
37c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/factor.h"     /*I "petscpc.h"  I*/
485317021SBarry Smith 
585317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
685317021SBarry Smith 
785317021SBarry Smith EXTERN_C_BEGIN
885317021SBarry Smith #undef __FUNCT__
985317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor"
1085317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
1185317021SBarry Smith {
1285317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
1385317021SBarry Smith 
1485317021SBarry Smith   PetscFunctionBegin;
1585317021SBarry Smith   ilu->info.zeropivot = z;
1685317021SBarry Smith   PetscFunctionReturn(0);
1785317021SBarry Smith }
1885317021SBarry Smith EXTERN_C_END
1985317021SBarry Smith 
2085317021SBarry Smith EXTERN_C_BEGIN
2185317021SBarry Smith #undef __FUNCT__
22d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType_Factor"
23d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
24d90ac83dSHong Zhang {
25d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
26d90ac83dSHong Zhang 
27d90ac83dSHong Zhang   PetscFunctionBegin;
28d90ac83dSHong Zhang   if (shifttype == (MatFactorShiftType)PETSC_DECIDE){
29f4db908eSBarry Smith     dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
30d90ac83dSHong Zhang   } else {
31f4db908eSBarry Smith     dir->info.shifttype = (PetscReal) shifttype;
32d90ac83dSHong Zhang   }
33d90ac83dSHong Zhang   PetscFunctionReturn(0);
34d90ac83dSHong Zhang }
35d90ac83dSHong Zhang 
36d90ac83dSHong Zhang #undef __FUNCT__
37d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount_Factor"
38d90ac83dSHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
39d90ac83dSHong Zhang {
40d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor*)pc->data;
41d90ac83dSHong Zhang 
42d90ac83dSHong Zhang   PetscFunctionBegin;
43d90ac83dSHong Zhang   if (shiftamount == (PetscReal) PETSC_DECIDE){
44d90ac83dSHong Zhang     dir->info.shiftamount = 1.e-12;
45d90ac83dSHong Zhang   } else {
46d90ac83dSHong Zhang     dir->info.shiftamount = shiftamount;
47d90ac83dSHong Zhang   }
48d90ac83dSHong Zhang   PetscFunctionReturn(0);
49d90ac83dSHong Zhang }
50d90ac83dSHong Zhang EXTERN_C_END
51d90ac83dSHong Zhang 
52d90ac83dSHong Zhang EXTERN_C_BEGIN
53d90ac83dSHong Zhang #undef __FUNCT__
54b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_Factor"
55b7c853c4SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
5685317021SBarry Smith {
5785317021SBarry Smith   PC_Factor         *ilu = (PC_Factor*)pc->data;
5885317021SBarry Smith 
5985317021SBarry Smith   PetscFunctionBegin;
6085317021SBarry 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)) {
6185317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
6285317021SBarry Smith   }
6385317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
6485317021SBarry Smith   ilu->info.dt      = dt;
6585317021SBarry Smith   ilu->info.dtcol   = dtcol;
6685317021SBarry Smith   ilu->info.dtcount = dtcount;
6785317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
6885317021SBarry Smith   PetscFunctionReturn(0);
6985317021SBarry Smith }
7085317021SBarry Smith EXTERN_C_END
7185317021SBarry Smith 
7285317021SBarry Smith EXTERN_C_BEGIN
7385317021SBarry Smith #undef __FUNCT__
7485317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor"
7585317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Factor(PC pc,PetscReal fill)
7685317021SBarry Smith {
7785317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
7885317021SBarry Smith 
7985317021SBarry Smith   PetscFunctionBegin;
8085317021SBarry Smith   dir->info.fill = fill;
8185317021SBarry Smith   PetscFunctionReturn(0);
8285317021SBarry Smith }
8385317021SBarry Smith EXTERN_C_END
8485317021SBarry Smith 
8585317021SBarry Smith EXTERN_C_BEGIN
8685317021SBarry Smith #undef __FUNCT__
8785317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
8885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
8985317021SBarry Smith {
9085317021SBarry Smith   PC_Factor      *dir = (PC_Factor*)pc->data;
9185317021SBarry Smith   PetscErrorCode ierr;
9285317021SBarry Smith   PetscTruth     flg;
9385317021SBarry Smith 
9485317021SBarry Smith   PetscFunctionBegin;
9585317021SBarry Smith   if (!pc->setupcalled) {
9685317021SBarry Smith      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
9785317021SBarry Smith      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
9885317021SBarry Smith   } else {
9985317021SBarry Smith     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
10085317021SBarry Smith     if (!flg) {
10185317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
10285317021SBarry Smith     }
10385317021SBarry Smith   }
10485317021SBarry Smith   PetscFunctionReturn(0);
10585317021SBarry Smith }
10685317021SBarry Smith EXTERN_C_END
10785317021SBarry Smith 
10885317021SBarry Smith EXTERN_C_BEGIN
10985317021SBarry Smith #undef __FUNCT__
11085317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor"
11185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_Factor(PC pc,PetscInt levels)
11285317021SBarry Smith {
11385317021SBarry Smith   PC_Factor      *ilu = (PC_Factor*)pc->data;
11485317021SBarry Smith 
11585317021SBarry Smith   PetscFunctionBegin;
11685317021SBarry Smith   if (!pc->setupcalled) {
11785317021SBarry Smith     ilu->info.levels = levels;
11885317021SBarry Smith   } else if (ilu->info.usedt || ilu->info.levels != levels) {
11985317021SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
12085317021SBarry Smith   }
12185317021SBarry Smith   PetscFunctionReturn(0);
12285317021SBarry Smith }
12385317021SBarry Smith EXTERN_C_END
12485317021SBarry Smith 
12585317021SBarry Smith EXTERN_C_BEGIN
12685317021SBarry Smith #undef __FUNCT__
12785317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
12885317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_Factor(PC pc)
12985317021SBarry Smith {
13085317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
13185317021SBarry Smith 
13285317021SBarry Smith   PetscFunctionBegin;
13385317021SBarry Smith   dir->info.diagonal_fill = 1;
13485317021SBarry Smith   PetscFunctionReturn(0);
13585317021SBarry Smith }
13685317021SBarry Smith EXTERN_C_END
13785317021SBarry Smith 
13885317021SBarry Smith 
13985317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
14085317021SBarry Smith 
14185317021SBarry Smith EXTERN_C_BEGIN
14285317021SBarry Smith #undef __FUNCT__
14385317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
14485317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot)
14585317021SBarry Smith {
14685317021SBarry Smith   PC_Factor *dir = (PC_Factor*)pc->data;
14785317021SBarry Smith 
14885317021SBarry Smith   PetscFunctionBegin;
14985317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
15085317021SBarry Smith   PetscFunctionReturn(0);
15185317021SBarry Smith }
15285317021SBarry Smith EXTERN_C_END
15385317021SBarry Smith 
15485317021SBarry Smith EXTERN_C_BEGIN
15585317021SBarry Smith #undef __FUNCT__
15685317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor"
15785317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatrix_Factor(PC pc,Mat *mat)
15885317021SBarry Smith {
15985317021SBarry Smith   PC_Factor *ilu = (PC_Factor*)pc->data;
16085317021SBarry Smith 
16185317021SBarry Smith   PetscFunctionBegin;
16285317021SBarry Smith   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
16385317021SBarry Smith   *mat = ilu->fact;
16485317021SBarry Smith   PetscFunctionReturn(0);
16585317021SBarry Smith }
16685317021SBarry Smith EXTERN_C_END
16785317021SBarry Smith 
16885317021SBarry Smith EXTERN_C_BEGIN
16985317021SBarry Smith #undef __FUNCT__
17085317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
17185317021SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
17285317021SBarry Smith {
17385317021SBarry Smith   PetscErrorCode ierr;
17485317021SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
17585317021SBarry Smith 
17685317021SBarry Smith   PetscFunctionBegin;
1777112b564SBarry Smith   if (lu->fact) {
17885317021SBarry Smith     const MatSolverPackage ltype;
17985317021SBarry Smith     PetscTruth             flg;
18085317021SBarry Smith     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
18185317021SBarry Smith     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
18285317021SBarry Smith     if (!flg) {
18385317021SBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
18485317021SBarry Smith     }
18585317021SBarry Smith   }
18685317021SBarry Smith   ierr = PetscStrfree(lu->solvertype);CHKERRQ(ierr);
18785317021SBarry Smith   ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
18885317021SBarry Smith   PetscFunctionReturn(0);
18985317021SBarry Smith }
19085317021SBarry Smith EXTERN_C_END
19185317021SBarry Smith 
19285317021SBarry Smith EXTERN_C_BEGIN
19385317021SBarry Smith #undef __FUNCT__
1947112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
1957112b564SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
1967112b564SBarry Smith {
1977112b564SBarry Smith   PC_Factor      *lu = (PC_Factor*)pc->data;
1987112b564SBarry Smith 
1997112b564SBarry Smith   PetscFunctionBegin;
2007112b564SBarry Smith   *stype = lu->solvertype;
2017112b564SBarry Smith   PetscFunctionReturn(0);
2027112b564SBarry Smith }
2037112b564SBarry Smith EXTERN_C_END
2047112b564SBarry Smith 
2057112b564SBarry Smith EXTERN_C_BEGIN
2067112b564SBarry Smith #undef __FUNCT__
2078ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot_Factor"
2088ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
20985317021SBarry Smith {
21085317021SBarry Smith   PC_Factor       *dir = (PC_Factor*)pc->data;
21185317021SBarry Smith 
21285317021SBarry Smith   PetscFunctionBegin;
21385317021SBarry Smith   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
21485317021SBarry Smith   dir->info.dtcol = dtcol;
21585317021SBarry Smith   PetscFunctionReturn(0);
21685317021SBarry Smith }
21785317021SBarry Smith EXTERN_C_END
2188ff23777SHong Zhang 
2198ff23777SHong Zhang EXTERN_C_BEGIN
2208ff23777SHong Zhang #undef __FUNCT__
2218ff23777SHong Zhang #define __FUNCT__ "PCSetFromOptions_Factor"
2228ff23777SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCSetFromOptions_Factor(PC pc)
2238ff23777SHong Zhang {
2248ff23777SHong Zhang   PC_Factor       *factor = (PC_Factor*)pc->data;
2258ff23777SHong Zhang   PetscErrorCode  ierr;
2268ff23777SHong Zhang   PetscTruth      flg = PETSC_FALSE,set;
2278ff23777SHong Zhang   char            tname[256], solvertype[64];
2288ff23777SHong Zhang   PetscFList      ordlist;
229*018dd85eSSatish Balay   PetscEnum       etmp;
2308ff23777SHong Zhang 
2318ff23777SHong Zhang   PetscFunctionBegin;
2328ff23777SHong Zhang   if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2338ff23777SHong Zhang   ierr = PetscOptionsTruth("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2348ff23777SHong Zhang   if (flg) {
2358ff23777SHong Zhang     ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
2368ff23777SHong Zhang   }
2378ff23777SHong 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);
2388ff23777SHong Zhang 
239d90ac83dSHong Zhang   ierr = PetscOptionsEnum("-pc_factor_shift_type","Shift added to diagonal","PCFactorSetShiftType",
240*018dd85eSSatish Balay                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr);
241*018dd85eSSatish Balay   if (flg) {
242*018dd85eSSatish Balay     ((PC_Factor*)factor)->info.shifttype = (PetscReal)(etmp);
243*018dd85eSSatish Balay   }
244d90ac83dSHong 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);
245d90ac83dSHong Zhang 
2468ff23777SHong 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);
2478ff23777SHong 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);
2488ff23777SHong Zhang 
2498ff23777SHong Zhang   flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
2508ff23777SHong Zhang   ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
2518ff23777SHong Zhang   if (set) {
2528ff23777SHong Zhang     ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
2538ff23777SHong Zhang   }
2548ff23777SHong Zhang 
2558ff23777SHong Zhang   flg  = PETSC_FALSE;
2568ff23777SHong Zhang   ierr = PetscOptionsTruth("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2578ff23777SHong Zhang   if (flg) {
2588ff23777SHong Zhang     ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
2598ff23777SHong Zhang   }
2608ff23777SHong Zhang   flg  = PETSC_FALSE;
2618ff23777SHong Zhang   ierr = PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2628ff23777SHong Zhang   if (flg) {
2638ff23777SHong Zhang     ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
2648ff23777SHong Zhang   }
2658ff23777SHong Zhang 
2668ff23777SHong Zhang   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
2678ff23777SHong Zhang   ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr);
2688ff23777SHong Zhang   if (flg) {
2698ff23777SHong Zhang     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
2708ff23777SHong Zhang   }
2718ff23777SHong Zhang 
2728ff23777SHong Zhang   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
2738ff23777SHong Zhang   ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
2748ff23777SHong Zhang   if (flg) {
2758ff23777SHong Zhang     ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
2768ff23777SHong Zhang   }
2778ff23777SHong Zhang   PetscFunctionReturn(0);
2788ff23777SHong Zhang }
2798ff23777SHong Zhang EXTERN_C_END
280914a5d51SHong Zhang 
281914a5d51SHong Zhang EXTERN_C_BEGIN
282914a5d51SHong Zhang #include "private/matimpl.h"
283914a5d51SHong Zhang #undef __FUNCT__
284914a5d51SHong Zhang #define __FUNCT__ "PCView_Factor"
285914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
286914a5d51SHong Zhang {
287914a5d51SHong Zhang   PC_Factor       *factor = (PC_Factor*)pc->data;
288914a5d51SHong Zhang   PetscErrorCode  ierr;
289914a5d51SHong Zhang   PetscTruth      isstring,iascii;
290914a5d51SHong Zhang 
291914a5d51SHong Zhang   PetscFunctionBegin;
292914a5d51SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
293914a5d51SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
294914a5d51SHong Zhang   if (iascii) {
295914a5d51SHong Zhang     if (factor->fact->factor == MAT_FACTOR_ILU || factor->fact->factor == MAT_FACTOR_ICC){
296914a5d51SHong Zhang       if (factor->info.dt > 0) {
297914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  drop tolerance %G\n",factor->info.dt);CHKERRQ(ierr);
298914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr);
299914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);CHKERRQ(ierr);
300914a5d51SHong Zhang       } else if (factor->info.levels == 1) {
301914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
302914a5d51SHong Zhang       } else {
303914a5d51SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
304914a5d51SHong Zhang       }
305914a5d51SHong Zhang     }
306914a5d51SHong Zhang 
307914a5d51SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %G\n",factor->info.zeropivot);CHKERRQ(ierr);
308f4db908eSBarry Smith     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
309914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  using Manteuffel shift\n");CHKERRQ(ierr);
310914a5d51SHong Zhang     }
311f4db908eSBarry Smith     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) {
312914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);
313914a5d51SHong Zhang     }
314f4db908eSBarry Smith     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) {
315914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);
316914a5d51SHong Zhang     }
317914a5d51SHong Zhang 
318914a5d51SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr);
319914a5d51SHong Zhang 
320914a5d51SHong Zhang     if (factor->fact) {
321914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %G, needed %G\n",factor->fact->info.fill_ratio_given,factor->fact->info.fill_ratio_needed);CHKERRQ(ierr);
322914a5d51SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");CHKERRQ(ierr);
323914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
324914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
325914a5d51SHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
326914a5d51SHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
327914a5d51SHong Zhang       ierr = MatView(factor->fact,viewer);CHKERRQ(ierr);
328914a5d51SHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
329914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
330914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
331914a5d51SHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
332914a5d51SHong Zhang     }
333914a5d51SHong Zhang 
334914a5d51SHong Zhang   } else if (isstring) {
335914a5d51SHong Zhang     if (factor->fact->factor == MAT_FACTOR_ILU || factor->fact->factor == MAT_FACTOR_ICC){
336914a5d51SHong Zhang       ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
337914a5d51SHong Zhang     }
338914a5d51SHong Zhang   } else {
339914a5d51SHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC_Factor",((PetscObject)viewer)->type_name);
340914a5d51SHong Zhang   }
341914a5d51SHong Zhang   PetscFunctionReturn(0);
342914a5d51SHong Zhang }
343914a5d51SHong Zhang EXTERN_C_END
344