185317021SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 385317021SBarry Smith 485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 585317021SBarry Smith 63ca39a21SBarry Smith PetscErrorCode PCFactorSetUpMatSolverType_Factor(PC pc) 7f8260c8fSBarry Smith { 8f8260c8fSBarry Smith PC_Factor *icc = (PC_Factor*)pc->data; 9f8260c8fSBarry Smith 10f8260c8fSBarry Smith PetscFunctionBegin; 11*28b400f6SJacob 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()"); 12f8260c8fSBarry Smith if (!pc->setupcalled && !((PC_Factor*)icc)->fact) { 135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact)); 14f8260c8fSBarry Smith } 15f8260c8fSBarry Smith PetscFunctionReturn(0); 16f8260c8fSBarry Smith } 17f8260c8fSBarry Smith 187087cfbeSBarry Smith PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc,PetscReal z) 1985317021SBarry Smith { 2085317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 2185317021SBarry Smith 2285317021SBarry Smith PetscFunctionBegin; 2385317021SBarry Smith ilu->info.zeropivot = z; 2485317021SBarry Smith PetscFunctionReturn(0); 2585317021SBarry Smith } 2685317021SBarry Smith 277087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype) 28d90ac83dSHong Zhang { 29d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 30d90ac83dSHong Zhang 31d90ac83dSHong Zhang PetscFunctionBegin; 322fa5cd67SKarl Rupp if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 332fa5cd67SKarl Rupp else { 34f4db908eSBarry Smith dir->info.shifttype = (PetscReal) shifttype; 359dbfc187SHong Zhang if ((shifttype == MAT_SHIFT_NONZERO || shifttype == MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) { 369dbfc187SHong Zhang dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ 376cc3b3c1SHong Zhang } 38d90ac83dSHong Zhang } 39d90ac83dSHong Zhang PetscFunctionReturn(0); 40d90ac83dSHong Zhang } 41d90ac83dSHong Zhang 427087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount) 43d90ac83dSHong Zhang { 44d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 45d90ac83dSHong Zhang 46d90ac83dSHong Zhang PetscFunctionBegin; 472fa5cd67SKarl Rupp if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 482fa5cd67SKarl Rupp else dir->info.shiftamount = shiftamount; 49d90ac83dSHong Zhang PetscFunctionReturn(0); 50d90ac83dSHong Zhang } 51d90ac83dSHong Zhang 527087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 5385317021SBarry Smith { 5485317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 5585317021SBarry Smith 5685317021SBarry Smith PetscFunctionBegin; 5785317021SBarry 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)) { 58ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use"); 5985317021SBarry Smith } 6085317021SBarry Smith ilu->info.usedt = PETSC_TRUE; 6185317021SBarry Smith ilu->info.dt = dt; 6285317021SBarry Smith ilu->info.dtcol = dtcol; 6385317021SBarry Smith ilu->info.dtcount = dtcount; 6485317021SBarry Smith ilu->info.fill = PETSC_DEFAULT; 6585317021SBarry Smith PetscFunctionReturn(0); 6685317021SBarry Smith } 6785317021SBarry Smith 687087cfbeSBarry Smith PetscErrorCode PCFactorSetFill_Factor(PC pc,PetscReal fill) 6985317021SBarry Smith { 7085317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 7185317021SBarry Smith 7285317021SBarry Smith PetscFunctionBegin; 7385317021SBarry Smith dir->info.fill = fill; 7485317021SBarry Smith PetscFunctionReturn(0); 7585317021SBarry Smith } 7685317021SBarry Smith 7719fd82e9SBarry Smith PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering) 7885317021SBarry Smith { 7985317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 80ace3abfcSBarry Smith PetscBool flg; 8185317021SBarry Smith 8285317021SBarry Smith PetscFunctionBegin; 8385317021SBarry Smith if (!pc->setupcalled) { 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(dir->ordering)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(ordering,(char**)&dir->ordering)); 8685317021SBarry Smith } else { 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(dir->ordering,ordering,&flg)); 88*28b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use"); 8985317021SBarry Smith } 9085317021SBarry Smith PetscFunctionReturn(0); 9185317021SBarry Smith } 9285317021SBarry Smith 932591b318SToby Isaac PetscErrorCode PCFactorGetLevels_Factor(PC pc,PetscInt *levels) 942591b318SToby Isaac { 952591b318SToby Isaac PC_Factor *ilu = (PC_Factor*)pc->data; 962591b318SToby Isaac 972591b318SToby Isaac PetscFunctionBegin; 982591b318SToby Isaac *levels = ilu->info.levels; 992591b318SToby Isaac PetscFunctionReturn(0); 1002591b318SToby Isaac } 1012591b318SToby Isaac 102c7f610a1SBarry Smith PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc,PetscReal *pivot) 103c7f610a1SBarry Smith { 104c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 105c7f610a1SBarry Smith 106c7f610a1SBarry Smith PetscFunctionBegin; 107c7f610a1SBarry Smith *pivot = ilu->info.zeropivot; 108c7f610a1SBarry Smith PetscFunctionReturn(0); 109c7f610a1SBarry Smith } 110c7f610a1SBarry Smith 111c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc,PetscReal *shift) 112c7f610a1SBarry Smith { 113c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 114c7f610a1SBarry Smith 115c7f610a1SBarry Smith PetscFunctionBegin; 116c7f610a1SBarry Smith *shift = ilu->info.shiftamount; 117c7f610a1SBarry Smith PetscFunctionReturn(0); 118c7f610a1SBarry Smith } 119c7f610a1SBarry Smith 120c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftType_Factor(PC pc,MatFactorShiftType *type) 121c7f610a1SBarry Smith { 122c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 123c7f610a1SBarry Smith 124c7f610a1SBarry Smith PetscFunctionBegin; 125b729e602SBarry Smith *type = (MatFactorShiftType) (int) ilu->info.shifttype; 126c7f610a1SBarry Smith PetscFunctionReturn(0); 127c7f610a1SBarry Smith } 128c7f610a1SBarry Smith 1297087cfbeSBarry Smith PetscErrorCode PCFactorSetLevels_Factor(PC pc,PetscInt levels) 13085317021SBarry Smith { 13185317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 13285317021SBarry Smith 13385317021SBarry Smith PetscFunctionBegin; 1342fa5cd67SKarl Rupp if (!pc->setupcalled) ilu->info.levels = levels; 1352fa5cd67SKarl Rupp else if (ilu->info.levels != levels) { 1365f80ce2aSJacob Faibussowitsch CHKERRQ((*pc->ops->reset)(pc)); /* remove previous factored matrices */ 1375b9c68c7SBarry Smith pc->setupcalled = 0; /* force a complete rebuild of preconditioner factored matrices */ 1385b9c68c7SBarry Smith ilu->info.levels = levels; 1392c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(ilu->info.usedt,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt"); 14085317021SBarry Smith PetscFunctionReturn(0); 14185317021SBarry Smith } 14285317021SBarry Smith 14392e9c092SBarry Smith PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg) 14485317021SBarry Smith { 14585317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 14685317021SBarry Smith 14785317021SBarry Smith PetscFunctionBegin; 14892e9c092SBarry Smith dir->info.diagonal_fill = (PetscReal) flg; 14992e9c092SBarry Smith PetscFunctionReturn(0); 15092e9c092SBarry Smith } 15192e9c092SBarry Smith 15292e9c092SBarry Smith PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg) 15392e9c092SBarry Smith { 15492e9c092SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 15592e9c092SBarry Smith 15692e9c092SBarry Smith PetscFunctionBegin; 15792e9c092SBarry Smith *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE; 15885317021SBarry Smith PetscFunctionReturn(0); 15985317021SBarry Smith } 16085317021SBarry Smith 16185317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 16285317021SBarry Smith 1637087cfbeSBarry Smith PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot) 16485317021SBarry Smith { 16585317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 16685317021SBarry Smith 16785317021SBarry Smith PetscFunctionBegin; 16885317021SBarry Smith dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 16985317021SBarry Smith PetscFunctionReturn(0); 17085317021SBarry Smith } 17185317021SBarry Smith 1727087cfbeSBarry Smith PetscErrorCode PCFactorGetMatrix_Factor(PC pc,Mat *mat) 17385317021SBarry Smith { 17485317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 17585317021SBarry Smith 17685317021SBarry Smith PetscFunctionBegin; 177*28b400f6SJacob Faibussowitsch PetscCheck(ilu->fact,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 17885317021SBarry Smith *mat = ilu->fact; 17985317021SBarry Smith PetscFunctionReturn(0); 18085317021SBarry Smith } 18185317021SBarry Smith 182978beb7fSPierre Jolivet /* allow access to preallocation information */ 183978beb7fSPierre Jolivet #include <petsc/private/matimpl.h> 184978beb7fSPierre Jolivet 185ea799195SBarry Smith PetscErrorCode PCFactorSetMatSolverType_Factor(PC pc,MatSolverType stype) 18685317021SBarry Smith { 18785317021SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 18885317021SBarry Smith 18985317021SBarry Smith PetscFunctionBegin; 190978beb7fSPierre Jolivet if (lu->fact && lu->fact->assembled) { 191ea799195SBarry Smith MatSolverType ltype; 192ace3abfcSBarry Smith PetscBool flg; 1935f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorGetSolverType(lu->fact,<ype)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcmp(stype,ltype,&flg)); 195*28b400f6SJacob 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); 19600c67f3bSHong Zhang } 19700c67f3bSHong Zhang 1985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(lu->solvertype)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrallocpy(stype,&lu->solvertype)); 20085317021SBarry Smith PetscFunctionReturn(0); 20185317021SBarry Smith } 20285317021SBarry Smith 203ea799195SBarry Smith PetscErrorCode PCFactorGetMatSolverType_Factor(PC pc,MatSolverType *stype) 2047112b564SBarry Smith { 2057112b564SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 2067112b564SBarry Smith 2077112b564SBarry Smith PetscFunctionBegin; 2087112b564SBarry Smith *stype = lu->solvertype; 2097112b564SBarry Smith PetscFunctionReturn(0); 2107112b564SBarry Smith } 2117112b564SBarry Smith 2127087cfbeSBarry Smith PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol) 21385317021SBarry Smith { 21485317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 21585317021SBarry Smith 21685317021SBarry Smith PetscFunctionBegin; 2172c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 21885317021SBarry Smith dir->info.dtcol = dtcol; 21985317021SBarry Smith PetscFunctionReturn(0); 22085317021SBarry Smith } 2218ff23777SHong Zhang 2224416b707SBarry Smith PetscErrorCode PCSetFromOptions_Factor(PetscOptionItems *PetscOptionsObject,PC pc) 2238ff23777SHong Zhang { 2248ff23777SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 2258afaa268SBarry Smith PetscBool flg,set; 2268ff23777SHong Zhang char tname[256], solvertype[64]; 227140e18c1SBarry Smith PetscFunctionList ordlist; 228018dd85eSSatish Balay PetscEnum etmp; 2298e37d05fSBarry Smith PetscBool inplace; 2308ff23777SHong Zhang 2318ff23777SHong Zhang PetscFunctionBegin; 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorGetUseInPlace(pc,&inplace)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set)); 2348e37d05fSBarry Smith if (set) { 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetUseInPlace(pc,flg)); 2368ff23777SHong Zhang } 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,NULL)); 2388ff23777SHong Zhang 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg)); 240018dd85eSSatish Balay if (flg) { 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetShiftType(pc,(MatFactorShiftType)etmp)); 242018dd85eSSatish Balay } 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,NULL)); 244d90ac83dSHong Zhang 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,NULL)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2478ff23777SHong Zhang 2485f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2498ff23777SHong Zhang if (set) { 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetPivotInBlocks(pc,flg)); 2518ff23777SHong Zhang } 2528ff23777SHong Zhang 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set)); 2548afaa268SBarry Smith if (set) { 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetReuseFill(pc,flg)); 2568ff23777SHong Zhang } 2575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set)); 2588afaa268SBarry Smith if (set) { 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetReuseOrdering(pc,flg)); 2608ff23777SHong Zhang } 2618ff23777SHong Zhang 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsDeprecated("-pc_factor_mat_solver_package","-pc_factor_mat_solver_type","3.9",NULL)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-pc_factor_mat_solver_type","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,sizeof(solvertype),&flg)); 2648ff23777SHong Zhang if (flg) { 2655f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetMatSolverType(pc,solvertype)); 2668ff23777SHong Zhang } 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetDefaultOrdering_Factor(pc)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrderingList(&ordlist)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,sizeof(tname),&flg)); 2704ac6704cSBarry Smith if (flg) { 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorSetMatOrderingType(pc,tname)); 2724ac6704cSBarry Smith } 2738ff23777SHong Zhang PetscFunctionReturn(0); 2748ff23777SHong Zhang } 275914a5d51SHong Zhang 276914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer) 277914a5d51SHong Zhang { 278914a5d51SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 2794ac6704cSBarry Smith PetscBool isstring,iascii,canuseordering; 2804ac6704cSBarry Smith MatInfo info; 2819bd791bbSBarry Smith MatOrderingType ordering; 282914a5d51SHong Zhang 283914a5d51SHong Zhang PetscFunctionBegin; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 286914a5d51SHong Zhang if (iascii) { 2873d1c1ea0SBarry Smith if (factor->inplace) { 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," in-place factorization\n")); 2893d1c1ea0SBarry Smith } else { 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," out-of-place factorization\n")); 2913d1c1ea0SBarry Smith } 2923d1c1ea0SBarry Smith 2935f80ce2aSJacob Faibussowitsch if (factor->reusefill) CHKERRQ(PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n")); 2945f80ce2aSJacob Faibussowitsch if (factor->reuseordering) CHKERRQ(PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n")); 295879e8a4dSBarry Smith if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) { 296914a5d51SHong Zhang if (factor->info.dt > 0) { 2975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," drop tolerance %g\n",(double)factor->info.dt)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," max nonzeros per row %D\n",factor->info.dtcount)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," column permutation tolerance %g\n",(double)factor->info.dtcol)); 300914a5d51SHong Zhang } else if (factor->info.levels == 1) { 3015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," %D level of fill\n",(PetscInt)factor->info.levels)); 302914a5d51SHong Zhang } else { 3035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," %D levels of fill\n",(PetscInt)factor->info.levels)); 304914a5d51SHong Zhang } 305914a5d51SHong Zhang } 306914a5d51SHong Zhang 3075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," tolerance for zero pivot %g\n",(double)factor->info.zeropivot)); 3085e9742b9SJed Brown if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */ 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype])); 310914a5d51SHong Zhang } 311914a5d51SHong Zhang 3124ac6704cSBarry Smith if (factor->fact) { 3135f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorGetCanUseOrdering(factor->fact,&canuseordering)); 3144ac6704cSBarry Smith if (!canuseordering) ordering = MATORDERINGEXTERNAL; 3154ac6704cSBarry Smith else ordering = factor->ordering; 3165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ordering)); 317978beb7fSPierre Jolivet if (!factor->fact->assembled) { 3185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," matrix solver type: %s\n",factor->fact->solvertype)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," matrix not yet factored; no additional information available\n")); 320978beb7fSPierre Jolivet } else { 3215f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(factor->fact,MAT_LOCAL,&info)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," Factored matrix follows:\n")); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(factor->fact,viewer)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 333914a5d51SHong Zhang } 334978beb7fSPierre Jolivet } 335914a5d51SHong Zhang 336914a5d51SHong Zhang } else if (isstring) { 3376335c336SBarry Smith MatFactorType t; 3385f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactorType(factor->fact,&t)); 3396335c336SBarry Smith if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) { 3405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering)); 341914a5d51SHong Zhang } 342914a5d51SHong Zhang } 343914a5d51SHong Zhang PetscFunctionReturn(0); 344914a5d51SHong Zhang } 345