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 PetscErrorCode ierr; 10f8260c8fSBarry Smith 11f8260c8fSBarry Smith PetscFunctionBegin; 12*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!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()"); 13f8260c8fSBarry Smith if (!pc->setupcalled && !((PC_Factor*)icc)->fact) { 14f8260c8fSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 15f8260c8fSBarry Smith } 16f8260c8fSBarry Smith PetscFunctionReturn(0); 17f8260c8fSBarry Smith } 18f8260c8fSBarry Smith 197087cfbeSBarry Smith PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc,PetscReal z) 2085317021SBarry Smith { 2185317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 2285317021SBarry Smith 2385317021SBarry Smith PetscFunctionBegin; 2485317021SBarry Smith ilu->info.zeropivot = z; 2585317021SBarry Smith PetscFunctionReturn(0); 2685317021SBarry Smith } 2785317021SBarry Smith 287087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype) 29d90ac83dSHong Zhang { 30d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 31d90ac83dSHong Zhang 32d90ac83dSHong Zhang PetscFunctionBegin; 332fa5cd67SKarl Rupp if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 342fa5cd67SKarl Rupp else { 35f4db908eSBarry Smith dir->info.shifttype = (PetscReal) shifttype; 369dbfc187SHong Zhang if ((shifttype == MAT_SHIFT_NONZERO || shifttype == MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) { 379dbfc187SHong Zhang dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ 386cc3b3c1SHong Zhang } 39d90ac83dSHong Zhang } 40d90ac83dSHong Zhang PetscFunctionReturn(0); 41d90ac83dSHong Zhang } 42d90ac83dSHong Zhang 437087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount) 44d90ac83dSHong Zhang { 45d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 46d90ac83dSHong Zhang 47d90ac83dSHong Zhang PetscFunctionBegin; 482fa5cd67SKarl Rupp if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 492fa5cd67SKarl Rupp else dir->info.shiftamount = shiftamount; 50d90ac83dSHong Zhang PetscFunctionReturn(0); 51d90ac83dSHong Zhang } 52d90ac83dSHong Zhang 537087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 5485317021SBarry Smith { 5585317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 5685317021SBarry Smith 5785317021SBarry Smith PetscFunctionBegin; 5885317021SBarry 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)) { 59ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use"); 6085317021SBarry Smith } 6185317021SBarry Smith ilu->info.usedt = PETSC_TRUE; 6285317021SBarry Smith ilu->info.dt = dt; 6385317021SBarry Smith ilu->info.dtcol = dtcol; 6485317021SBarry Smith ilu->info.dtcount = dtcount; 6585317021SBarry Smith ilu->info.fill = PETSC_DEFAULT; 6685317021SBarry Smith PetscFunctionReturn(0); 6785317021SBarry Smith } 6885317021SBarry Smith 697087cfbeSBarry Smith PetscErrorCode PCFactorSetFill_Factor(PC pc,PetscReal fill) 7085317021SBarry Smith { 7185317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 7285317021SBarry Smith 7385317021SBarry Smith PetscFunctionBegin; 7485317021SBarry Smith dir->info.fill = fill; 7585317021SBarry Smith PetscFunctionReturn(0); 7685317021SBarry Smith } 7785317021SBarry Smith 7819fd82e9SBarry Smith PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering) 7985317021SBarry Smith { 8085317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 8185317021SBarry Smith PetscErrorCode ierr; 82ace3abfcSBarry Smith PetscBool flg; 8385317021SBarry Smith 8485317021SBarry Smith PetscFunctionBegin; 8585317021SBarry Smith if (!pc->setupcalled) { 86503cfb0cSBarry Smith ierr = PetscFree(dir->ordering);CHKERRQ(ierr); 8719fd82e9SBarry Smith ierr = PetscStrallocpy(ordering,(char**)&dir->ordering);CHKERRQ(ierr); 8885317021SBarry Smith } else { 8985317021SBarry Smith ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 90*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use"); 9185317021SBarry Smith } 9285317021SBarry Smith PetscFunctionReturn(0); 9385317021SBarry Smith } 9485317021SBarry Smith 952591b318SToby Isaac PetscErrorCode PCFactorGetLevels_Factor(PC pc,PetscInt *levels) 962591b318SToby Isaac { 972591b318SToby Isaac PC_Factor *ilu = (PC_Factor*)pc->data; 982591b318SToby Isaac 992591b318SToby Isaac PetscFunctionBegin; 1002591b318SToby Isaac *levels = ilu->info.levels; 1012591b318SToby Isaac PetscFunctionReturn(0); 1022591b318SToby Isaac } 1032591b318SToby Isaac 104c7f610a1SBarry Smith PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc,PetscReal *pivot) 105c7f610a1SBarry Smith { 106c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 107c7f610a1SBarry Smith 108c7f610a1SBarry Smith PetscFunctionBegin; 109c7f610a1SBarry Smith *pivot = ilu->info.zeropivot; 110c7f610a1SBarry Smith PetscFunctionReturn(0); 111c7f610a1SBarry Smith } 112c7f610a1SBarry Smith 113c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc,PetscReal *shift) 114c7f610a1SBarry Smith { 115c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 116c7f610a1SBarry Smith 117c7f610a1SBarry Smith PetscFunctionBegin; 118c7f610a1SBarry Smith *shift = ilu->info.shiftamount; 119c7f610a1SBarry Smith PetscFunctionReturn(0); 120c7f610a1SBarry Smith } 121c7f610a1SBarry Smith 122c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftType_Factor(PC pc,MatFactorShiftType *type) 123c7f610a1SBarry Smith { 124c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 125c7f610a1SBarry Smith 126c7f610a1SBarry Smith PetscFunctionBegin; 127b729e602SBarry Smith *type = (MatFactorShiftType) (int) ilu->info.shifttype; 128c7f610a1SBarry Smith PetscFunctionReturn(0); 129c7f610a1SBarry Smith } 130c7f610a1SBarry Smith 1317087cfbeSBarry Smith PetscErrorCode PCFactorSetLevels_Factor(PC pc,PetscInt levels) 13285317021SBarry Smith { 13385317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 1345b9c68c7SBarry Smith PetscErrorCode ierr; 13585317021SBarry Smith 13685317021SBarry Smith PetscFunctionBegin; 1372fa5cd67SKarl Rupp if (!pc->setupcalled) ilu->info.levels = levels; 1382fa5cd67SKarl Rupp else if (ilu->info.levels != levels) { 1395b9c68c7SBarry Smith ierr = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */ 1405b9c68c7SBarry Smith pc->setupcalled = 0; /* force a complete rebuild of preconditioner factored matrices */ 1415b9c68c7SBarry Smith ilu->info.levels = levels; 142*2c71b3e2SJacob Faibussowitsch } else PetscCheckFalse(ilu->info.usedt,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt"); 14385317021SBarry Smith PetscFunctionReturn(0); 14485317021SBarry Smith } 14585317021SBarry Smith 14692e9c092SBarry Smith PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg) 14785317021SBarry Smith { 14885317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 14985317021SBarry Smith 15085317021SBarry Smith PetscFunctionBegin; 15192e9c092SBarry Smith dir->info.diagonal_fill = (PetscReal) flg; 15292e9c092SBarry Smith PetscFunctionReturn(0); 15392e9c092SBarry Smith } 15492e9c092SBarry Smith 15592e9c092SBarry Smith PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg) 15692e9c092SBarry Smith { 15792e9c092SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 15892e9c092SBarry Smith 15992e9c092SBarry Smith PetscFunctionBegin; 16092e9c092SBarry Smith *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE; 16185317021SBarry Smith PetscFunctionReturn(0); 16285317021SBarry Smith } 16385317021SBarry Smith 16485317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 16585317021SBarry Smith 1667087cfbeSBarry Smith PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot) 16785317021SBarry Smith { 16885317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 16985317021SBarry Smith 17085317021SBarry Smith PetscFunctionBegin; 17185317021SBarry Smith dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 17285317021SBarry Smith PetscFunctionReturn(0); 17385317021SBarry Smith } 17485317021SBarry Smith 1757087cfbeSBarry Smith PetscErrorCode PCFactorGetMatrix_Factor(PC pc,Mat *mat) 17685317021SBarry Smith { 17785317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 17885317021SBarry Smith 17985317021SBarry Smith PetscFunctionBegin; 180*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!ilu->fact,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 18185317021SBarry Smith *mat = ilu->fact; 18285317021SBarry Smith PetscFunctionReturn(0); 18385317021SBarry Smith } 18485317021SBarry Smith 185978beb7fSPierre Jolivet /* allow access to preallocation information */ 186978beb7fSPierre Jolivet #include <petsc/private/matimpl.h> 187978beb7fSPierre Jolivet 188ea799195SBarry Smith PetscErrorCode PCFactorSetMatSolverType_Factor(PC pc,MatSolverType stype) 18985317021SBarry Smith { 19085317021SBarry Smith PetscErrorCode ierr; 19185317021SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 19285317021SBarry Smith 19385317021SBarry Smith PetscFunctionBegin; 194978beb7fSPierre Jolivet if (lu->fact && lu->fact->assembled) { 195ea799195SBarry Smith MatSolverType ltype; 196ace3abfcSBarry Smith PetscBool flg; 1973ca39a21SBarry Smith ierr = MatFactorGetSolverType(lu->fact,<ype);CHKERRQ(ierr); 19885317021SBarry Smith ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 199*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!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); 20000c67f3bSHong Zhang } 20100c67f3bSHong Zhang 202503cfb0cSBarry Smith ierr = PetscFree(lu->solvertype);CHKERRQ(ierr); 20385317021SBarry Smith ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr); 20485317021SBarry Smith PetscFunctionReturn(0); 20585317021SBarry Smith } 20685317021SBarry Smith 207ea799195SBarry Smith PetscErrorCode PCFactorGetMatSolverType_Factor(PC pc,MatSolverType *stype) 2087112b564SBarry Smith { 2097112b564SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 2107112b564SBarry Smith 2117112b564SBarry Smith PetscFunctionBegin; 2127112b564SBarry Smith *stype = lu->solvertype; 2137112b564SBarry Smith PetscFunctionReturn(0); 2147112b564SBarry Smith } 2157112b564SBarry Smith 2167087cfbeSBarry Smith PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol) 21785317021SBarry Smith { 21885317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 21985317021SBarry Smith 22085317021SBarry Smith PetscFunctionBegin; 221*2c71b3e2SJacob 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); 22285317021SBarry Smith dir->info.dtcol = dtcol; 22385317021SBarry Smith PetscFunctionReturn(0); 22485317021SBarry Smith } 2258ff23777SHong Zhang 2264416b707SBarry Smith PetscErrorCode PCSetFromOptions_Factor(PetscOptionItems *PetscOptionsObject,PC pc) 2278ff23777SHong Zhang { 2288ff23777SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 2298ff23777SHong Zhang PetscErrorCode ierr; 2308afaa268SBarry Smith PetscBool flg,set; 2318ff23777SHong Zhang char tname[256], solvertype[64]; 232140e18c1SBarry Smith PetscFunctionList ordlist; 233018dd85eSSatish Balay PetscEnum etmp; 2348e37d05fSBarry Smith PetscBool inplace; 2358ff23777SHong Zhang 2368ff23777SHong Zhang PetscFunctionBegin; 2378e37d05fSBarry Smith ierr = PCFactorGetUseInPlace(pc,&inplace);CHKERRQ(ierr); 23892e9c092SBarry Smith ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set);CHKERRQ(ierr); 2398e37d05fSBarry Smith if (set) { 2408e37d05fSBarry Smith ierr = PCFactorSetUseInPlace(pc,flg);CHKERRQ(ierr); 2418ff23777SHong Zhang } 2428afaa268SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,NULL);CHKERRQ(ierr); 2438ff23777SHong Zhang 2448afaa268SBarry Smith ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr); 245018dd85eSSatish Balay if (flg) { 2466cc3b3c1SHong Zhang ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr); 247018dd85eSSatish Balay } 2480a545947SLisandro Dalcin ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,NULL);CHKERRQ(ierr); 249d90ac83dSHong Zhang 2500a545947SLisandro Dalcin ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,NULL);CHKERRQ(ierr); 2518ff23777SHong 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); 2528ff23777SHong Zhang 2538afaa268SBarry Smith ierr = 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);CHKERRQ(ierr); 2548ff23777SHong Zhang if (set) { 2558ff23777SHong Zhang ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 2568ff23777SHong Zhang } 2578ff23777SHong Zhang 2588afaa268SBarry Smith ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 2598afaa268SBarry Smith if (set) { 2608afaa268SBarry Smith ierr = PCFactorSetReuseFill(pc,flg);CHKERRQ(ierr); 2618ff23777SHong Zhang } 2628afaa268SBarry Smith ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 2638afaa268SBarry Smith if (set) { 2648afaa268SBarry Smith ierr = PCFactorSetReuseOrdering(pc,flg);CHKERRQ(ierr); 2658ff23777SHong Zhang } 2668ff23777SHong Zhang 2679f3a6782SPatrick Sanan ierr = PetscOptionsDeprecated("-pc_factor_mat_solver_package","-pc_factor_mat_solver_type","3.9",NULL);CHKERRQ(ierr); 268589a23caSBarry Smith ierr = PetscOptionsString("-pc_factor_mat_solver_type","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,sizeof(solvertype),&flg);CHKERRQ(ierr); 2698ff23777SHong Zhang if (flg) { 2703ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solvertype);CHKERRQ(ierr); 2718ff23777SHong Zhang } 2724ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 2734ac6704cSBarry Smith ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 2744ac6704cSBarry Smith ierr = PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,sizeof(tname),&flg);CHKERRQ(ierr); 2754ac6704cSBarry Smith if (flg) { 2764ac6704cSBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 2774ac6704cSBarry Smith } 2788ff23777SHong Zhang PetscFunctionReturn(0); 2798ff23777SHong Zhang } 280914a5d51SHong Zhang 281914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer) 282914a5d51SHong Zhang { 283914a5d51SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 284914a5d51SHong Zhang PetscErrorCode ierr; 2854ac6704cSBarry Smith PetscBool isstring,iascii,canuseordering; 2864ac6704cSBarry Smith MatInfo info; 2879bd791bbSBarry Smith MatOrderingType ordering; 288914a5d51SHong Zhang 289914a5d51SHong Zhang PetscFunctionBegin; 290251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 291251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 292914a5d51SHong Zhang if (iascii) { 2933d1c1ea0SBarry Smith if (factor->inplace) { 2943d1c1ea0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr); 2953d1c1ea0SBarry Smith } else { 2963d1c1ea0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr); 2973d1c1ea0SBarry Smith } 2983d1c1ea0SBarry Smith 2993d1c1ea0SBarry Smith if (factor->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 3003d1c1ea0SBarry Smith if (factor->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 301879e8a4dSBarry Smith if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) { 302914a5d51SHong Zhang if (factor->info.dt > 0) { 30357622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," drop tolerance %g\n",(double)factor->info.dt);CHKERRQ(ierr); 304914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr); 30557622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," column permutation tolerance %g\n",(double)factor->info.dtcol);CHKERRQ(ierr); 306914a5d51SHong Zhang } else if (factor->info.levels == 1) { 307914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 308914a5d51SHong Zhang } else { 309914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 310914a5d51SHong Zhang } 311914a5d51SHong Zhang } 312914a5d51SHong Zhang 31357622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerance for zero pivot %g\n",(double)factor->info.zeropivot);CHKERRQ(ierr); 3145e9742b9SJed Brown if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */ 3155e9742b9SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);CHKERRQ(ierr); 316914a5d51SHong Zhang } 317914a5d51SHong Zhang 3184ac6704cSBarry Smith if (factor->fact) { 319f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(factor->fact,&canuseordering);CHKERRQ(ierr); 3204ac6704cSBarry Smith if (!canuseordering) ordering = MATORDERINGEXTERNAL; 3214ac6704cSBarry Smith else ordering = factor->ordering; 3229bd791bbSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ordering);CHKERRQ(ierr); 323978beb7fSPierre Jolivet if (!factor->fact->assembled) { 324978beb7fSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," matrix solver type: %s\n",factor->fact->solvertype);CHKERRQ(ierr); 325978beb7fSPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," matrix not yet factored; no additional information available\n");CHKERRQ(ierr); 326978beb7fSPierre Jolivet } else { 3276335c336SBarry Smith ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 32857622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);CHKERRQ(ierr); 329914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows:\n");CHKERRQ(ierr); 330914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 331914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 332914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 333914a5d51SHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 334914a5d51SHong Zhang ierr = MatView(factor->fact,viewer);CHKERRQ(ierr); 335914a5d51SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 336914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 337914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 338914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 339914a5d51SHong Zhang } 340978beb7fSPierre Jolivet } 341914a5d51SHong Zhang 342914a5d51SHong Zhang } else if (isstring) { 3436335c336SBarry Smith MatFactorType t; 3446335c336SBarry Smith ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr); 3456335c336SBarry Smith if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) { 346c3b366b1Sprj- ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr); 347914a5d51SHong Zhang } 348914a5d51SHong Zhang } 349914a5d51SHong Zhang PetscFunctionReturn(0); 350914a5d51SHong Zhang } 351