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 */ 8*7c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/factor.h" /*I "petscpc.h" I*/ 99b54502bSHong Zhang 109b54502bSHong Zhang typedef struct { 11075768bcSBarry Smith PC_Factor hdr; 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 PetscTruth reuseordering; /* reuses previous reordering computed */ 169b54502bSHong Zhang PetscTruth reusefill; /* reuse fill from previous Cholesky */ 179b54502bSHong Zhang } PC_Cholesky; 189b54502bSHong Zhang 199b54502bSHong Zhang EXTERN_C_BEGIN 209b54502bSHong Zhang #undef __FUNCT__ 212401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_Cholesky" 222401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag) 239b54502bSHong Zhang { 249b54502bSHong Zhang PC_Cholesky *lu; 259b54502bSHong Zhang 269b54502bSHong Zhang PetscFunctionBegin; 279b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 289b54502bSHong Zhang lu->reuseordering = flag; 299b54502bSHong Zhang PetscFunctionReturn(0); 309b54502bSHong Zhang } 319b54502bSHong Zhang EXTERN_C_END 329b54502bSHong Zhang 339b54502bSHong Zhang EXTERN_C_BEGIN 349b54502bSHong Zhang #undef __FUNCT__ 352401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_Cholesky" 362401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag) 379b54502bSHong Zhang { 389b54502bSHong Zhang PC_Cholesky *lu; 399b54502bSHong Zhang 409b54502bSHong Zhang PetscFunctionBegin; 419b54502bSHong Zhang lu = (PC_Cholesky*)pc->data; 429b54502bSHong Zhang lu->reusefill = flag; 439b54502bSHong Zhang PetscFunctionReturn(0); 449b54502bSHong Zhang } 459b54502bSHong Zhang EXTERN_C_END 469b54502bSHong Zhang 479b54502bSHong Zhang #undef __FUNCT__ 489b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_Cholesky" 499b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_Cholesky(PC pc) 509b54502bSHong Zhang { 519b54502bSHong Zhang PC_Cholesky *lu = (PC_Cholesky*)pc->data; 529b54502bSHong Zhang PetscErrorCode ierr; 539b54502bSHong Zhang PetscTruth flg; 545c9eb25fSBarry Smith char tname[256], solvertype[64]; 559b54502bSHong Zhang PetscFList ordlist; 569b54502bSHong Zhang 579b54502bSHong Zhang PetscFunctionBegin; 589b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 599b54502bSHong Zhang ierr = PetscOptionsHead("Cholesky options");CHKERRQ(ierr); 602401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 619b54502bSHong Zhang if (flg) { 622401956bSBarry Smith ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 639b54502bSHong Zhang } 64075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr); 659b54502bSHong Zhang 662401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 679b54502bSHong Zhang if (flg) { 682401956bSBarry Smith ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 699b54502bSHong Zhang } 702401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 719b54502bSHong Zhang if (flg) { 722401956bSBarry Smith ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 739b54502bSHong Zhang } 749b54502bSHong Zhang 759b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 76075768bcSBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr); 779b54502bSHong Zhang if (flg) { 78e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 799b54502bSHong Zhang } 805c9eb25fSBarry Smith 815c9eb25fSBarry Smith /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 82075768bcSBarry Smith ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific Cholesky solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 835c9eb25fSBarry Smith if (flg) { 84c7393fdbSBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 855c9eb25fSBarry Smith } 865c9eb25fSBarry Smith 879f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 889b54502bSHong Zhang if (flg) { 89afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 909b54502bSHong Zhang } 91075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);CHKERRQ(ierr); 929f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 939b54502bSHong Zhang if (flg) { 94afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 959b54502bSHong Zhang } 96075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)lu)->info.zeropivot,&((PC_Factor*)lu)->info.zeropivot,0);CHKERRQ(ierr); 979b54502bSHong Zhang 989b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 999b54502bSHong Zhang PetscFunctionReturn(0); 1009b54502bSHong Zhang } 1019b54502bSHong Zhang 1029b54502bSHong Zhang #undef __FUNCT__ 1039b54502bSHong Zhang #define __FUNCT__ "PCView_Cholesky" 1049b54502bSHong Zhang static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer) 1059b54502bSHong Zhang { 1069b54502bSHong Zhang PC_Cholesky *lu = (PC_Cholesky*)pc->data; 1079b54502bSHong Zhang PetscErrorCode ierr; 1089b54502bSHong Zhang PetscTruth iascii,isstring; 1099b54502bSHong Zhang 1109b54502bSHong Zhang PetscFunctionBegin; 1119b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1129b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1139b54502bSHong Zhang if (iascii) { 1149b54502bSHong Zhang 1159b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");CHKERRQ(ierr);} 1169b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");CHKERRQ(ierr);} 117075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr); 118075768bcSBarry Smith if (((PC_Factor*)lu)->fact) { 119f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cholesky: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 120f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 121f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 122f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 123f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 124f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 125075768bcSBarry Smith ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr); 1269b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 127f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 128f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 129f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1309b54502bSHong Zhang } 1319b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 1329b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1339b54502bSHong Zhang } else if (isstring) { 134075768bcSBarry Smith ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 1359b54502bSHong Zhang } else { 1362401956bSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name); 1379b54502bSHong Zhang } 1389b54502bSHong Zhang PetscFunctionReturn(0); 1399b54502bSHong Zhang } 1409b54502bSHong Zhang 1419b54502bSHong Zhang 1429b54502bSHong Zhang #undef __FUNCT__ 1439b54502bSHong Zhang #define __FUNCT__ "PCSetUp_Cholesky" 1449b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 1459b54502bSHong Zhang { 1469b54502bSHong Zhang PetscErrorCode ierr; 1479b54502bSHong Zhang PetscTruth flg; 1489b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1499b54502bSHong Zhang 1509b54502bSHong Zhang PetscFunctionBegin; 151075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1529b54502bSHong Zhang 1539b54502bSHong Zhang if (dir->inplace) { 1549b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 1559b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 1569b54502bSHong Zhang dir->row = 0; 1579b54502bSHong Zhang } 1589b54502bSHong Zhang if (dir->col) { 1599b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1609b54502bSHong Zhang dir->col = 0; 1619b54502bSHong Zhang } 162075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1639b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 1649b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 1659b54502bSHong Zhang dir->col=0; 1669b54502bSHong Zhang } 167efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 168075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 169075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1709b54502bSHong Zhang } else { 1719b54502bSHong Zhang MatInfo info; 1729b54502bSHong Zhang if (!pc->setupcalled) { 173075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1749bfd6278SHong Zhang /* check if dir->row == dir->col */ 1759bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 1769bfd6278SHong Zhang if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 1779bfd6278SHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 1789b54502bSHong Zhang dir->col=0; 1799bfd6278SHong Zhang 1807adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 1819b54502bSHong Zhang if (flg) { 1829b54502bSHong Zhang PetscReal tol = 1.e-10; 1837adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 1849b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 1859b54502bSHong Zhang } 186efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 187075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 188075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 189075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1909b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 191075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1929b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1939b54502bSHong Zhang if (!dir->reuseordering) { 1949b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 1959b54502bSHong Zhang ierr = ISDestroy(dir->row);CHKERRQ(ierr); 1969b54502bSHong Zhang dir->row = 0; 1979b54502bSHong Zhang } 1989b54502bSHong Zhang if (dir->col) { 1999b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2009b54502bSHong Zhang dir->col =0; 2019b54502bSHong Zhang } 202075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2039b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 2049b54502bSHong Zhang ierr = ISDestroy(dir->col);CHKERRQ(ierr); 2059b54502bSHong Zhang dir->col=0; 2069b54502bSHong Zhang } 2077adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);CHKERRQ(ierr); 2089b54502bSHong Zhang if (flg) { 2099b54502bSHong Zhang PetscReal tol = 1.e-10; 2107adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);CHKERRQ(ierr); 2119b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 2129b54502bSHong Zhang } 213efee365bSSatish Balay if (dir->row) {ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr);} 2149b54502bSHong Zhang } 215075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 216075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 217075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 218075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2199b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 220075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 2219b54502bSHong Zhang } 222075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 2239b54502bSHong Zhang } 2249b54502bSHong Zhang PetscFunctionReturn(0); 2259b54502bSHong Zhang } 2269b54502bSHong Zhang 2279b54502bSHong Zhang #undef __FUNCT__ 2289b54502bSHong Zhang #define __FUNCT__ "PCDestroy_Cholesky" 2299b54502bSHong Zhang static PetscErrorCode PCDestroy_Cholesky(PC pc) 2309b54502bSHong Zhang { 2319b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2329b54502bSHong Zhang PetscErrorCode ierr; 2339b54502bSHong Zhang 2349b54502bSHong Zhang PetscFunctionBegin; 235075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 2369b54502bSHong Zhang if (dir->row) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2379b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 238075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 239075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 2409b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 2419b54502bSHong Zhang PetscFunctionReturn(0); 2429b54502bSHong Zhang } 2439b54502bSHong Zhang 2449b54502bSHong Zhang #undef __FUNCT__ 2459b54502bSHong Zhang #define __FUNCT__ "PCApply_Cholesky" 2469b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 2479b54502bSHong Zhang { 2489b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2499b54502bSHong Zhang PetscErrorCode ierr; 2509b54502bSHong Zhang 2519b54502bSHong Zhang PetscFunctionBegin; 2529b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 253075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2549b54502bSHong Zhang PetscFunctionReturn(0); 2559b54502bSHong Zhang } 2569b54502bSHong Zhang 2579b54502bSHong Zhang #undef __FUNCT__ 2589b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_Cholesky" 2599b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 2609b54502bSHong Zhang { 2619b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2629b54502bSHong Zhang PetscErrorCode ierr; 2639b54502bSHong Zhang 2649b54502bSHong Zhang PetscFunctionBegin; 2659b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 266075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2679b54502bSHong Zhang PetscFunctionReturn(0); 2689b54502bSHong Zhang } 2699b54502bSHong Zhang 2709b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2719b54502bSHong Zhang 2729b54502bSHong Zhang EXTERN_C_BEGIN 2739b54502bSHong Zhang #undef __FUNCT__ 2742401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_Cholesky" 2752401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc) 2769b54502bSHong Zhang { 2779b54502bSHong Zhang PC_Cholesky *dir; 2789b54502bSHong Zhang 2799b54502bSHong Zhang PetscFunctionBegin; 2809b54502bSHong Zhang dir = (PC_Cholesky*)pc->data; 2819b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2829b54502bSHong Zhang PetscFunctionReturn(0); 2839b54502bSHong Zhang } 2849b54502bSHong Zhang EXTERN_C_END 2859b54502bSHong Zhang 2869b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2879b54502bSHong Zhang 2889b54502bSHong Zhang #undef __FUNCT__ 2892401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering" 2909b54502bSHong Zhang /*@ 2912401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2929b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2939b54502bSHong Zhang following factors. 2949b54502bSHong Zhang 2959b54502bSHong Zhang Collective on PC 2969b54502bSHong Zhang 2979b54502bSHong Zhang Input Parameters: 2989b54502bSHong Zhang + pc - the preconditioner context 2999b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 3009b54502bSHong Zhang 3019b54502bSHong Zhang Options Database Key: 3022401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 3039b54502bSHong Zhang 3049b54502bSHong Zhang Level: intermediate 3059b54502bSHong Zhang 3069b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 3079b54502bSHong Zhang 3082401956bSBarry Smith .seealso: PCFactorSetReuseFill() 3099b54502bSHong Zhang @*/ 3102401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag) 3119b54502bSHong Zhang { 3129b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 3139b54502bSHong Zhang 3149b54502bSHong Zhang PetscFunctionBegin; 3159b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3162401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 3179b54502bSHong Zhang if (f) { 3189b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 3199b54502bSHong Zhang } 3209b54502bSHong Zhang PetscFunctionReturn(0); 3219b54502bSHong Zhang } 3229b54502bSHong Zhang 3239b54502bSHong Zhang 3249b54502bSHong Zhang /*MC 32596fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 3269b54502bSHong Zhang 3279b54502bSHong Zhang Options Database Keys: 3282401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 329c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 3302401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 33155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 3322401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 3332401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 334f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 335f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 336f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 3379b54502bSHong Zhang 3389b54502bSHong Zhang Notes: Not all options work for all matrix formats 3399b54502bSHong Zhang 3409b54502bSHong Zhang Level: beginner 3419b54502bSHong Zhang 3429b54502bSHong Zhang Concepts: Cholesky factorization, direct solver 3439b54502bSHong Zhang 3449b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3459b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3469b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3479b54502bSHong Zhang 3489b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 349a4fd02acSBarry Smith PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 35055ba2a51SBarry Smith PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 351e5a9bf91SBarry Smith PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 3529b54502bSHong Zhang 3539b54502bSHong Zhang M*/ 3549b54502bSHong Zhang 3559b54502bSHong Zhang EXTERN_C_BEGIN 3569b54502bSHong Zhang #undef __FUNCT__ 3579b54502bSHong Zhang #define __FUNCT__ "PCCreate_Cholesky" 358dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc) 3599b54502bSHong Zhang { 3609b54502bSHong Zhang PetscErrorCode ierr; 3619b54502bSHong Zhang PC_Cholesky *dir; 3629b54502bSHong Zhang 3639b54502bSHong Zhang PetscFunctionBegin; 36438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Cholesky,&dir);CHKERRQ(ierr); 3659b54502bSHong Zhang 366075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3679b54502bSHong Zhang dir->inplace = PETSC_FALSE; 368075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 369075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 370075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = 0.0; 371075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 0.0; /* false */ 372075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3739b54502bSHong Zhang dir->col = 0; 3749b54502bSHong Zhang dir->row = 0; 375075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 376075768bcSBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3779b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3789b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3799b54502bSHong Zhang pc->data = (void*)dir; 3809b54502bSHong Zhang 3819b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 3829b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3839b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3849b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3859b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 3869b54502bSHong Zhang pc->ops->view = PCView_Cholesky; 3879b54502bSHong Zhang pc->ops->applyrichardson = 0; 38885317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3899b54502bSHong Zhang 39085317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 39185317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 3927112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3937112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 39485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 39585317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 39685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor", 39785317021SBarry Smith PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr); 39885317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor", 39985317021SBarry Smith PCFactorSetShiftPd_Factor);CHKERRQ(ierr); 400afaefe49SHong Zhang 40185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 40285317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 4032401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky", 4042401956bSBarry Smith PCFactorSetUseInPlace_Cholesky);CHKERRQ(ierr); 40585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 40685317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 4072401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky", 4082401956bSBarry Smith PCFactorSetReuseOrdering_Cholesky);CHKERRQ(ierr); 4092401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky", 4102401956bSBarry Smith PCFactorSetReuseFill_Cholesky);CHKERRQ(ierr); 4119b54502bSHong Zhang PetscFunctionReturn(0); 4129b54502bSHong Zhang } 4139b54502bSHong Zhang EXTERN_C_END 414