xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 1511cd715a1f0c8d257549c5ebe5cee9c6feed4d)
185317021SBarry Smith 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h"  I*/
385317021SBarry Smith 
485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
585317021SBarry Smith 
69371c9d4SSatish Balay PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC pc) {
7f8260c8fSBarry Smith   PC_Factor *icc = (PC_Factor *)pc->data;
8f8260c8fSBarry Smith 
9f8260c8fSBarry Smith   PetscFunctionBegin;
1028b400f6SJacob Faibussowitsch   PetscCheck(pc->pmat, 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()");
1148a46eb9SPierre Jolivet   if (!pc->setupcalled && !((PC_Factor *)icc)->fact) PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)icc)->solvertype, ((PC_Factor *)icc)->factortype, &((PC_Factor *)icc)->fact));
12f8260c8fSBarry Smith   PetscFunctionReturn(0);
13f8260c8fSBarry Smith }
14f8260c8fSBarry Smith 
159371c9d4SSatish Balay PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc, PetscReal z) {
1685317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
1785317021SBarry Smith 
1885317021SBarry Smith   PetscFunctionBegin;
1985317021SBarry Smith   ilu->info.zeropivot = z;
2085317021SBarry Smith   PetscFunctionReturn(0);
2185317021SBarry Smith }
2285317021SBarry Smith 
239371c9d4SSatish Balay PetscErrorCode PCFactorSetShiftType_Factor(PC pc, MatFactorShiftType shifttype) {
24d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor *)pc->data;
25d90ac83dSHong Zhang 
26d90ac83dSHong Zhang   PetscFunctionBegin;
272fa5cd67SKarl Rupp   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
282fa5cd67SKarl Rupp   else {
29f4db908eSBarry Smith     dir->info.shifttype = (PetscReal)shifttype;
309371c9d4SSatish Balay     if ((shifttype == MAT_SHIFT_NONZERO || shifttype == MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) { dir->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ }
31d90ac83dSHong Zhang   }
32d90ac83dSHong Zhang   PetscFunctionReturn(0);
33d90ac83dSHong Zhang }
34d90ac83dSHong Zhang 
359371c9d4SSatish Balay PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc, PetscReal shiftamount) {
36d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor *)pc->data;
37d90ac83dSHong Zhang 
38d90ac83dSHong Zhang   PetscFunctionBegin;
392fa5cd67SKarl Rupp   if (shiftamount == (PetscReal)PETSC_DECIDE) dir->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON;
402fa5cd67SKarl Rupp   else dir->info.shiftamount = shiftamount;
41d90ac83dSHong Zhang   PetscFunctionReturn(0);
42d90ac83dSHong Zhang }
43d90ac83dSHong Zhang 
449371c9d4SSatish Balay PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount) {
4585317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
4685317021SBarry Smith 
4785317021SBarry Smith   PetscFunctionBegin;
4885317021SBarry 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)) {
49ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change tolerance after use");
5085317021SBarry Smith   }
5185317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
5285317021SBarry Smith   ilu->info.dt      = dt;
5385317021SBarry Smith   ilu->info.dtcol   = dtcol;
5485317021SBarry Smith   ilu->info.dtcount = dtcount;
5585317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
5685317021SBarry Smith   PetscFunctionReturn(0);
5785317021SBarry Smith }
5885317021SBarry Smith 
599371c9d4SSatish Balay PetscErrorCode PCFactorSetFill_Factor(PC pc, PetscReal fill) {
6085317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
6185317021SBarry Smith 
6285317021SBarry Smith   PetscFunctionBegin;
6385317021SBarry Smith   dir->info.fill = fill;
6485317021SBarry Smith   PetscFunctionReturn(0);
6585317021SBarry Smith }
6685317021SBarry Smith 
679371c9d4SSatish Balay PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc, MatOrderingType ordering) {
6885317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
69ace3abfcSBarry Smith   PetscBool  flg;
7085317021SBarry Smith 
7185317021SBarry Smith   PetscFunctionBegin;
7285317021SBarry Smith   if (!pc->setupcalled) {
739566063dSJacob Faibussowitsch     PetscCall(PetscFree(dir->ordering));
749566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(ordering, (char **)&dir->ordering));
7585317021SBarry Smith   } else {
769566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(dir->ordering, ordering, &flg));
7728b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change ordering after use");
7885317021SBarry Smith   }
7985317021SBarry Smith   PetscFunctionReturn(0);
8085317021SBarry Smith }
8185317021SBarry Smith 
829371c9d4SSatish Balay PetscErrorCode PCFactorGetLevels_Factor(PC pc, PetscInt *levels) {
832591b318SToby Isaac   PC_Factor *ilu = (PC_Factor *)pc->data;
842591b318SToby Isaac 
852591b318SToby Isaac   PetscFunctionBegin;
862591b318SToby Isaac   *levels = ilu->info.levels;
872591b318SToby Isaac   PetscFunctionReturn(0);
882591b318SToby Isaac }
892591b318SToby Isaac 
909371c9d4SSatish Balay PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc, PetscReal *pivot) {
91c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
92c7f610a1SBarry Smith 
93c7f610a1SBarry Smith   PetscFunctionBegin;
94c7f610a1SBarry Smith   *pivot = ilu->info.zeropivot;
95c7f610a1SBarry Smith   PetscFunctionReturn(0);
96c7f610a1SBarry Smith }
97c7f610a1SBarry Smith 
989371c9d4SSatish Balay PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc, PetscReal *shift) {
99c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
100c7f610a1SBarry Smith 
101c7f610a1SBarry Smith   PetscFunctionBegin;
102c7f610a1SBarry Smith   *shift = ilu->info.shiftamount;
103c7f610a1SBarry Smith   PetscFunctionReturn(0);
104c7f610a1SBarry Smith }
105c7f610a1SBarry Smith 
1069371c9d4SSatish Balay PetscErrorCode PCFactorGetShiftType_Factor(PC pc, MatFactorShiftType *type) {
107c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
108c7f610a1SBarry Smith 
109c7f610a1SBarry Smith   PetscFunctionBegin;
110b729e602SBarry Smith   *type = (MatFactorShiftType)(int)ilu->info.shifttype;
111c7f610a1SBarry Smith   PetscFunctionReturn(0);
112c7f610a1SBarry Smith }
113c7f610a1SBarry Smith 
1149371c9d4SSatish Balay PetscErrorCode PCFactorSetLevels_Factor(PC pc, PetscInt levels) {
11585317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
11685317021SBarry Smith 
11785317021SBarry Smith   PetscFunctionBegin;
1182fa5cd67SKarl Rupp   if (!pc->setupcalled) ilu->info.levels = levels;
1192fa5cd67SKarl Rupp   else if (ilu->info.levels != levels) {
120dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, reset); /* remove previous factored matrices */
1215b9c68c7SBarry Smith     pc->setupcalled  = 0;          /* force a complete rebuild of preconditioner factored matrices */
1225b9c68c7SBarry Smith     ilu->info.levels = levels;
1232472a847SBarry Smith   } else PetscCheck(!ilu->info.usedt, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change levels after use with ILUdt");
12485317021SBarry Smith   PetscFunctionReturn(0);
12585317021SBarry Smith }
12685317021SBarry Smith 
1279371c9d4SSatish Balay PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc, PetscBool flg) {
12885317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
12985317021SBarry Smith 
13085317021SBarry Smith   PetscFunctionBegin;
13192e9c092SBarry Smith   dir->info.diagonal_fill = (PetscReal)flg;
13292e9c092SBarry Smith   PetscFunctionReturn(0);
13392e9c092SBarry Smith }
13492e9c092SBarry Smith 
1359371c9d4SSatish Balay PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc, PetscBool *flg) {
13692e9c092SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
13792e9c092SBarry Smith 
13892e9c092SBarry Smith   PetscFunctionBegin;
13992e9c092SBarry Smith   *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE;
14085317021SBarry Smith   PetscFunctionReturn(0);
14185317021SBarry Smith }
14285317021SBarry Smith 
14385317021SBarry Smith /* ------------------------------------------------------------------------------------------*/
14485317021SBarry Smith 
1459371c9d4SSatish Balay PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc, PetscBool pivot) {
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 
1539371c9d4SSatish Balay PetscErrorCode PCFactorGetMatrix_Factor(PC pc, Mat *mat) {
15485317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
15585317021SBarry Smith 
15685317021SBarry Smith   PetscFunctionBegin;
15728b400f6SJacob Faibussowitsch   PetscCheck(ilu->fact, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
15885317021SBarry Smith   *mat = ilu->fact;
15985317021SBarry Smith   PetscFunctionReturn(0);
16085317021SBarry Smith }
16185317021SBarry Smith 
162978beb7fSPierre Jolivet /* allow access to preallocation information */
163978beb7fSPierre Jolivet #include <petsc/private/matimpl.h>
164978beb7fSPierre Jolivet 
1659371c9d4SSatish Balay PetscErrorCode PCFactorSetMatSolverType_Factor(PC pc, MatSolverType stype) {
16685317021SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
16785317021SBarry Smith 
16885317021SBarry Smith   PetscFunctionBegin;
169978beb7fSPierre Jolivet   if (lu->fact && lu->fact->assembled) {
170ea799195SBarry Smith     MatSolverType ltype;
171ace3abfcSBarry Smith     PetscBool     flg;
1729566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(lu->fact, &ltype));
1739566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(stype, ltype, &flg));
17428b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change solver matrix package from %s to %s after PC has been setup or used", ltype, stype);
17500c67f3bSHong Zhang   }
17600c67f3bSHong Zhang 
1779566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->solvertype));
1789566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(stype, &lu->solvertype));
17985317021SBarry Smith   PetscFunctionReturn(0);
18085317021SBarry Smith }
18185317021SBarry Smith 
1829371c9d4SSatish Balay PetscErrorCode PCFactorGetMatSolverType_Factor(PC pc, MatSolverType *stype) {
1837112b564SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
1847112b564SBarry Smith 
1857112b564SBarry Smith   PetscFunctionBegin;
1867112b564SBarry Smith   *stype = lu->solvertype;
1877112b564SBarry Smith   PetscFunctionReturn(0);
1887112b564SBarry Smith }
1897112b564SBarry Smith 
1909371c9d4SSatish Balay PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc, PetscReal dtcol) {
19185317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
19285317021SBarry Smith 
19385317021SBarry Smith   PetscFunctionBegin;
1942472a847SBarry Smith   PetscCheck(dtcol >= 0.0 && dtcol <= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Column pivot tolerance is %g must be between 0 and 1", (double)dtcol);
19585317021SBarry Smith   dir->info.dtcol = dtcol;
19685317021SBarry Smith   PetscFunctionReturn(0);
19785317021SBarry Smith }
1988ff23777SHong Zhang 
1999371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_Factor(PC pc, PetscOptionItems *PetscOptionsObject) {
2008ff23777SHong Zhang   PC_Factor        *factor = (PC_Factor *)pc->data;
2018afaa268SBarry Smith   PetscBool         flg, set;
2028ff23777SHong Zhang   char              tname[256], solvertype[64];
203140e18c1SBarry Smith   PetscFunctionList ordlist;
204018dd85eSSatish Balay   PetscEnum         etmp;
2058e37d05fSBarry Smith   PetscBool         inplace;
2068ff23777SHong Zhang 
2078ff23777SHong Zhang   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(PCFactorGetUseInPlace(pc, &inplace));
2099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_in_place", "Form factored matrix in the same memory as the matrix", "PCFactorSetUseInPlace", inplace, &flg, &set));
2101baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetUseInPlace(pc, flg));
2119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_fill", "Expected non-zeros in factored matrix", "PCFactorSetFill", ((PC_Factor *)factor)->info.fill, &((PC_Factor *)factor)->info.fill, NULL));
2128ff23777SHong Zhang 
2139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_factor_shift_type", "Type of shift to add to diagonal", "PCFactorSetShiftType", MatFactorShiftTypes, (PetscEnum)(int)((PC_Factor *)factor)->info.shifttype, &etmp, &flg));
2141baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetShiftType(pc, (MatFactorShiftType)etmp));
2159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_shift_amount", "Shift added to diagonal", "PCFactorSetShiftAmount", ((PC_Factor *)factor)->info.shiftamount, &((PC_Factor *)factor)->info.shiftamount, NULL));
216d90ac83dSHong Zhang 
2179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_zeropivot", "Pivot is considered zero if less than", "PCFactorSetZeroPivot", ((PC_Factor *)factor)->info.zeropivot, &((PC_Factor *)factor)->info.zeropivot, NULL));
2189566063dSJacob Faibussowitsch   PetscCall(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));
2198ff23777SHong Zhang 
2209566063dSJacob Faibussowitsch   PetscCall(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));
2211baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetPivotInBlocks(pc, flg));
2228ff23777SHong Zhang 
2239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_reuse_fill", "Use fill from previous factorization", "PCFactorSetReuseFill", PETSC_FALSE, &flg, &set));
2241baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetReuseFill(pc, flg));
2259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_reuse_ordering", "Reuse ordering from previous factorization", "PCFactorSetReuseOrdering", PETSC_FALSE, &flg, &set));
2261baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetReuseOrdering(pc, flg));
2278ff23777SHong Zhang 
2289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-pc_factor_mat_solver_package", "-pc_factor_mat_solver_type", "3.9", NULL));
2299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-pc_factor_mat_solver_type", "Specific direct solver to use", "MatGetFactor", ((PC_Factor *)factor)->solvertype, solvertype, sizeof(solvertype), &flg));
2301baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetMatSolverType(pc, solvertype));
2319566063dSJacob Faibussowitsch   PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
2329566063dSJacob Faibussowitsch   PetscCall(MatGetOrderingList(&ordlist));
2339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_factor_mat_ordering_type", "Reordering to reduce nonzeros in factored matrix", "PCFactorSetMatOrderingType", ordlist, ((PC_Factor *)factor)->ordering, tname, sizeof(tname), &flg));
2341baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetMatOrderingType(pc, tname));
2358ff23777SHong Zhang   PetscFunctionReturn(0);
2368ff23777SHong Zhang }
237914a5d51SHong Zhang 
2389371c9d4SSatish Balay PetscErrorCode PCView_Factor(PC pc, PetscViewer viewer) {
239914a5d51SHong Zhang   PC_Factor        *factor = (PC_Factor *)pc->data;
2404ac6704cSBarry Smith   PetscBool         isstring, iascii, canuseordering;
2414ac6704cSBarry Smith   MatInfo           info;
2429bd791bbSBarry Smith   MatOrderingType   ordering;
243*1511cd71SPierre Jolivet   PetscViewerFormat format;
244914a5d51SHong Zhang 
245914a5d51SHong Zhang   PetscFunctionBegin;
2469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
248914a5d51SHong Zhang   if (iascii) {
2493d1c1ea0SBarry Smith     if (factor->inplace) {
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  in-place factorization\n"));
2513d1c1ea0SBarry Smith     } else {
2529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  out-of-place factorization\n"));
2533d1c1ea0SBarry Smith     }
2543d1c1ea0SBarry Smith 
2559566063dSJacob Faibussowitsch     if (factor->reusefill) PetscCall(PetscViewerASCIIPrintf(viewer, "  Reusing fill from past factorization\n"));
2569566063dSJacob Faibussowitsch     if (factor->reuseordering) PetscCall(PetscViewerASCIIPrintf(viewer, "  Reusing reordering from past factorization\n"));
257879e8a4dSBarry Smith     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
258914a5d51SHong Zhang       if (factor->info.dt > 0) {
2599566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  drop tolerance %g\n", (double)factor->info.dt));
26063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  max nonzeros per row %" PetscInt_FMT "\n", (PetscInt)factor->info.dtcount));
2619566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  column permutation tolerance %g\n", (double)factor->info.dtcol));
262914a5d51SHong Zhang       } else if (factor->info.levels == 1) {
26363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " level of fill\n", (PetscInt)factor->info.levels));
264914a5d51SHong Zhang       } else {
26563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " levels of fill\n", (PetscInt)factor->info.levels));
266914a5d51SHong Zhang       }
267914a5d51SHong Zhang     }
268914a5d51SHong Zhang 
2699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerance for zero pivot %g\n", (double)factor->info.zeropivot));
2705e9742b9SJed Brown     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
2719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  using %s [%s]\n", MatFactorShiftTypesDetail[(int)factor->info.shifttype], MatFactorShiftTypes[(int)factor->info.shifttype]));
272914a5d51SHong Zhang     }
273914a5d51SHong Zhang 
2744ac6704cSBarry Smith     if (factor->fact) {
2759566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(factor->fact, &canuseordering));
2764ac6704cSBarry Smith       if (!canuseordering) ordering = MATORDERINGEXTERNAL;
2774ac6704cSBarry Smith       else ordering = factor->ordering;
2789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix ordering: %s\n", ordering));
279978beb7fSPierre Jolivet       if (!factor->fact->assembled) {
2809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix solver type: %s\n", factor->fact->solvertype));
2819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix not yet factored; no additional information available\n"));
282978beb7fSPierre Jolivet       } else {
2839566063dSJacob Faibussowitsch         PetscCall(MatGetInfo(factor->fact, MAT_LOCAL, &info));
2849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  factor fill ratio given %g, needed %g\n", (double)info.fill_ratio_given, (double)info.fill_ratio_needed));
2859566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Factored matrix follows:\n"));
2869566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2879566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2889566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
289*1511cd71SPierre Jolivet         PetscCall(PetscViewerGetFormat(viewer, &format));
290*1511cd71SPierre Jolivet         PetscCall(PetscViewerPushFormat(viewer, format != PETSC_VIEWER_ASCII_INFO_DETAIL ? PETSC_VIEWER_ASCII_INFO : PETSC_VIEWER_ASCII_INFO_DETAIL));
2919566063dSJacob Faibussowitsch         PetscCall(MatView(factor->fact, viewer));
2929566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
2939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
2949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
2959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
296914a5d51SHong Zhang       }
297978beb7fSPierre Jolivet     }
298914a5d51SHong Zhang 
299914a5d51SHong Zhang   } else if (isstring) {
3006335c336SBarry Smith     MatFactorType t;
3019566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(factor->fact, &t));
30248a46eb9SPierre Jolivet     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) PetscCall(PetscViewerStringSPrintf(viewer, " lvls=%" PetscInt_FMT ",order=%s", (PetscInt)factor->info.levels, factor->ordering));
303914a5d51SHong Zhang   }
304914a5d51SHong Zhang   PetscFunctionReturn(0);
305914a5d51SHong Zhang }
306