1dba47a55SKris Buschelman #define PETSCKSP_DLL 29b54502bSHong Zhang 39b54502bSHong Zhang /* 49b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 59b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 69b54502bSHong Zhang a direct solver. 79b54502bSHong Zhang */ 86356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 99b54502bSHong Zhang 109b54502bSHong Zhang typedef struct { 119b54502bSHong Zhang Mat fact; /* factored matrix */ 129b54502bSHong Zhang PetscReal actualfill; /* actual fill in factor */ 139b54502bSHong Zhang PetscTruth inplace; /* flag indicating in-place factorization */ 149b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 159b54502bSHong Zhang MatOrderingType ordering; /* matrix ordering */ 169b54502bSHong Zhang PetscTruth reuseordering; /* reuses previous reordering computed */ 179b54502bSHong Zhang PetscTruth reusefill; /* reuse fill from previous Cholesky */ 189b54502bSHong Zhang MatFactorInfo info; 19*5c9eb25fSBarry Smith MatSolverType solvertype; 209b54502bSHong Zhang } PC_Cholesky; 219b54502bSHong Zhang 229b54502bSHong Zhang EXTERN_C_BEGIN 239b54502bSHong Zhang #undef __FUNCT__ 24afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_Cholesky" 25dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Cholesky(PC pc,PetscReal z) 26afaefe49SHong Zhang { 27afaefe49SHong Zhang PC_Cholesky *ch; 28afaefe49SHong Zhang 29afaefe49SHong Zhang PetscFunctionBegin; 30afaefe49SHong Zhang ch = (PC_Cholesky*)pc->data; 31afaefe49SHong Zhang ch->info.zeropivot = z; 32afaefe49SHong Zhang PetscFunctionReturn(0); 33afaefe49SHong Zhang } 34afaefe49SHong Zhang EXTERN_C_END 35afaefe49SHong Zhang 36afaefe49SHong Zhang EXTERN_C_BEGIN 37afaefe49SHong Zhang #undef __FUNCT__ 38afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_Cholesky" 39dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Cholesky(PC pc,PetscReal shift) 40afaefe49SHong Zhang { 41afaefe49SHong Zhang PC_Cholesky *dir; 42afaefe49SHong Zhang 43afaefe49SHong Zhang PetscFunctionBegin; 44afaefe49SHong Zhang dir = (PC_Cholesky*)pc->data; 45afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 46afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 47afaefe49SHong Zhang } else { 48afaefe49SHong Zhang dir->info.shiftnz = shift; 49afaefe49SHong Zhang } 50afaefe49SHong Zhang PetscFunctionReturn(0); 51afaefe49SHong Zhang } 52afaefe49SHong Zhang EXTERN_C_END 53afaefe49SHong Zhang 54afaefe49SHong Zhang EXTERN_C_BEGIN 55afaefe49SHong Zhang #undef __FUNCT__ 56afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_Cholesky" 57dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Cholesky(PC pc,PetscTruth shift) 58afaefe49SHong Zhang { 59afaefe49SHong Zhang PC_Cholesky *dir; 60afaefe49SHong Zhang 61afaefe49SHong Zhang PetscFunctionBegin; 62afaefe49SHong Zhang dir = (PC_Cholesky*)pc->data; 63fbf22428SSatish Balay if (shift) { 64fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 65fbf22428SSatish Balay dir->info.shiftpd = 1.0; 66fbf22428SSatish Balay } else { 67fbf22428SSatish Balay dir->info.shiftpd = 0.0; 68fbf22428SSatish Balay } 69afaefe49SHong Zhang PetscFunctionReturn(0); 70afaefe49SHong Zhang } 71afaefe49SHong Zhang EXTERN_C_END 72afaefe49SHong Zhang 73afaefe49SHong Zhang EXTERN_C_BEGIN 74afaefe49SHong Zhang #undef __FUNCT__ 752401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 762401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag) 779b54502bSHong Zhang { 789b54502bSHong Zhang PC_Cholesky *lu; 799b54502bSHong Zhang 809b54502bSHong Zhang PetscFunctionBegin; 819b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 829b54502bSHong Zhang lu->reuseordering = flag; 839b54502bSHong Zhang PetscFunctionReturn(0); 849b54502bSHong Zhang } 859b54502bSHong Zhang EXTERN_C_END 869b54502bSHong Zhang 879b54502bSHong Zhang EXTERN_C_BEGIN 889b54502bSHong Zhang #undef __FUNCT__ 892401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 902401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag) 919b54502bSHong Zhang { 929b54502bSHong Zhang PC_Cholesky *lu; 939b54502bSHong Zhang 949b54502bSHong Zhang PetscFunctionBegin; 959b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 969b54502bSHong Zhang lu->reusefill = flag; 979b54502bSHong Zhang PetscFunctionReturn(0); 989b54502bSHong Zhang } 999b54502bSHong Zhang EXTERN_C_END 1009b54502bSHong Zhang 1019b54502bSHong Zhang #undef __FUNCT__ 1029b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 1039b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 1049b54502bSHong Zhang { 1059b54502bSHong Zhang PC_Cholesky *lu = (PC_Cholesky*)pc->data; 1069b54502bSHong Zhang PetscErrorCode ierr; 1079b54502bSHong Zhang PetscTruth flg; 108*5c9eb25fSBarry Smith char tname[256], solvertype[64]; 1099b54502bSHong Zhang PetscFList ordlist; 1109b54502bSHong Zhang 1119b54502bSHong Zhang PetscFunctionBegin; 1129b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1139b54502bSHong Zhang ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 1142401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 1159b54502bSHong Zhang if (flg) { 1162401956bSBarry Smith ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 1179b54502bSHong Zhang } 11855ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 1199b54502bSHong Zhang 1202401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 1219b54502bSHong Zhang if (flg) { 1222401956bSBarry Smith ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1239b54502bSHong Zhang } 1242401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 1259b54502bSHong Zhang if (flg) { 1262401956bSBarry Smith ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1279b54502bSHong Zhang } 1289b54502bSHong Zhang 1299b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 130e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 1319b54502bSHong Zhang if (flg) { 132e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 1339b54502bSHong Zhang } 134*5c9eb25fSBarry Smith 135*5c9eb25fSBarry Smith /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 136*5c9eb25fSBarry Smith ierr = PetscOptionsString("-pc_mat_solver_type","Specific Cholesky solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 137*5c9eb25fSBarry Smith if (flg) { 138*5c9eb25fSBarry Smith ierr = PetscStrallocpy(solvertype,&lu->solvertype);CHKERRQ(ierr); 139*5c9eb25fSBarry Smith } 140*5c9eb25fSBarry Smith 1419f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 1429b54502bSHong Zhang if (flg) { 143afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1449b54502bSHong Zhang } 1459f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 1469f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 1479b54502bSHong Zhang if (flg) { 148afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 1499b54502bSHong Zhang } 150ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 1519b54502bSHong Zhang 1529b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1539b54502bSHong Zhang PetscFunctionReturn(0); 1549b54502bSHong Zhang } 1559b54502bSHong Zhang 1569b54502bSHong Zhang #undef __FUNCT__ 1579b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 1589b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 1599b54502bSHong Zhang { 1609b54502bSHong Zhang PC_Cholesky *lu = (PC_Cholesky*)pc->data; 1619b54502bSHong Zhang PetscErrorCode ierr; 1629b54502bSHong Zhang PetscTruth iascii,isstring; 1639b54502bSHong Zhang 1649b54502bSHong Zhang PetscFunctionBegin; 1659b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1669b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1679b54502bSHong Zhang if (iascii) { 1689b54502bSHong Zhang 1699b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr);} 1709b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr);} 1719b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 1729b54502bSHong Zhang if (lu->fact) { 173f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cholesky: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 174f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 175f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 176f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 177f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 178f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 1799b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 1809b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 181f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 182f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 183f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1849b54502bSHong Zhang } 1859b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 1869b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1879b54502bSHong Zhang } else if (isstring) { 1889b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 1899b54502bSHong Zhang } else { 1902401956bSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name); 1919b54502bSHong Zhang } 1929b54502bSHong Zhang PetscFunctionReturn(0); 1939b54502bSHong Zhang } 1949b54502bSHong Zhang 1959b54502bSHong Zhang #undef __FUNCT__ 196a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_Cholesky" 197a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_Cholesky(PC pc,Mat *mat) 1989b54502bSHong Zhang { 1999b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2009b54502bSHong Zhang 2019b54502bSHong Zhang PetscFunctionBegin; 2029b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 2039b54502bSHong Zhang *mat = dir->fact; 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang #undef __FUNCT__ 2089b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 2099b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 2109b54502bSHong Zhang { 2119b54502bSHong Zhang PetscErrorCode ierr; 2129b54502bSHong Zhang PetscTruth flg; 2139b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2149b54502bSHong Zhang 2159b54502bSHong Zhang PetscFunctionBegin; 2169b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 2179b54502bSHong Zhang 2189b54502bSHong Zhang if (dir->inplace) { 2199b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 2209b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 2219b54502bSHong Zhang dir->row = 0; 2229b54502bSHong Zhang } 2239b54502bSHong Zhang if (dir->col) { 2249b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2259b54502bSHong Zhang dir->col = 0; 2269b54502bSHong Zhang } 2279b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2289b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 2299b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2309b54502bSHong Zhang dir->col=0; 2319b54502bSHong Zhang } 232efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 2339b54502bSHong Zhang ierr = MatCholeskyFactor(pc->pmat,dir->row,&dir->info);CHKERRQ(ierr); 2349b54502bSHong Zhang dir->fact = pc->pmat; 2359b54502bSHong Zhang } else { 2369b54502bSHong Zhang MatInfo info; 2379b54502bSHong Zhang if (!pc->setupcalled) { 2389b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2399bfd6278SHong Zhang /* check if dir->row == dir->col */ 2409bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 2419bfd6278SHong Zhang if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 2429bfd6278SHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 2439b54502bSHong Zhang dir->col=0; 2449bfd6278SHong Zhang 2457adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 2469b54502bSHong Zhang if (flg) { 2479b54502bSHong Zhang PetscReal tol = 1.e-10; 2487adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 2499b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 2509b54502bSHong Zhang } 251efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 252*5c9eb25fSBarry Smith ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_CHOLESKY,&dir->fact);CHKERRQ(ierr); 2539b54502bSHong Zhang ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr); 2549b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2559b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 25652e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2579b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2589b54502bSHong Zhang if (!dir->reuseordering) { 2599b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 2609b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 2619b54502bSHong Zhang dir->row = 0; 2629b54502bSHong Zhang } 2639b54502bSHong Zhang if (dir->col) { 2649b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2659b54502bSHong Zhang dir->col =0; 2669b54502bSHong Zhang } 2679b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2689b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 2699b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2709b54502bSHong Zhang dir->col=0; 2719b54502bSHong Zhang } 2727adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 2739b54502bSHong Zhang if (flg) { 2749b54502bSHong Zhang PetscReal tol = 1.e-10; 2757adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 2769b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 2779b54502bSHong Zhang } 278efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 2799b54502bSHong Zhang } 2809b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 281*5c9eb25fSBarry Smith ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_CHOLESKY,&dir->fact);CHKERRQ(ierr); 2829b54502bSHong Zhang ierr = MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);CHKERRQ(ierr); 2839b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2849b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 28552e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2869b54502bSHong Zhang } 2879b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 2889b54502bSHong Zhang } 2899b54502bSHong Zhang PetscFunctionReturn(0); 2909b54502bSHong Zhang } 2919b54502bSHong Zhang 2929b54502bSHong Zhang #undef __FUNCT__ 2939b54502bSHong Zhang #define __FUNCT__ "PCDestroy_Cholesky" 2949b54502bSHong Zhang static PetscErrorCode PCDestroy_Cholesky(PC pc) 2959b54502bSHong Zhang { 2969b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2979b54502bSHong Zhang PetscErrorCode ierr; 2989b54502bSHong Zhang 2999b54502bSHong Zhang PetscFunctionBegin; 3009b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 3019b54502bSHong Zhang if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 3029b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 3039b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 304*5c9eb25fSBarry Smith ierr = PetscStrfree(dir->solvertype);CHKERRQ(ierr); 3059b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 3069b54502bSHong Zhang PetscFunctionReturn(0); 3079b54502bSHong Zhang } 3089b54502bSHong Zhang 3099b54502bSHong Zhang #undef __FUNCT__ 3109b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 3119b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 3129b54502bSHong Zhang { 3139b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 3149b54502bSHong Zhang PetscErrorCode ierr; 3159b54502bSHong Zhang 3169b54502bSHong Zhang PetscFunctionBegin; 3179b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 3189b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 3199b54502bSHong Zhang PetscFunctionReturn(0); 3209b54502bSHong Zhang } 3219b54502bSHong Zhang 3229b54502bSHong Zhang #undef __FUNCT__ 3239b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 3249b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 3259b54502bSHong Zhang { 3269b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 3279b54502bSHong Zhang PetscErrorCode ierr; 3289b54502bSHong Zhang 3299b54502bSHong Zhang PetscFunctionBegin; 3309b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 3319b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 3329b54502bSHong Zhang PetscFunctionReturn(0); 3339b54502bSHong Zhang } 3349b54502bSHong Zhang 3359b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3369b54502bSHong Zhang 3379b54502bSHong Zhang EXTERN_C_BEGIN 3389b54502bSHong Zhang #undef __FUNCT__ 33955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_Cholesky" 34055ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill) 3419b54502bSHong Zhang { 3429b54502bSHong Zhang PC_Cholesky *dir; 3439b54502bSHong Zhang 3449b54502bSHong Zhang PetscFunctionBegin; 3459b54502bSHong Zhang dir = (PC_Cholesky*)pc->data; 3469b54502bSHong Zhang dir->info.fill = fill; 3479b54502bSHong Zhang PetscFunctionReturn(0); 3489b54502bSHong Zhang } 3499b54502bSHong Zhang EXTERN_C_END 3509b54502bSHong Zhang 3519b54502bSHong Zhang EXTERN_C_BEGIN 3529b54502bSHong Zhang #undef __FUNCT__ 3532401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 3542401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc) 3559b54502bSHong Zhang { 3569b54502bSHong Zhang PC_Cholesky *dir; 3579b54502bSHong Zhang 3589b54502bSHong Zhang PetscFunctionBegin; 3599b54502bSHong Zhang dir = (PC_Cholesky*)pc->data; 3609b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3619b54502bSHong Zhang PetscFunctionReturn(0); 3629b54502bSHong Zhang } 3639b54502bSHong Zhang EXTERN_C_END 3649b54502bSHong Zhang 3659b54502bSHong Zhang EXTERN_C_BEGIN 3669b54502bSHong Zhang #undef __FUNCT__ 367e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_Cholesky" 368e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_Cholesky(PC pc,MatOrderingType ordering) 3699b54502bSHong Zhang { 3709b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 3719b54502bSHong Zhang PetscErrorCode ierr; 3729b54502bSHong Zhang 3739b54502bSHong Zhang PetscFunctionBegin; 3749b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3759b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 3769b54502bSHong Zhang PetscFunctionReturn(0); 3779b54502bSHong Zhang } 3789b54502bSHong Zhang EXTERN_C_END 3799b54502bSHong Zhang 3809b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3819b54502bSHong Zhang 3829b54502bSHong Zhang #undef __FUNCT__ 3832401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 3849b54502bSHong Zhang /*@ 3852401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 3869b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 3879b54502bSHong Zhang following factors. 3889b54502bSHong Zhang 3899b54502bSHong Zhang Collective on PC 3909b54502bSHong Zhang 3919b54502bSHong Zhang Input Parameters: 3929b54502bSHong Zhang + pc - the preconditioner context 3939b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 3949b54502bSHong Zhang 3959b54502bSHong Zhang Options Database Key: 3962401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 3979b54502bSHong Zhang 3989b54502bSHong Zhang Level: intermediate 3999b54502bSHong Zhang 4009b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 4019b54502bSHong Zhang 4022401956bSBarry Smith .seealso: PCFactorSetReuseFill() 4039b54502bSHong Zhang @*/ 4042401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag) 4059b54502bSHong Zhang { 4069b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4079b54502bSHong Zhang 4089b54502bSHong Zhang PetscFunctionBegin; 4099b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4102401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 4119b54502bSHong Zhang if (f) { 4129b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 4139b54502bSHong Zhang } 4149b54502bSHong Zhang PetscFunctionReturn(0); 4159b54502bSHong Zhang } 4169b54502bSHong Zhang 4179b54502bSHong Zhang #undef __FUNCT__ 4182401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill" 4199b54502bSHong Zhang /*@ 4202401956bSBarry Smith PCFactorSetReuseFill - When matrices with same different nonzero structure are factored, 4212401956bSBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 4229b54502bSHong Zhang 4239b54502bSHong Zhang Collective on PC 4249b54502bSHong Zhang 4259b54502bSHong Zhang Input Parameters: 4269b54502bSHong Zhang + pc - the preconditioner context 4279b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 4289b54502bSHong Zhang 4299b54502bSHong Zhang Options Database Key: 4302401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 4319b54502bSHong Zhang 4329b54502bSHong Zhang Level: intermediate 4339b54502bSHong Zhang 4349b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, Cholesky 4359b54502bSHong Zhang 4362401956bSBarry Smith .seealso: PCFactorSetReuseOrdering() 4379b54502bSHong Zhang @*/ 4382401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag) 4399b54502bSHong Zhang { 4409b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4419b54502bSHong Zhang 4429b54502bSHong Zhang PetscFunctionBegin; 4439b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,2); 4442401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 4459b54502bSHong Zhang if (f) { 4469b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 4479b54502bSHong Zhang } 4489b54502bSHong Zhang PetscFunctionReturn(0); 4499b54502bSHong Zhang } 4509b54502bSHong Zhang 4519b54502bSHong Zhang /*MC 45296fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 4539b54502bSHong Zhang 4549b54502bSHong Zhang Options Database Keys: 4552401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 4562401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 45755ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 4582401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 4592401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 460f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 461f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 462f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 4639b54502bSHong Zhang 4649b54502bSHong Zhang Notes: Not all options work for all matrix formats 4659b54502bSHong Zhang 4669b54502bSHong Zhang Level: beginner 4679b54502bSHong Zhang 4689b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 4699b54502bSHong Zhang 4709b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 4719b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 4729b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 4739b54502bSHong Zhang 4749b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 475a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 47655ba2a51SBarry Smith PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 477e5a9bf91SBarry Smith PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 4789b54502bSHong Zhang 4799b54502bSHong Zhang M*/ 4809b54502bSHong Zhang 4819b54502bSHong Zhang EXTERN_C_BEGIN 4829b54502bSHong Zhang #undef __FUNCT__ 4839b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 484dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc) 4859b54502bSHong Zhang { 4869b54502bSHong Zhang PetscErrorCode ierr; 4879b54502bSHong Zhang PC_Cholesky *dir; 4889b54502bSHong Zhang 4899b54502bSHong Zhang PetscFunctionBegin; 49038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 4919b54502bSHong Zhang 4929b54502bSHong Zhang dir->fact = 0; 4939b54502bSHong Zhang dir->inplace = PETSC_FALSE; 4949b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 4959b54502bSHong Zhang dir->info.fill = 5.0; 4960a29876aSHong Zhang dir->info.shiftnz = 0.0; 497fbf22428SSatish Balay dir->info.shiftpd = 0.0; /* false */ 4989b54502bSHong Zhang dir->info.shift_fraction = 0.0; 4999b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 5009b54502bSHong Zhang dir->col = 0; 5019b54502bSHong Zhang dir->row = 0; 5029b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 503*5c9eb25fSBarry Smith ierr = PetscStrallocpy("petsc",&dir->solvertype);CHKERRQ(ierr); 5049b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 5059b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 5069b54502bSHong Zhang pc->data = (void*)dir; 5079b54502bSHong Zhang 5089b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 5099b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 5109b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 5119b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 5129b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 5139b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 5149b54502bSHong Zhang pc->ops->applyrichardson = 0; 515a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Cholesky; 5169b54502bSHong Zhang 517afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky", 518afaefe49SHong Zhang PCFactorSetZeroPivot_Cholesky);CHKERRQ(ierr); 519afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky", 520afaefe49SHong Zhang PCFactorSetShiftNonzero_Cholesky);CHKERRQ(ierr); 521afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky", 522afaefe49SHong Zhang PCFactorSetShiftPd_Cholesky);CHKERRQ(ierr); 523afaefe49SHong Zhang 52455ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky", 52555ba2a51SBarry Smith PCFactorSetFill_Cholesky);CHKERRQ(ierr); 5262401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky", 5272401956bSBarry Smith PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 528e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Cholesky", 529e5a9bf91SBarry Smith PCFactorSetMatOrderingType_Cholesky);CHKERRQ(ierr); 5302401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky", 5312401956bSBarry Smith PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 5322401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky", 5332401956bSBarry Smith PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 5349b54502bSHong Zhang PetscFunctionReturn(0); 5359b54502bSHong Zhang } 5369b54502bSHong Zhang EXTERN_C_END 537