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 #undef __FUNCT__ 8f8260c8fSBarry Smith #define __FUNCT__ "PCFactorSetUpMatSolverPackage_Factor" 9f8260c8fSBarry Smith PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc) 10f8260c8fSBarry Smith { 11f8260c8fSBarry Smith PC_Factor *icc = (PC_Factor*)pc->data; 12f8260c8fSBarry Smith PetscErrorCode ierr; 13f8260c8fSBarry Smith 14f8260c8fSBarry Smith PetscFunctionBegin; 15ad5697c6SBarry Smith if (!pc->pmat) SETERRQ(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()"); 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 2285317021SBarry Smith #undef __FUNCT__ 2385317021SBarry Smith #define __FUNCT__ "PCFactorSetZeroPivot_Factor" 247087cfbeSBarry Smith PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc,PetscReal z) 2585317021SBarry Smith { 2685317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 2785317021SBarry Smith 2885317021SBarry Smith PetscFunctionBegin; 2985317021SBarry Smith ilu->info.zeropivot = z; 3085317021SBarry Smith PetscFunctionReturn(0); 3185317021SBarry Smith } 3285317021SBarry Smith 3385317021SBarry Smith #undef __FUNCT__ 34d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftType_Factor" 357087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype) 36d90ac83dSHong Zhang { 37d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 38d90ac83dSHong Zhang 39d90ac83dSHong Zhang PetscFunctionBegin; 402fa5cd67SKarl Rupp if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 412fa5cd67SKarl Rupp else { 42f4db908eSBarry Smith dir->info.shifttype = (PetscReal) shifttype; 439dbfc187SHong Zhang if ((shifttype == MAT_SHIFT_NONZERO || shifttype == MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) { 449dbfc187SHong Zhang dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */ 456cc3b3c1SHong Zhang } 46d90ac83dSHong Zhang } 47d90ac83dSHong Zhang PetscFunctionReturn(0); 48d90ac83dSHong Zhang } 49d90ac83dSHong Zhang 50d90ac83dSHong Zhang #undef __FUNCT__ 51d90ac83dSHong Zhang #define __FUNCT__ "PCFactorSetShiftAmount_Factor" 527087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount) 53d90ac83dSHong Zhang { 54d90ac83dSHong Zhang PC_Factor *dir = (PC_Factor*)pc->data; 55d90ac83dSHong Zhang 56d90ac83dSHong Zhang PetscFunctionBegin; 572fa5cd67SKarl Rupp if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 582fa5cd67SKarl Rupp else dir->info.shiftamount = shiftamount; 59d90ac83dSHong Zhang PetscFunctionReturn(0); 60d90ac83dSHong Zhang } 61d90ac83dSHong Zhang 62d90ac83dSHong Zhang #undef __FUNCT__ 63b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_Factor" 647087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 6585317021SBarry Smith { 6685317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 6785317021SBarry Smith 6885317021SBarry Smith PetscFunctionBegin; 6985317021SBarry 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)) { 70ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use"); 7185317021SBarry Smith } 7285317021SBarry Smith ilu->info.usedt = PETSC_TRUE; 7385317021SBarry Smith ilu->info.dt = dt; 7485317021SBarry Smith ilu->info.dtcol = dtcol; 7585317021SBarry Smith ilu->info.dtcount = dtcount; 7685317021SBarry Smith ilu->info.fill = PETSC_DEFAULT; 7785317021SBarry Smith PetscFunctionReturn(0); 7885317021SBarry Smith } 7985317021SBarry Smith 8085317021SBarry Smith #undef __FUNCT__ 8185317021SBarry Smith #define __FUNCT__ "PCFactorSetFill_Factor" 827087cfbeSBarry Smith PetscErrorCode PCFactorSetFill_Factor(PC pc,PetscReal fill) 8385317021SBarry Smith { 8485317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 8585317021SBarry Smith 8685317021SBarry Smith PetscFunctionBegin; 8785317021SBarry Smith dir->info.fill = fill; 8885317021SBarry Smith PetscFunctionReturn(0); 8985317021SBarry Smith } 9085317021SBarry Smith 9185317021SBarry Smith #undef __FUNCT__ 9285317021SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Factor" 9319fd82e9SBarry Smith PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering) 9485317021SBarry Smith { 9585317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 9685317021SBarry Smith PetscErrorCode ierr; 97ace3abfcSBarry Smith PetscBool flg; 9885317021SBarry Smith 9985317021SBarry Smith PetscFunctionBegin; 10085317021SBarry Smith if (!pc->setupcalled) { 101503cfb0cSBarry Smith ierr = PetscFree(dir->ordering);CHKERRQ(ierr); 10219fd82e9SBarry Smith ierr = PetscStrallocpy(ordering,(char**)&dir->ordering);CHKERRQ(ierr); 10385317021SBarry Smith } else { 10485317021SBarry Smith ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 105ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use"); 10685317021SBarry Smith } 10785317021SBarry Smith PetscFunctionReturn(0); 10885317021SBarry Smith } 10985317021SBarry Smith 11085317021SBarry Smith #undef __FUNCT__ 1112591b318SToby Isaac #define __FUNCT__ "PCFactorGetLevels_Factor" 1122591b318SToby Isaac PetscErrorCode PCFactorGetLevels_Factor(PC pc,PetscInt *levels) 1132591b318SToby Isaac { 1142591b318SToby Isaac PC_Factor *ilu = (PC_Factor*)pc->data; 1152591b318SToby Isaac 1162591b318SToby Isaac PetscFunctionBegin; 1172591b318SToby Isaac *levels = ilu->info.levels; 1182591b318SToby Isaac PetscFunctionReturn(0); 1192591b318SToby Isaac } 1202591b318SToby Isaac 1212591b318SToby Isaac #undef __FUNCT__ 122c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetZeroPivot_Factor" 123c7f610a1SBarry Smith PetscErrorCode PCFactorGetZeroPivot_Factor(PC pc,PetscReal *pivot) 124c7f610a1SBarry Smith { 125c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 126c7f610a1SBarry Smith 127c7f610a1SBarry Smith PetscFunctionBegin; 128c7f610a1SBarry Smith *pivot = ilu->info.zeropivot; 129c7f610a1SBarry Smith PetscFunctionReturn(0); 130c7f610a1SBarry Smith } 131c7f610a1SBarry Smith 132c7f610a1SBarry Smith #undef __FUNCT__ 133c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetShiftAmount_Factor" 134c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftAmount_Factor(PC pc,PetscReal *shift) 135c7f610a1SBarry Smith { 136c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 137c7f610a1SBarry Smith 138c7f610a1SBarry Smith PetscFunctionBegin; 139c7f610a1SBarry Smith *shift = ilu->info.shiftamount; 140c7f610a1SBarry Smith PetscFunctionReturn(0); 141c7f610a1SBarry Smith } 142c7f610a1SBarry Smith 143c7f610a1SBarry Smith #undef __FUNCT__ 144c7f610a1SBarry Smith #define __FUNCT__ "PCFactorGetShiftTyoe_Factor" 145c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftType_Factor(PC pc,MatFactorShiftType *type) 146c7f610a1SBarry Smith { 147c7f610a1SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 148c7f610a1SBarry Smith 149c7f610a1SBarry Smith PetscFunctionBegin; 150c7f610a1SBarry Smith *type = (MatFactorShiftType) ilu->info.shifttype; 151c7f610a1SBarry Smith PetscFunctionReturn(0); 152c7f610a1SBarry Smith } 153c7f610a1SBarry Smith 154c7f610a1SBarry Smith #undef __FUNCT__ 15585317021SBarry Smith #define __FUNCT__ "PCFactorSetLevels_Factor" 1567087cfbeSBarry Smith PetscErrorCode PCFactorSetLevels_Factor(PC pc,PetscInt levels) 15785317021SBarry Smith { 15885317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 1595b9c68c7SBarry Smith PetscErrorCode ierr; 16085317021SBarry Smith 16185317021SBarry Smith PetscFunctionBegin; 1622fa5cd67SKarl Rupp if (!pc->setupcalled) ilu->info.levels = levels; 1632fa5cd67SKarl Rupp else if (ilu->info.levels != levels) { 1645b9c68c7SBarry Smith ierr = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */ 1655b9c68c7SBarry Smith pc->setupcalled = 0; /* force a complete rebuild of preconditioner factored matrices */ 1665b9c68c7SBarry Smith ilu->info.levels = levels; 167ce94432eSBarry Smith } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt"); 16885317021SBarry Smith PetscFunctionReturn(0); 16985317021SBarry Smith } 17085317021SBarry Smith 17185317021SBarry Smith #undef __FUNCT__ 17285317021SBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor" 17392e9c092SBarry Smith PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg) 17485317021SBarry Smith { 17585317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 17685317021SBarry Smith 17785317021SBarry Smith PetscFunctionBegin; 17892e9c092SBarry Smith dir->info.diagonal_fill = (PetscReal) flg; 17992e9c092SBarry Smith PetscFunctionReturn(0); 18092e9c092SBarry Smith } 18192e9c092SBarry Smith 18292e9c092SBarry Smith #undef __FUNCT__ 18392e9c092SBarry Smith #define __FUNCT__ "PCFactorGetAllowDiagonalFill_Factor" 18492e9c092SBarry Smith PetscErrorCode PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg) 18592e9c092SBarry Smith { 18692e9c092SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 18792e9c092SBarry Smith 18892e9c092SBarry Smith PetscFunctionBegin; 18992e9c092SBarry Smith *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE; 19085317021SBarry Smith PetscFunctionReturn(0); 19185317021SBarry Smith } 19285317021SBarry Smith 19385317021SBarry Smith /* ------------------------------------------------------------------------------------------*/ 19485317021SBarry Smith 19585317021SBarry Smith #undef __FUNCT__ 19685317021SBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor" 1977087cfbeSBarry Smith PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot) 19885317021SBarry Smith { 19985317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 20085317021SBarry Smith 20185317021SBarry Smith PetscFunctionBegin; 20285317021SBarry Smith dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 20385317021SBarry Smith PetscFunctionReturn(0); 20485317021SBarry Smith } 20585317021SBarry Smith 20685317021SBarry Smith #undef __FUNCT__ 20785317021SBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Factor" 2087087cfbeSBarry Smith PetscErrorCode PCFactorGetMatrix_Factor(PC pc,Mat *mat) 20985317021SBarry Smith { 21085317021SBarry Smith PC_Factor *ilu = (PC_Factor*)pc->data; 21185317021SBarry Smith 21285317021SBarry Smith PetscFunctionBegin; 213ce94432eSBarry Smith if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 21485317021SBarry Smith *mat = ilu->fact; 21585317021SBarry Smith PetscFunctionReturn(0); 21685317021SBarry Smith } 21785317021SBarry Smith 21885317021SBarry Smith #undef __FUNCT__ 21985317021SBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor" 2207087cfbeSBarry Smith PetscErrorCode PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype) 22185317021SBarry Smith { 22285317021SBarry Smith PetscErrorCode ierr; 22385317021SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 22485317021SBarry Smith 22585317021SBarry Smith PetscFunctionBegin; 2267112b564SBarry Smith if (lu->fact) { 22785317021SBarry Smith const MatSolverPackage ltype; 228ace3abfcSBarry Smith PetscBool flg; 22985317021SBarry Smith ierr = MatFactorGetSolverPackage(lu->fact,<ype);CHKERRQ(ierr); 23085317021SBarry Smith ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 231ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used"); 23200c67f3bSHong Zhang } 23300c67f3bSHong Zhang 234503cfb0cSBarry Smith ierr = PetscFree(lu->solvertype);CHKERRQ(ierr); 23585317021SBarry Smith ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr); 23685317021SBarry Smith PetscFunctionReturn(0); 23785317021SBarry Smith } 23885317021SBarry Smith 23985317021SBarry Smith #undef __FUNCT__ 2407112b564SBarry Smith #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor" 2417087cfbeSBarry Smith PetscErrorCode PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype) 2427112b564SBarry Smith { 2437112b564SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 2447112b564SBarry Smith 2457112b564SBarry Smith PetscFunctionBegin; 2467112b564SBarry Smith *stype = lu->solvertype; 2477112b564SBarry Smith PetscFunctionReturn(0); 2487112b564SBarry Smith } 2497112b564SBarry Smith 2507112b564SBarry Smith #undef __FUNCT__ 2518ff23777SHong Zhang #define __FUNCT__ "PCFactorSetColumnPivot_Factor" 2527087cfbeSBarry Smith PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol) 25385317021SBarry Smith { 25485317021SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 25585317021SBarry Smith 25685317021SBarry Smith PetscFunctionBegin; 25757622a8eSBarry Smith if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",(double)dtcol); 25885317021SBarry Smith dir->info.dtcol = dtcol; 25985317021SBarry Smith PetscFunctionReturn(0); 26085317021SBarry Smith } 2618ff23777SHong Zhang 2628ff23777SHong Zhang #undef __FUNCT__ 2638ff23777SHong Zhang #define __FUNCT__ "PCSetFromOptions_Factor" 2644416b707SBarry Smith PetscErrorCode PCSetFromOptions_Factor(PetscOptionItems *PetscOptionsObject,PC pc) 2658ff23777SHong Zhang { 2668ff23777SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 2678ff23777SHong Zhang PetscErrorCode ierr; 2688afaa268SBarry Smith PetscBool flg,set; 2698ff23777SHong Zhang char tname[256], solvertype[64]; 270140e18c1SBarry Smith PetscFunctionList ordlist; 271018dd85eSSatish Balay PetscEnum etmp; 2728e37d05fSBarry Smith PetscBool inplace; 2738ff23777SHong Zhang 2748ff23777SHong Zhang PetscFunctionBegin; 2758e37d05fSBarry Smith ierr = PCFactorGetUseInPlace(pc,&inplace);CHKERRQ(ierr); 27692e9c092SBarry Smith ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set);CHKERRQ(ierr); 2778e37d05fSBarry Smith if (set) { 2788e37d05fSBarry Smith ierr = PCFactorSetUseInPlace(pc,flg);CHKERRQ(ierr); 2798ff23777SHong Zhang } 2808afaa268SBarry 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); 2818ff23777SHong Zhang 2828afaa268SBarry 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); 283018dd85eSSatish Balay if (flg) { 2846cc3b3c1SHong Zhang ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr); 285018dd85eSSatish Balay } 286d90ac83dSHong 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); 287d90ac83dSHong Zhang 2888ff23777SHong 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); 2898ff23777SHong 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); 2908ff23777SHong Zhang 2918afaa268SBarry 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); 2928ff23777SHong Zhang if (set) { 2938ff23777SHong Zhang ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 2948ff23777SHong Zhang } 2958ff23777SHong Zhang 2968afaa268SBarry Smith ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 2978afaa268SBarry Smith if (set) { 2988afaa268SBarry Smith ierr = PCFactorSetReuseFill(pc,flg);CHKERRQ(ierr); 2998ff23777SHong Zhang } 3008afaa268SBarry Smith ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 3018afaa268SBarry Smith if (set) { 3028afaa268SBarry Smith ierr = PCFactorSetReuseOrdering(pc,flg);CHKERRQ(ierr); 3038ff23777SHong Zhang } 3048ff23777SHong Zhang 3058ff23777SHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 306a264d7a6SBarry Smith ierr = PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr); 3078ff23777SHong Zhang if (flg) { 3088ff23777SHong Zhang ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 3098ff23777SHong Zhang } 3108ff23777SHong Zhang 3118ff23777SHong Zhang /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 3128ff23777SHong Zhang ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 3138ff23777SHong Zhang if (flg) { 3148ff23777SHong Zhang ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 3158ff23777SHong Zhang } 3168ff23777SHong Zhang PetscFunctionReturn(0); 3178ff23777SHong Zhang } 318914a5d51SHong Zhang 319914a5d51SHong Zhang #undef __FUNCT__ 320914a5d51SHong Zhang #define __FUNCT__ "PCView_Factor" 321914a5d51SHong Zhang PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer) 322914a5d51SHong Zhang { 323914a5d51SHong Zhang PC_Factor *factor = (PC_Factor*)pc->data; 324914a5d51SHong Zhang PetscErrorCode ierr; 325ace3abfcSBarry Smith PetscBool isstring,iascii; 326914a5d51SHong Zhang 327914a5d51SHong Zhang PetscFunctionBegin; 328251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 329251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 330914a5d51SHong Zhang if (iascii) { 331*3d1c1ea0SBarry Smith if (factor->inplace) { 332*3d1c1ea0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr); 333*3d1c1ea0SBarry Smith } else { 334*3d1c1ea0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr); 335*3d1c1ea0SBarry Smith } 336*3d1c1ea0SBarry Smith 337*3d1c1ea0SBarry Smith if (factor->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 338*3d1c1ea0SBarry Smith if (factor->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 339879e8a4dSBarry Smith if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) { 340914a5d51SHong Zhang if (factor->info.dt > 0) { 34157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," drop tolerance %g\n",(double)factor->info.dt);CHKERRQ(ierr); 342914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr); 34357622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," column permutation tolerance %g\n",(double)factor->info.dtcol);CHKERRQ(ierr); 344914a5d51SHong Zhang } else if (factor->info.levels == 1) { 345914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 346914a5d51SHong Zhang } else { 347914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr); 348914a5d51SHong Zhang } 349914a5d51SHong Zhang } 350914a5d51SHong Zhang 35157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerance for zero pivot %g\n",(double)factor->info.zeropivot);CHKERRQ(ierr); 3525e9742b9SJed Brown if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */ 3535e9742b9SJed Brown ierr = PetscViewerASCIIPrintf(viewer," using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);CHKERRQ(ierr); 354914a5d51SHong Zhang } 355914a5d51SHong Zhang 356914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr); 357914a5d51SHong Zhang 358914a5d51SHong Zhang if (factor->fact) { 3596335c336SBarry Smith MatInfo info; 3606335c336SBarry Smith ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 36157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);CHKERRQ(ierr); 362914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows:\n");CHKERRQ(ierr); 363914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 364914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 365914a5d51SHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 366914a5d51SHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 367914a5d51SHong Zhang ierr = MatView(factor->fact,viewer);CHKERRQ(ierr); 368914a5d51SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 369914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 370914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 371914a5d51SHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 372914a5d51SHong Zhang } 373914a5d51SHong Zhang 374914a5d51SHong Zhang } else if (isstring) { 3756335c336SBarry Smith MatFactorType t; 3766335c336SBarry Smith ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr); 3776335c336SBarry Smith if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) { 378914a5d51SHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 379914a5d51SHong Zhang } 380914a5d51SHong Zhang } 381914a5d51SHong Zhang PetscFunctionReturn(0); 382914a5d51SHong Zhang } 383