185317021SBarry Smith 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 385317021SBarry Smith 485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 585317021SBarry Smith 6f8260c8fSBarry Smith 7f8260c8fSBarry Smith EXTERN_C_BEGIN 8f8260c8fSBarry Smith #undef __FUNCT__ 9f8260c8fSBarry Smith #define __FUNCT__ "PCFactorSetUpMatSolverPackage_Factor" 10f8260c8fSBarry Smith PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc) 11f8260c8fSBarry Smith { 12f8260c8fSBarry Smith PC_Factor *icc = (PC_Factor*)pc->data; 13f8260c8fSBarry Smith PetscErrorCode ierr; 14f8260c8fSBarry Smith 15f8260c8fSBarry Smith PetscFunctionBegin; 16f8260c8fSBarry Smith if (!pc->setupcalled && !((PC_Factor*)icc)->fact){ 17f8260c8fSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr); 18f8260c8fSBarry Smith } 19f8260c8fSBarry Smith PetscFunctionReturn(0); 20f8260c8fSBarry Smith } 21f8260c8fSBarry Smith EXTERN_C_END 22f8260c8fSBarry Smith 2385317021SBarry Smith EXTERN_C_BEGIN 2485317021SBarry Smith #undef __FUNCT__ 2585317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor" 267087cfbeSBarry Smith PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc,PetscReal z) 2785317021SBarry Smith { 2885317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 2985317021SBarry Smith 3085317021SBarry Smith PetscFunctionBegin; 3185317021SBarry Smith ilu->info.zeropivot = z; 3285317021SBarry Smith PetscFunctionReturn(0); 3385317021SBarry Smith } 3485317021SBarry Smith EXTERN_C_END 3585317021SBarry Smith 3685317021SBarry Smith EXTERN_C_BEGIN 3785317021SBarry Smith #undef __FUNCT__ 38d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType_Factor" 397087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype) 40d90ac83dSHong Zhang { 41d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 42d90ac83dSHong Zhang 43d90ac83dSHong Zhang PetscFunctionBegin; 44d90ac83dSHong Zhang if (shifttype == (MatFactorShiftType)PETSC_DECIDE){ 45f4db908eSBarry Smith dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 46d90ac83dSHong Zhang } else { 47f4db908eSBarry Smith dir->info.shifttype = (PetscReal) shifttype; 486cc3b3c1SHong Zhang if (shifttype == MAT_SHIFT_NONZERO && dir->info.shiftamount == 0.0){ 496cc3b3c1SHong Zhang dir->info.shiftamount = 1.e-12; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ 506cc3b3c1SHong Zhang } 51d90ac83dSHong Zhang } 52d90ac83dSHong Zhang PetscFunctionReturn(0); 53d90ac83dSHong Zhang } 54d90ac83dSHong Zhang 55d90ac83dSHong Zhang #undef __FUNCT__ 56d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount_Factor" 577087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount) 58d90ac83dSHong Zhang { 59d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 60d90ac83dSHong Zhang 61d90ac83dSHong Zhang PetscFunctionBegin; 62d90ac83dSHong Zhang if (shiftamount == (PetscReal) PETSC_DECIDE){ 63d90ac83dSHong Zhang dir->info.shiftamount = 1.e-12; 64d90ac83dSHong Zhang } else { 65d90ac83dSHong Zhang dir->info.shiftamount = shiftamount; 66d90ac83dSHong Zhang } 67d90ac83dSHong Zhang PetscFunctionReturn(0); 68d90ac83dSHong Zhang } 69d90ac83dSHong Zhang EXTERN_C_END 70d90ac83dSHong Zhang 71d90ac83dSHong Zhang EXTERN_C_BEGIN 72d90ac83dSHong Zhang #undef __FUNCT__ 73b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_Factor" 747087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 7585317021SBarry Smith { 7685317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 7785317021SBarry Smith 7885317021SBarry Smith PetscFunctionBegin; 7985317021SBarry 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)) { 80e7e72b3dSBarry Smith SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use"); 8185317021SBarry Smith } 8285317021SBarry Smith ilu->info.usedt = PETSC_TRUE; 8385317021SBarry Smith ilu->info.dt = dt; 8485317021SBarry Smith ilu->info.dtcol = dtcol; 8585317021SBarry Smith ilu->info.dtcount = dtcount; 8685317021SBarry Smith ilu->info.fill = PETSC_DEFAULT; 8785317021SBarry Smith PetscFunctionReturn(0); 8885317021SBarry Smith } 8985317021SBarry Smith EXTERN_C_END 9085317021SBarry Smith 9185317021SBarry Smith EXTERN_C_BEGIN 9285317021SBarry Smith #undef __FUNCT__ 9385317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor" 947087cfbeSBarry Smith PetscErrorCode PCFactorSetFill_Factor(PC pc,PetscReal fill) 9585317021SBarry Smith { 9685317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 9785317021SBarry Smith 9885317021SBarry Smith PetscFunctionBegin; 9985317021SBarry Smith dir->info.fill = fill; 10085317021SBarry Smith PetscFunctionReturn(0); 10185317021SBarry Smith } 10285317021SBarry Smith EXTERN_C_END 10385317021SBarry Smith 10485317021SBarry Smith EXTERN_C_BEGIN 10585317021SBarry Smith #undef __FUNCT__ 10685317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor" 1077087cfbeSBarry Smith PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering) 10885317021SBarry Smith { 10985317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 11085317021SBarry Smith PetscErrorCode ierr; 111ace3abfcSBarry Smith PetscBool flg; 11285317021SBarry Smith 11385317021SBarry Smith PetscFunctionBegin; 11485317021SBarry Smith if (!pc->setupcalled) { 115503cfb0cSBarry Smith ierr = PetscFree(dir->ordering);CHKERRQ(ierr); 11685317021SBarry Smith ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 11785317021SBarry Smith } else { 11885317021SBarry Smith ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 119e7e72b3dSBarry Smith if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering 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__ "PCFactorSetLevels_Factor" 1287087cfbeSBarry Smith PetscErrorCode PCFactorSetLevels_Factor(PC pc,PetscInt levels) 12985317021SBarry Smith { 13085317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 131*5b9c68c7SBarry Smith PetscErrorCode ierr; 13285317021SBarry Smith 13385317021SBarry Smith PetscFunctionBegin; 13485317021SBarry Smith if (!pc->setupcalled) { 13585317021SBarry Smith ilu->info.levels = levels; 136*5b9c68c7SBarry Smith } else if (ilu->info.levels != levels) { 137*5b9c68c7SBarry Smith ierr = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */ 138*5b9c68c7SBarry Smith pc->setupcalled = 0; /* force a complete rebuild of preconditioner factored matrices */ 139*5b9c68c7SBarry Smith ilu->info.levels = levels; 140*5b9c68c7SBarry Smith } else if (ilu->info.usedt) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt"); 14185317021SBarry Smith PetscFunctionReturn(0); 14285317021SBarry Smith } 14385317021SBarry Smith EXTERN_C_END 14485317021SBarry Smith 14585317021SBarry Smith EXTERN_C_BEGIN 14685317021SBarry Smith #undef __FUNCT__ 14785317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor" 1487087cfbeSBarry Smith PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc) 14985317021SBarry Smith { 15085317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 15185317021SBarry Smith 15285317021SBarry Smith PetscFunctionBegin; 15385317021SBarry Smith dir->info.diagonal_fill = 1; 15485317021SBarry Smith PetscFunctionReturn(0); 15585317021SBarry Smith } 15685317021SBarry Smith EXTERN_C_END 15785317021SBarry Smith 15885317021SBarry Smith 15985317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 16085317021SBarry Smith 16185317021SBarry Smith EXTERN_C_BEGIN 16285317021SBarry Smith #undef __FUNCT__ 16385317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor" 1647087cfbeSBarry Smith PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot) 16585317021SBarry Smith { 16685317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 16785317021SBarry Smith 16885317021SBarry Smith PetscFunctionBegin; 16985317021SBarry Smith dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 17085317021SBarry Smith PetscFunctionReturn(0); 17185317021SBarry Smith } 17285317021SBarry Smith EXTERN_C_END 17385317021SBarry Smith 17485317021SBarry Smith #undef __FUNCT__ 17585317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor" 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; 181e7e72b3dSBarry Smith if (!ilu->fact) SETERRQ(((PetscObject)pc)->comm,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 18685317021SBarry Smith EXTERN_C_BEGIN 18785317021SBarry Smith #undef __FUNCT__ 18885317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor" 1897087cfbeSBarry Smith PetscErrorCode PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype) 19085317021SBarry Smith { 19185317021SBarry Smith PetscErrorCode ierr; 19285317021SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 19385317021SBarry Smith 19485317021SBarry Smith PetscFunctionBegin; 1957112b564SBarry Smith if (lu->fact) { 19685317021SBarry Smith const MatSolverPackage ltype; 197ace3abfcSBarry Smith PetscBool flg; 19885317021SBarry Smith ierr = MatFactorGetSolverPackage(lu->fact,<ype);CHKERRQ(ierr); 19985317021SBarry Smith ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 200e7e72b3dSBarry Smith if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used"); 201a1f19f5aSHong Zhang } else { 202503cfb0cSBarry Smith ierr = PetscFree(lu->solvertype);CHKERRQ(ierr); 20385317021SBarry Smith ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr); 204f8260c8fSBarry Smith } 20585317021SBarry Smith PetscFunctionReturn(0); 20685317021SBarry Smith } 20785317021SBarry Smith EXTERN_C_END 20885317021SBarry Smith 20985317021SBarry Smith EXTERN_C_BEGIN 21085317021SBarry Smith #undef __FUNCT__ 2117112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor" 2127087cfbeSBarry Smith PetscErrorCode PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype) 2137112b564SBarry Smith { 2147112b564SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 2157112b564SBarry Smith 2167112b564SBarry Smith PetscFunctionBegin; 2177112b564SBarry Smith *stype = lu->solvertype; 2187112b564SBarry Smith PetscFunctionReturn(0); 2197112b564SBarry Smith } 2207112b564SBarry Smith EXTERN_C_END 2217112b564SBarry Smith 2227112b564SBarry Smith EXTERN_C_BEGIN 2237112b564SBarry Smith #undef __FUNCT__ 2248ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot_Factor" 2257087cfbeSBarry Smith PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol) 22685317021SBarry Smith { 22785317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 22885317021SBarry Smith 22985317021SBarry Smith PetscFunctionBegin; 23065e19b50SBarry Smith if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 23185317021SBarry Smith dir->info.dtcol = dtcol; 23285317021SBarry Smith PetscFunctionReturn(0); 23385317021SBarry Smith } 23485317021SBarry Smith EXTERN_C_END 2358ff23777SHong Zhang 2368ff23777SHong Zhang EXTERN_C_BEGIN 2378ff23777SHong Zhang #undef __FUNCT__ 2388ff23777SHong Zhang #define __FUNCT__ "PCSetFromOptions_Factor" 2397087cfbeSBarry Smith PetscErrorCode PCSetFromOptions_Factor(PC pc) 2408ff23777SHong Zhang { 2418ff23777SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 2428ff23777SHong Zhang PetscErrorCode ierr; 243ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE,set; 2448ff23777SHong Zhang char tname[256], solvertype[64]; 2458ff23777SHong Zhang PetscFList ordlist; 246018dd85eSSatish Balay PetscEnum etmp; 2478ff23777SHong Zhang 2488ff23777SHong Zhang PetscFunctionBegin; 2498ff23777SHong Zhang if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 250acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 2518ff23777SHong Zhang if (flg) { 2528ff23777SHong Zhang ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 2538ff23777SHong Zhang } 2548ff23777SHong 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); 2558ff23777SHong Zhang 2561ae3d29cSBarry Smith ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType", 257018dd85eSSatish Balay MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr); 258018dd85eSSatish Balay if (flg) { 2596cc3b3c1SHong Zhang ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr); 260018dd85eSSatish Balay } 261d90ac83dSHong 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); 262d90ac83dSHong Zhang 2638ff23777SHong 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); 2648ff23777SHong 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); 2658ff23777SHong Zhang 2668ff23777SHong Zhang flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 267acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 2688ff23777SHong Zhang if (set) { 2698ff23777SHong Zhang ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 2708ff23777SHong Zhang } 2718ff23777SHong Zhang 2728ff23777SHong Zhang flg = PETSC_FALSE; 273acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 2748ff23777SHong Zhang if (flg) { 2758ff23777SHong Zhang ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 2768ff23777SHong Zhang } 2778ff23777SHong Zhang flg = PETSC_FALSE; 278acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 2798ff23777SHong Zhang if (flg) { 2808ff23777SHong Zhang ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 2818ff23777SHong Zhang } 2828ff23777SHong Zhang 2838ff23777SHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 2848ff23777SHong 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); 2858ff23777SHong Zhang if (flg) { 2868ff23777SHong Zhang ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 2878ff23777SHong Zhang } 2888ff23777SHong Zhang 2898ff23777SHong Zhang /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 2908ff23777SHong Zhang ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 2918ff23777SHong Zhang if (flg) { 2928ff23777SHong Zhang ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 2938ff23777SHong Zhang } 2948ff23777SHong Zhang PetscFunctionReturn(0); 2958ff23777SHong Zhang } 2968ff23777SHong Zhang EXTERN_C_END 297914a5d51SHong Zhang 298914a5d51SHong Zhang EXTERN_C_BEGIN 299914a5d51SHong Zhang #undef __FUNCT__ 300914a5d51SHong Zhang #define __FUNCT__ "PCView_Factor" 301914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer) 302914a5d51SHong Zhang { 303914a5d51SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 304914a5d51SHong Zhang PetscErrorCode ierr; 305ace3abfcSBarry Smith PetscBool isstring,iascii; 306914a5d51SHong Zhang 307914a5d51SHong Zhang PetscFunctionBegin; 3082692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 3092692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 310914a5d51SHong Zhang if (iascii) { 311879e8a4dSBarry Smith if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC){ 312914a5d51SHong Zhang if (factor->info.dt > 0) { 313914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," drop tolerance %G\n",factor->info.dt);CHKERRQ(ierr); 314914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr); 315914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);CHKERRQ(ierr); 316914a5d51SHong Zhang } else if (factor->info.levels == 1) { 317914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 318914a5d51SHong Zhang } else { 319914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 320914a5d51SHong Zhang } 321914a5d51SHong Zhang } 322914a5d51SHong Zhang 323914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," tolerance for zero pivot %G\n",factor->info.zeropivot);CHKERRQ(ierr); 324f4db908eSBarry Smith if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 325914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," using Manteuffel shift\n");CHKERRQ(ierr); 326914a5d51SHong Zhang } 327f4db908eSBarry Smith if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) { 328914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr); 329914a5d51SHong Zhang } 330f4db908eSBarry Smith if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) { 331914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr); 332914a5d51SHong Zhang } 333914a5d51SHong Zhang 334914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr); 335914a5d51SHong Zhang 336914a5d51SHong Zhang if (factor->fact) { 3376335c336SBarry Smith MatInfo info; 3386335c336SBarry Smith ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 3396335c336SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," factor fill ratio given %G, needed %G\n",info.fill_ratio_given,info.fill_ratio_needed);CHKERRQ(ierr); 340914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows:\n");CHKERRQ(ierr); 341914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 342914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 343914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 344914a5d51SHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 345914a5d51SHong Zhang ierr = MatView(factor->fact,viewer);CHKERRQ(ierr); 346914a5d51SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 347914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 348914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 349914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 350914a5d51SHong Zhang } 351914a5d51SHong Zhang 352914a5d51SHong Zhang } else if (isstring) { 3536335c336SBarry Smith MatFactorType t; 3546335c336SBarry Smith ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr); 3556335c336SBarry Smith if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC){ 356914a5d51SHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 357914a5d51SHong Zhang } 358914a5d51SHong Zhang } else { 35965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC_Factor",((PetscObject)viewer)->type_name); 360914a5d51SHong Zhang } 361914a5d51SHong Zhang PetscFunctionReturn(0); 362914a5d51SHong Zhang } 363914a5d51SHong Zhang EXTERN_C_END 364