xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
185317021SBarry Smith 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h"  I*/
385317021SBarry Smith 
4d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC pc)
5d71ae5a4SJacob Faibussowitsch {
6f8260c8fSBarry Smith   PC_Factor *icc = (PC_Factor *)pc->data;
7f8260c8fSBarry Smith 
8f8260c8fSBarry Smith   PetscFunctionBegin;
928b400f6SJacob 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()");
1048a46eb9SPierre Jolivet   if (!pc->setupcalled && !((PC_Factor *)icc)->fact) PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)icc)->solvertype, ((PC_Factor *)icc)->factortype, &((PC_Factor *)icc)->fact));
11*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12f8260c8fSBarry Smith }
13f8260c8fSBarry Smith 
14d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc, PetscReal z)
15d71ae5a4SJacob Faibussowitsch {
1685317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
1785317021SBarry Smith 
1885317021SBarry Smith   PetscFunctionBegin;
1985317021SBarry Smith   ilu->info.zeropivot = z;
20*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2185317021SBarry Smith }
2285317021SBarry Smith 
23d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType_Factor(PC pc, MatFactorShiftType shifttype)
24d71ae5a4SJacob Faibussowitsch {
25d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor *)pc->data;
26d90ac83dSHong Zhang 
27d90ac83dSHong Zhang   PetscFunctionBegin;
282fa5cd67SKarl Rupp   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
292fa5cd67SKarl Rupp   else {
30f4db908eSBarry Smith     dir->info.shifttype = (PetscReal)shifttype;
319371c9d4SSatish 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 */ }
32d90ac83dSHong Zhang   }
33*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34d90ac83dSHong Zhang }
35d90ac83dSHong Zhang 
36d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc, PetscReal shiftamount)
37d71ae5a4SJacob Faibussowitsch {
38d90ac83dSHong Zhang   PC_Factor *dir = (PC_Factor *)pc->data;
39d90ac83dSHong Zhang 
40d90ac83dSHong Zhang   PetscFunctionBegin;
412fa5cd67SKarl Rupp   if (shiftamount == (PetscReal)PETSC_DECIDE) dir->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON;
422fa5cd67SKarl Rupp   else dir->info.shiftamount = shiftamount;
43*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44d90ac83dSHong Zhang }
45d90ac83dSHong Zhang 
46d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc, PetscReal dt, PetscReal dtcol, PetscInt dtcount)
47d71ae5a4SJacob Faibussowitsch {
4885317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
4985317021SBarry Smith 
5085317021SBarry Smith   PetscFunctionBegin;
5185317021SBarry 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)) {
52ce94432eSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change tolerance after use");
5385317021SBarry Smith   }
5485317021SBarry Smith   ilu->info.usedt   = PETSC_TRUE;
5585317021SBarry Smith   ilu->info.dt      = dt;
5685317021SBarry Smith   ilu->info.dtcol   = dtcol;
5785317021SBarry Smith   ilu->info.dtcount = dtcount;
5885317021SBarry Smith   ilu->info.fill    = PETSC_DEFAULT;
59*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6085317021SBarry Smith }
6185317021SBarry Smith 
62d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill_Factor(PC pc, PetscReal fill)
63d71ae5a4SJacob Faibussowitsch {
6485317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
6585317021SBarry Smith 
6685317021SBarry Smith   PetscFunctionBegin;
6785317021SBarry Smith   dir->info.fill = fill;
68*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6985317021SBarry Smith }
7085317021SBarry Smith 
71d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc, MatOrderingType ordering)
72d71ae5a4SJacob Faibussowitsch {
7385317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
74ace3abfcSBarry Smith   PetscBool  flg;
7585317021SBarry Smith 
7685317021SBarry Smith   PetscFunctionBegin;
7785317021SBarry Smith   if (!pc->setupcalled) {
789566063dSJacob Faibussowitsch     PetscCall(PetscFree(dir->ordering));
799566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(ordering, (char **)&dir->ordering));
8085317021SBarry Smith   } else {
819566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(dir->ordering, ordering, &flg));
8228b400f6SJacob Faibussowitsch     PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change ordering after use");
8385317021SBarry Smith   }
84*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8585317021SBarry Smith }
8685317021SBarry Smith 
87d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels_Factor(PC pc, PetscInt *levels)
88d71ae5a4SJacob Faibussowitsch {
892591b318SToby Isaac   PC_Factor *ilu = (PC_Factor *)pc->data;
902591b318SToby Isaac 
912591b318SToby Isaac   PetscFunctionBegin;
922591b318SToby Isaac   *levels = ilu->info.levels;
93*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
942591b318SToby Isaac }
952591b318SToby Isaac 
96d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc, PetscReal *pivot)
97d71ae5a4SJacob Faibussowitsch {
98c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
99c7f610a1SBarry Smith 
100c7f610a1SBarry Smith   PetscFunctionBegin;
101c7f610a1SBarry Smith   *pivot = ilu->info.zeropivot;
102*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103c7f610a1SBarry Smith }
104c7f610a1SBarry Smith 
105d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc, PetscReal *shift)
106d71ae5a4SJacob Faibussowitsch {
107c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
108c7f610a1SBarry Smith 
109c7f610a1SBarry Smith   PetscFunctionBegin;
110c7f610a1SBarry Smith   *shift = ilu->info.shiftamount;
111*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112c7f610a1SBarry Smith }
113c7f610a1SBarry Smith 
114d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType_Factor(PC pc, MatFactorShiftType *type)
115d71ae5a4SJacob Faibussowitsch {
116c7f610a1SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
117c7f610a1SBarry Smith 
118c7f610a1SBarry Smith   PetscFunctionBegin;
119b729e602SBarry Smith   *type = (MatFactorShiftType)(int)ilu->info.shifttype;
120*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121c7f610a1SBarry Smith }
122c7f610a1SBarry Smith 
123d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels_Factor(PC pc, PetscInt levels)
124d71ae5a4SJacob Faibussowitsch {
12585317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
12685317021SBarry Smith 
12785317021SBarry Smith   PetscFunctionBegin;
1282fa5cd67SKarl Rupp   if (!pc->setupcalled) ilu->info.levels = levels;
1292fa5cd67SKarl Rupp   else if (ilu->info.levels != levels) {
130dbbe0bcdSBarry Smith     PetscUseTypeMethod(pc, reset); /* remove previous factored matrices */
1315b9c68c7SBarry Smith     pc->setupcalled  = 0;          /* force a complete rebuild of preconditioner factored matrices */
1325b9c68c7SBarry Smith     ilu->info.levels = levels;
1332472a847SBarry Smith   } else PetscCheck(!ilu->info.usedt, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Cannot change levels after use with ILUdt");
134*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13585317021SBarry Smith }
13685317021SBarry Smith 
137d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc, PetscBool flg)
138d71ae5a4SJacob Faibussowitsch {
13985317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
14085317021SBarry Smith 
14185317021SBarry Smith   PetscFunctionBegin;
14292e9c092SBarry Smith   dir->info.diagonal_fill = (PetscReal)flg;
143*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14492e9c092SBarry Smith }
14592e9c092SBarry Smith 
146d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc, PetscBool *flg)
147d71ae5a4SJacob Faibussowitsch {
14892e9c092SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
14992e9c092SBarry Smith 
15092e9c092SBarry Smith   PetscFunctionBegin;
15192e9c092SBarry Smith   *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE;
152*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15385317021SBarry Smith }
15485317021SBarry Smith 
155d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc, PetscBool pivot)
156d71ae5a4SJacob Faibussowitsch {
15785317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
15885317021SBarry Smith 
15985317021SBarry Smith   PetscFunctionBegin;
16085317021SBarry Smith   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
161*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16285317021SBarry Smith }
16385317021SBarry Smith 
164d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix_Factor(PC pc, Mat *mat)
165d71ae5a4SJacob Faibussowitsch {
16685317021SBarry Smith   PC_Factor *ilu = (PC_Factor *)pc->data;
16785317021SBarry Smith 
16885317021SBarry Smith   PetscFunctionBegin;
16928b400f6SJacob Faibussowitsch   PetscCheck(ilu->fact, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
17085317021SBarry Smith   *mat = ilu->fact;
171*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17285317021SBarry Smith }
17385317021SBarry Smith 
174978beb7fSPierre Jolivet /* allow access to preallocation information */
175978beb7fSPierre Jolivet #include <petsc/private/matimpl.h>
176978beb7fSPierre Jolivet 
177d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType_Factor(PC pc, MatSolverType stype)
178d71ae5a4SJacob Faibussowitsch {
17985317021SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
18085317021SBarry Smith 
18185317021SBarry Smith   PetscFunctionBegin;
182978beb7fSPierre Jolivet   if (lu->fact && lu->fact->assembled) {
183ea799195SBarry Smith     MatSolverType ltype;
184ace3abfcSBarry Smith     PetscBool     flg;
1859566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(lu->fact, &ltype));
1869566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(stype, ltype, &flg));
18728b400f6SJacob 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);
18800c67f3bSHong Zhang   }
18900c67f3bSHong Zhang 
1909566063dSJacob Faibussowitsch   PetscCall(PetscFree(lu->solvertype));
1919566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(stype, &lu->solvertype));
192*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19385317021SBarry Smith }
19485317021SBarry Smith 
195d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType_Factor(PC pc, MatSolverType *stype)
196d71ae5a4SJacob Faibussowitsch {
1977112b564SBarry Smith   PC_Factor *lu = (PC_Factor *)pc->data;
1987112b564SBarry Smith 
1997112b564SBarry Smith   PetscFunctionBegin;
2007112b564SBarry Smith   *stype = lu->solvertype;
201*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2027112b564SBarry Smith }
2037112b564SBarry Smith 
204d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc, PetscReal dtcol)
205d71ae5a4SJacob Faibussowitsch {
20685317021SBarry Smith   PC_Factor *dir = (PC_Factor *)pc->data;
20785317021SBarry Smith 
20885317021SBarry Smith   PetscFunctionBegin;
2092472a847SBarry 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);
21085317021SBarry Smith   dir->info.dtcol = dtcol;
211*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21285317021SBarry Smith }
2138ff23777SHong Zhang 
214d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_Factor(PC pc, PetscOptionItems *PetscOptionsObject)
215d71ae5a4SJacob Faibussowitsch {
2168ff23777SHong Zhang   PC_Factor        *factor = (PC_Factor *)pc->data;
2178afaa268SBarry Smith   PetscBool         flg, set;
2188ff23777SHong Zhang   char              tname[256], solvertype[64];
219140e18c1SBarry Smith   PetscFunctionList ordlist;
220018dd85eSSatish Balay   PetscEnum         etmp;
2218e37d05fSBarry Smith   PetscBool         inplace;
2228ff23777SHong Zhang 
2238ff23777SHong Zhang   PetscFunctionBegin;
2249566063dSJacob Faibussowitsch   PetscCall(PCFactorGetUseInPlace(pc, &inplace));
2259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_in_place", "Form factored matrix in the same memory as the matrix", "PCFactorSetUseInPlace", inplace, &flg, &set));
2261baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetUseInPlace(pc, flg));
2279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_fill", "Expected non-zeros in factored matrix", "PCFactorSetFill", ((PC_Factor *)factor)->info.fill, &((PC_Factor *)factor)->info.fill, NULL));
2288ff23777SHong Zhang 
2299566063dSJacob 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));
2301baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetShiftType(pc, (MatFactorShiftType)etmp));
2319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_shift_amount", "Shift added to diagonal", "PCFactorSetShiftAmount", ((PC_Factor *)factor)->info.shiftamount, &((PC_Factor *)factor)->info.shiftamount, NULL));
232d90ac83dSHong Zhang 
2339566063dSJacob 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));
2349566063dSJacob 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));
2358ff23777SHong Zhang 
2369566063dSJacob 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));
2371baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetPivotInBlocks(pc, flg));
2388ff23777SHong Zhang 
2399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_reuse_fill", "Use fill from previous factorization", "PCFactorSetReuseFill", PETSC_FALSE, &flg, &set));
2401baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetReuseFill(pc, flg));
2419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_factor_reuse_ordering", "Reuse ordering from previous factorization", "PCFactorSetReuseOrdering", PETSC_FALSE, &flg, &set));
2421baa6e33SBarry Smith   if (set) PetscCall(PCFactorSetReuseOrdering(pc, flg));
2438ff23777SHong Zhang 
2449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsDeprecated("-pc_factor_mat_solver_package", "-pc_factor_mat_solver_type", "3.9", NULL));
2459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-pc_factor_mat_solver_type", "Specific direct solver to use", "MatGetFactor", ((PC_Factor *)factor)->solvertype, solvertype, sizeof(solvertype), &flg));
2461baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetMatSolverType(pc, solvertype));
2479566063dSJacob Faibussowitsch   PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
2489566063dSJacob Faibussowitsch   PetscCall(MatGetOrderingList(&ordlist));
2499566063dSJacob 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));
2501baa6e33SBarry Smith   if (flg) PetscCall(PCFactorSetMatOrderingType(pc, tname));
251*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2528ff23777SHong Zhang }
253914a5d51SHong Zhang 
254d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_Factor(PC pc, PetscViewer viewer)
255d71ae5a4SJacob Faibussowitsch {
256914a5d51SHong Zhang   PC_Factor        *factor = (PC_Factor *)pc->data;
2574ac6704cSBarry Smith   PetscBool         isstring, iascii, canuseordering;
2584ac6704cSBarry Smith   MatInfo           info;
2599bd791bbSBarry Smith   MatOrderingType   ordering;
2601511cd71SPierre Jolivet   PetscViewerFormat format;
261914a5d51SHong Zhang 
262914a5d51SHong Zhang   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
265914a5d51SHong Zhang   if (iascii) {
2663d1c1ea0SBarry Smith     if (factor->inplace) {
2679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  in-place factorization\n"));
2683d1c1ea0SBarry Smith     } else {
2699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  out-of-place factorization\n"));
2703d1c1ea0SBarry Smith     }
2713d1c1ea0SBarry Smith 
2729566063dSJacob Faibussowitsch     if (factor->reusefill) PetscCall(PetscViewerASCIIPrintf(viewer, "  Reusing fill from past factorization\n"));
2739566063dSJacob Faibussowitsch     if (factor->reuseordering) PetscCall(PetscViewerASCIIPrintf(viewer, "  Reusing reordering from past factorization\n"));
274879e8a4dSBarry Smith     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
275914a5d51SHong Zhang       if (factor->info.dt > 0) {
2769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  drop tolerance %g\n", (double)factor->info.dt));
27763a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  max nonzeros per row %" PetscInt_FMT "\n", (PetscInt)factor->info.dtcount));
2789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  column permutation tolerance %g\n", (double)factor->info.dtcol));
279914a5d51SHong Zhang       } else if (factor->info.levels == 1) {
28063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " level of fill\n", (PetscInt)factor->info.levels));
281914a5d51SHong Zhang       } else {
28263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  %" PetscInt_FMT " levels of fill\n", (PetscInt)factor->info.levels));
283914a5d51SHong Zhang       }
284914a5d51SHong Zhang     }
285914a5d51SHong Zhang 
2869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerance for zero pivot %g\n", (double)factor->info.zeropivot));
2875e9742b9SJed Brown     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
2889566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  using %s [%s]\n", MatFactorShiftTypesDetail[(int)factor->info.shifttype], MatFactorShiftTypes[(int)factor->info.shifttype]));
289914a5d51SHong Zhang     }
290914a5d51SHong Zhang 
2914ac6704cSBarry Smith     if (factor->fact) {
2929566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(factor->fact, &canuseordering));
2934ac6704cSBarry Smith       if (!canuseordering) ordering = MATORDERINGEXTERNAL;
2944ac6704cSBarry Smith       else ordering = factor->ordering;
2959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix ordering: %s\n", ordering));
296978beb7fSPierre Jolivet       if (!factor->fact->assembled) {
2979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix solver type: %s\n", factor->fact->solvertype));
2989566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  matrix not yet factored; no additional information available\n"));
299978beb7fSPierre Jolivet       } else {
3009566063dSJacob Faibussowitsch         PetscCall(MatGetInfo(factor->fact, MAT_LOCAL, &info));
3019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  factor fill ratio given %g, needed %g\n", (double)info.fill_ratio_given, (double)info.fill_ratio_needed));
3029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    Factored matrix follows:\n"));
3039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3061511cd71SPierre Jolivet         PetscCall(PetscViewerGetFormat(viewer, &format));
3071511cd71SPierre Jolivet         PetscCall(PetscViewerPushFormat(viewer, format != PETSC_VIEWER_ASCII_INFO_DETAIL ? PETSC_VIEWER_ASCII_INFO : PETSC_VIEWER_ASCII_INFO_DETAIL));
3089566063dSJacob Faibussowitsch         PetscCall(MatView(factor->fact, viewer));
3099566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(viewer));
3109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
313914a5d51SHong Zhang       }
314978beb7fSPierre Jolivet     }
315914a5d51SHong Zhang 
316914a5d51SHong Zhang   } else if (isstring) {
3176335c336SBarry Smith     MatFactorType t;
3189566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(factor->fact, &t));
31948a46eb9SPierre Jolivet     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) PetscCall(PetscViewerStringSPrintf(viewer, " lvls=%" PetscInt_FMT ",order=%s", (PetscInt)factor->info.levels, factor->ordering));
320914a5d51SHong Zhang   }
321*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
322914a5d51SHong Zhang }
323