1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 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 */ 8ee45ca4aSHong Zhang 96356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 10ee45ca4aSHong Zhang #include "src/ksp/pc/impls/factor/lu/lu.h" 119b54502bSHong Zhang 129b54502bSHong Zhang EXTERN_C_BEGIN 139b54502bSHong Zhang #undef __FUNCT__ 14c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_LU" 15c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_LU(PC pc,const MatSolverPackage stype) 16c62fe420SBarry Smith { 17c62fe420SBarry Smith PetscErrorCode ierr; 18c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 19c62fe420SBarry Smith 20c62fe420SBarry Smith PetscFunctionBegin; 21c62fe420SBarry Smith ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr); 22c62fe420SBarry Smith PetscFunctionReturn(0); 23c62fe420SBarry Smith } 24c62fe420SBarry Smith EXTERN_C_END 25c62fe420SBarry Smith 26c62fe420SBarry Smith EXTERN_C_BEGIN 27c62fe420SBarry Smith #undef __FUNCT__ 28afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU" 29dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z) 30afaefe49SHong Zhang { 31c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 32afaefe49SHong Zhang 33afaefe49SHong Zhang PetscFunctionBegin; 34afaefe49SHong Zhang lu->info.zeropivot = z; 35afaefe49SHong Zhang PetscFunctionReturn(0); 36afaefe49SHong Zhang } 37afaefe49SHong Zhang EXTERN_C_END 38afaefe49SHong Zhang 39afaefe49SHong Zhang EXTERN_C_BEGIN 40afaefe49SHong Zhang #undef __FUNCT__ 41afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU" 42dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift) 43afaefe49SHong Zhang { 44c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 45afaefe49SHong Zhang 46afaefe49SHong Zhang PetscFunctionBegin; 47afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 48afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 49afaefe49SHong Zhang } else { 50afaefe49SHong Zhang dir->info.shiftnz = shift; 51afaefe49SHong Zhang } 52afaefe49SHong Zhang PetscFunctionReturn(0); 53afaefe49SHong Zhang } 54afaefe49SHong Zhang EXTERN_C_END 55afaefe49SHong Zhang 56afaefe49SHong Zhang EXTERN_C_BEGIN 57afaefe49SHong Zhang #undef __FUNCT__ 58afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU" 59dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift) 60afaefe49SHong Zhang { 61c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 62afaefe49SHong Zhang 63afaefe49SHong Zhang PetscFunctionBegin; 64fbf22428SSatish Balay if (shift) { 65fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 66fbf22428SSatish Balay dir->info.shiftpd = 1.0; 67fbf22428SSatish Balay } else { 68fbf22428SSatish Balay dir->info.shiftpd = 0.0; 69fbf22428SSatish Balay } 70afaefe49SHong Zhang PetscFunctionReturn(0); 71afaefe49SHong Zhang } 72afaefe49SHong Zhang EXTERN_C_END 73afaefe49SHong Zhang 74afaefe49SHong Zhang EXTERN_C_BEGIN 75afaefe49SHong Zhang #undef __FUNCT__ 762401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 772401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 789b54502bSHong Zhang { 799b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 809b54502bSHong Zhang 819b54502bSHong Zhang PetscFunctionBegin; 829b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 839b54502bSHong Zhang if (z == PETSC_DECIDE) { 849b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 859b54502bSHong Zhang } else { 869b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 879b54502bSHong Zhang } 889b54502bSHong Zhang PetscFunctionReturn(0); 899b54502bSHong Zhang } 909b54502bSHong Zhang EXTERN_C_END 919b54502bSHong Zhang 929b54502bSHong Zhang EXTERN_C_BEGIN 939b54502bSHong Zhang #undef __FUNCT__ 942401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 952401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag) 969b54502bSHong Zhang { 97c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 989b54502bSHong Zhang 999b54502bSHong Zhang PetscFunctionBegin; 1009b54502bSHong Zhang lu->reuseordering = flag; 1019b54502bSHong Zhang PetscFunctionReturn(0); 1029b54502bSHong Zhang } 1039b54502bSHong Zhang EXTERN_C_END 1049b54502bSHong Zhang 1059b54502bSHong Zhang EXTERN_C_BEGIN 1069b54502bSHong Zhang #undef __FUNCT__ 1072401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 1082401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag) 1099b54502bSHong Zhang { 110c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 1119b54502bSHong Zhang 1129b54502bSHong Zhang PetscFunctionBegin; 1139b54502bSHong Zhang lu->reusefill = flag; 1149b54502bSHong Zhang PetscFunctionReturn(0); 1159b54502bSHong Zhang } 1169b54502bSHong Zhang EXTERN_C_END 1179b54502bSHong Zhang 1189b54502bSHong Zhang #undef __FUNCT__ 1199b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 1209b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 1219b54502bSHong Zhang { 1229b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1239b54502bSHong Zhang PetscErrorCode ierr; 1249b54502bSHong Zhang PetscTruth flg,set; 1259989ab13SBarry Smith char tname[256], solvertype[64]; 1269b54502bSHong Zhang PetscFList ordlist; 1279b54502bSHong Zhang PetscReal tol; 1289b54502bSHong Zhang 1299b54502bSHong Zhang PetscFunctionBegin; 1309b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1319b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 1322401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 1339b54502bSHong Zhang if (flg) { 1342401956bSBarry Smith ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 1359b54502bSHong Zhang } 13655ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 1379b54502bSHong Zhang 1389f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 1399b54502bSHong Zhang if (flg) { 140afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1419b54502bSHong Zhang } 1429f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 1439f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 1449b54502bSHong Zhang if (flg) { 145afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 1469b54502bSHong Zhang } 147ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 1489b54502bSHong Zhang 1492401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 1509b54502bSHong Zhang if (flg) { 1512401956bSBarry Smith ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1529b54502bSHong Zhang } 1532401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 1549b54502bSHong Zhang if (flg) { 1552401956bSBarry Smith ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1569b54502bSHong Zhang } 1579b54502bSHong Zhang 1589b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 159e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 1609b54502bSHong Zhang if (flg) { 161e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang 1645c9eb25fSBarry Smith /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 165c7393fdbSBarry Smith ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",lu->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 1665c9eb25fSBarry Smith if (flg) { 167c7393fdbSBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 1685c9eb25fSBarry Smith } 1695c9eb25fSBarry Smith 1702401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1719b54502bSHong Zhang if (flg) { 1729b54502bSHong Zhang tol = PETSC_DECIDE; 1732401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1742401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1759b54502bSHong Zhang } 1769b54502bSHong Zhang 1772401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 1789b54502bSHong Zhang 1799b54502bSHong Zhang flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 1802401956bSBarry Smith ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 1819b54502bSHong Zhang if (set) { 1822401956bSBarry Smith ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 1839b54502bSHong Zhang } 1849b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1859b54502bSHong Zhang PetscFunctionReturn(0); 1869b54502bSHong Zhang } 1879b54502bSHong Zhang 1889b54502bSHong Zhang #undef __FUNCT__ 1899b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 1909b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 1919b54502bSHong Zhang { 1929b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1939b54502bSHong Zhang PetscErrorCode ierr; 1949b54502bSHong Zhang PetscTruth iascii,isstring; 1959b54502bSHong Zhang 1969b54502bSHong Zhang PetscFunctionBegin; 1979b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1989b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1999b54502bSHong Zhang if (iascii) { 2009b54502bSHong Zhang 2019b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 2029b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 2039b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 204a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr); 2050a29876aSHong Zhang if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 2069b54502bSHong Zhang if (lu->fact) { 207f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 208f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 209f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 210f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 211f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 212f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 2139b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 2149b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 215f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 216f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 217f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2189b54502bSHong Zhang } 2199b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 2209b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 2219b54502bSHong Zhang } else if (isstring) { 2229b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 2239b54502bSHong Zhang } else { 2249b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 2259b54502bSHong Zhang } 2269b54502bSHong Zhang PetscFunctionReturn(0); 2279b54502bSHong Zhang } 2289b54502bSHong Zhang 2299b54502bSHong Zhang #undef __FUNCT__ 230a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU" 231a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat) 2329b54502bSHong Zhang { 2339b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2349b54502bSHong Zhang 2359b54502bSHong Zhang PetscFunctionBegin; 2369b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 2379b54502bSHong Zhang *mat = dir->fact; 2389b54502bSHong Zhang PetscFunctionReturn(0); 2399b54502bSHong Zhang } 2409b54502bSHong Zhang 2419b54502bSHong Zhang #undef __FUNCT__ 2429b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 2439b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 2449b54502bSHong Zhang { 2459b54502bSHong Zhang PetscErrorCode ierr; 2469b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2479b54502bSHong Zhang 2489b54502bSHong Zhang PetscFunctionBegin; 2499b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 2509b54502bSHong Zhang 2519b54502bSHong Zhang if (dir->inplace) { 2529b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2539b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2549b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 25503c60df9SBarry Smith if (dir->row) { 25603c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 25703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 25803c60df9SBarry Smith } 2599b54502bSHong Zhang ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 2609b54502bSHong Zhang dir->fact = pc->pmat; 2619b54502bSHong Zhang } else { 2629b54502bSHong Zhang MatInfo info; 2639b54502bSHong Zhang if (!pc->setupcalled) { 2649b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2659b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2669b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2679b54502bSHong Zhang } 26803c60df9SBarry Smith if (dir->row) { 26903c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 27003c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 27103c60df9SBarry Smith } 2725c9eb25fSBarry Smith ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr); 2739b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2749b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2759b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 27652e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2779b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2789b54502bSHong Zhang if (!dir->reuseordering) { 2799b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2809b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2819b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2829b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2839b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2849b54502bSHong Zhang } 28503c60df9SBarry Smith if (dir->row) { 28603c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 28703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 28803c60df9SBarry Smith } 2899b54502bSHong Zhang } 2909b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 2915c9eb25fSBarry Smith ierr = MatGetFactor(pc->pmat,dir->solvertype,MAT_FACTOR_LU,&dir->fact);CHKERRQ(ierr); 2929b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2939b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2949b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 29552e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2969b54502bSHong Zhang } 2979b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 2989b54502bSHong Zhang } 2999b54502bSHong Zhang PetscFunctionReturn(0); 3009b54502bSHong Zhang } 3019b54502bSHong Zhang 3029b54502bSHong Zhang #undef __FUNCT__ 3039b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 3049b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 3059b54502bSHong Zhang { 3069b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3079b54502bSHong Zhang PetscErrorCode ierr; 3089b54502bSHong Zhang 3099b54502bSHong Zhang PetscFunctionBegin; 3109b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 3119b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 3129b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 3139b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3145c9eb25fSBarry Smith ierr = PetscStrfree(dir->solvertype);CHKERRQ(ierr); 3159b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 3169b54502bSHong Zhang PetscFunctionReturn(0); 3179b54502bSHong Zhang } 3189b54502bSHong Zhang 3199b54502bSHong Zhang #undef __FUNCT__ 3209b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 3219b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 3229b54502bSHong Zhang { 3239b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3249b54502bSHong Zhang PetscErrorCode ierr; 3259b54502bSHong Zhang 3269b54502bSHong Zhang PetscFunctionBegin; 3279b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 3289b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 3299b54502bSHong Zhang PetscFunctionReturn(0); 3309b54502bSHong Zhang } 3319b54502bSHong Zhang 3329b54502bSHong Zhang #undef __FUNCT__ 3339b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 3349b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 3359b54502bSHong Zhang { 3369b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3379b54502bSHong Zhang PetscErrorCode ierr; 3389b54502bSHong Zhang 3399b54502bSHong Zhang PetscFunctionBegin; 3409b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 3419b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 3429b54502bSHong Zhang PetscFunctionReturn(0); 3439b54502bSHong Zhang } 3449b54502bSHong Zhang 3459b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3469b54502bSHong Zhang 3479b54502bSHong Zhang EXTERN_C_BEGIN 3489b54502bSHong Zhang #undef __FUNCT__ 34955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU" 35055ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 3519b54502bSHong Zhang { 352c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 3539b54502bSHong Zhang 3549b54502bSHong Zhang PetscFunctionBegin; 3559b54502bSHong Zhang dir->info.fill = fill; 3569b54502bSHong Zhang PetscFunctionReturn(0); 3579b54502bSHong Zhang } 3589b54502bSHong Zhang EXTERN_C_END 3599b54502bSHong Zhang 3609b54502bSHong Zhang EXTERN_C_BEGIN 3619b54502bSHong Zhang #undef __FUNCT__ 3622401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 3632401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 3649b54502bSHong Zhang { 365c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 3669b54502bSHong Zhang 3679b54502bSHong Zhang PetscFunctionBegin; 3689b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3699b54502bSHong Zhang PetscFunctionReturn(0); 3709b54502bSHong Zhang } 3719b54502bSHong Zhang EXTERN_C_END 3729b54502bSHong Zhang 3739b54502bSHong Zhang EXTERN_C_BEGIN 3749b54502bSHong Zhang #undef __FUNCT__ 375e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 376e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering) 3779b54502bSHong Zhang { 3789b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3799b54502bSHong Zhang PetscErrorCode ierr; 3809b54502bSHong Zhang 3819b54502bSHong Zhang PetscFunctionBegin; 3829b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3839b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 3849b54502bSHong Zhang PetscFunctionReturn(0); 3859b54502bSHong Zhang } 3869b54502bSHong Zhang EXTERN_C_END 3879b54502bSHong Zhang 3889b54502bSHong Zhang EXTERN_C_BEGIN 3899b54502bSHong Zhang #undef __FUNCT__ 3902401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU" 3912401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 3929b54502bSHong Zhang { 3939b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3949b54502bSHong Zhang 3959b54502bSHong Zhang PetscFunctionBegin; 396a83599f4SBarry Smith if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol); 3979b54502bSHong Zhang dir->info.dtcol = dtcol; 3989b54502bSHong Zhang PetscFunctionReturn(0); 3999b54502bSHong Zhang } 4009b54502bSHong Zhang EXTERN_C_END 4019b54502bSHong Zhang 4029b54502bSHong Zhang EXTERN_C_BEGIN 4039b54502bSHong Zhang #undef __FUNCT__ 4042401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 4052401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 4069b54502bSHong Zhang { 4079b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 4089b54502bSHong Zhang 4099b54502bSHong Zhang PetscFunctionBegin; 4109b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 4119b54502bSHong Zhang PetscFunctionReturn(0); 4129b54502bSHong Zhang } 4139b54502bSHong Zhang EXTERN_C_END 4149b54502bSHong Zhang 4159b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 4169b54502bSHong Zhang 4179b54502bSHong Zhang #undef __FUNCT__ 4182401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 4199b54502bSHong Zhang /*@ 4202401956bSBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 4219b54502bSHong Zhang 4229b54502bSHong Zhang Collective on PC 4239b54502bSHong Zhang 4249b54502bSHong Zhang Input Parameters: 4259b54502bSHong Zhang + pc - the preconditioner context 4269b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 4279b54502bSHong Zhang 4289b54502bSHong Zhang Options Database Key: 4292401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal 4309b54502bSHong Zhang 4319b54502bSHong Zhang Level: intermediate 4329b54502bSHong Zhang 4339b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4349b54502bSHong Zhang 43555ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 4369b54502bSHong Zhang @*/ 4372401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 4389b54502bSHong Zhang { 4399b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4409b54502bSHong Zhang 4419b54502bSHong Zhang PetscFunctionBegin; 4429b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4432401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 4449b54502bSHong Zhang if (f) { 4459b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 4469b54502bSHong Zhang } 4479b54502bSHong Zhang PetscFunctionReturn(0); 4489b54502bSHong Zhang } 4499b54502bSHong Zhang 4509b54502bSHong Zhang #undef __FUNCT__ 451c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage" 452c62fe420SBarry Smith /*@ 453c7393fdbSBarry Smith PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 454c62fe420SBarry Smith 455c62fe420SBarry Smith Collective on PC 456c62fe420SBarry Smith 457c62fe420SBarry Smith Input Parameters: 458c62fe420SBarry Smith + pc - the preconditioner context 459c62fe420SBarry Smith - stype - for example, spooles, superlu, superlu_d 460c62fe420SBarry Smith 461c62fe420SBarry Smith Options Database Key: 462c7393fdbSBarry Smith . -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps 463c62fe420SBarry Smith 464c62fe420SBarry Smith Level: intermediate 465c62fe420SBarry Smith 466c62fe420SBarry Smith Note: 467c62fe420SBarry Smith By default this will use the PETSc factorization if it exists 468c62fe420SBarry Smith 469c62fe420SBarry Smith .keywords: PC, set, factorization, direct, fill 470c62fe420SBarry Smith 471c7393fdbSBarry Smith .seealso: MatGetFactor(), MatSolverPackage 472c62fe420SBarry Smith 473c62fe420SBarry Smith @*/ 474c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 475c62fe420SBarry Smith { 476c7393fdbSBarry Smith PetscErrorCode ierr,(*f)(PC,const MatSolverPackage); 477c62fe420SBarry Smith 478c62fe420SBarry Smith PetscFunctionBegin; 479c62fe420SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 480c7393fdbSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr); 481c62fe420SBarry Smith if (f) { 482c62fe420SBarry Smith ierr = (*f)(pc,stype);CHKERRQ(ierr); 483c62fe420SBarry Smith } 484c62fe420SBarry Smith PetscFunctionReturn(0); 485c62fe420SBarry Smith } 486c62fe420SBarry Smith 487c62fe420SBarry Smith #undef __FUNCT__ 48855ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill" 4899b54502bSHong Zhang /*@ 49055ba2a51SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 4919b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 4929b54502bSHong Zhang 4939b54502bSHong Zhang Collective on PC 4949b54502bSHong Zhang 4959b54502bSHong Zhang Input Parameters: 4969b54502bSHong Zhang + pc - the preconditioner context 4979b54502bSHong Zhang - fill - amount of expected fill 4989b54502bSHong Zhang 4999b54502bSHong Zhang Options Database Key: 50055ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 5019b54502bSHong Zhang 5029b54502bSHong Zhang Level: intermediate 5039b54502bSHong Zhang 5049b54502bSHong Zhang Note: 5059b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 5066cf91177SBarry Smith fill to expect. By running with the option -info PETSc will print the 5079b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 5089b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 5099b54502bSHong Zhang 5109b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 5119b54502bSHong Zhang 5129b54502bSHong Zhang @*/ 51355ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 5149b54502bSHong Zhang { 5159b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5169b54502bSHong Zhang 5179b54502bSHong Zhang PetscFunctionBegin; 5189b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5199b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 52055ba2a51SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5219b54502bSHong Zhang if (f) { 5229b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 5239b54502bSHong Zhang } 5249b54502bSHong Zhang PetscFunctionReturn(0); 5259b54502bSHong Zhang } 5269b54502bSHong Zhang 5279b54502bSHong Zhang #undef __FUNCT__ 5282401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace" 5299b54502bSHong Zhang /*@ 5302401956bSBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 5319b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 5329b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 5339b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 5349b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 5359b54502bSHong Zhang space. 5369b54502bSHong Zhang 5379b54502bSHong Zhang Collective on PC 5389b54502bSHong Zhang 5399b54502bSHong Zhang Input Parameters: 5409b54502bSHong Zhang . pc - the preconditioner context 5419b54502bSHong Zhang 5429b54502bSHong Zhang Options Database Key: 5432401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 5449b54502bSHong Zhang 5459b54502bSHong Zhang Notes: 5462401956bSBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 5479b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 5489b54502bSHong Zhang a call to KSPSetOperators(). 5499b54502bSHong Zhang This is because the Krylov space methods require an application of the 5509b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 5519b54502bSHong Zhang been factored in-place, replacing the original matrix. 5529b54502bSHong Zhang 5539b54502bSHong Zhang Level: intermediate 5549b54502bSHong Zhang 5559b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 5569b54502bSHong Zhang 5579b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 5589b54502bSHong Zhang @*/ 5592401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 5609b54502bSHong Zhang { 5619b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 5629b54502bSHong Zhang 5639b54502bSHong Zhang PetscFunctionBegin; 5649b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5652401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 5669b54502bSHong Zhang if (f) { 5679b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 5689b54502bSHong Zhang } 5699b54502bSHong Zhang PetscFunctionReturn(0); 5709b54502bSHong Zhang } 5719b54502bSHong Zhang 5729b54502bSHong Zhang #undef __FUNCT__ 573e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType" 5749b54502bSHong Zhang /*@C 575e5a9bf91SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 5769b54502bSHong Zhang be used in the LU factorization. 5779b54502bSHong Zhang 5789b54502bSHong Zhang Collective on PC 5799b54502bSHong Zhang 5809b54502bSHong Zhang Input Parameters: 5819b54502bSHong Zhang + pc - the preconditioner context 5829b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 5839b54502bSHong Zhang 5849b54502bSHong Zhang Options Database Key: 5852401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 5869b54502bSHong Zhang 5879b54502bSHong Zhang Level: intermediate 5889b54502bSHong Zhang 5899b54502bSHong Zhang Notes: nested dissection is used by default 5909b54502bSHong Zhang 5912b6de112SBarry Smith For Cholesky and ICC and the SBAIJ format reorderings are not available, 5922b6de112SBarry Smith since only the upper triangular part of the matrix is stored. You can use the 5932b6de112SBarry Smith SeqAIJ format in this case to get reorderings. 5942b6de112SBarry Smith 5959b54502bSHong Zhang @*/ 596e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering) 5979b54502bSHong Zhang { 598e36ff184SSatish Balay PetscErrorCode ierr,(*f)(PC,const MatOrderingType); 5999b54502bSHong Zhang 6009b54502bSHong Zhang PetscFunctionBegin; 601e5a9bf91SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 6029b54502bSHong Zhang if (f) { 6039b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 6049b54502bSHong Zhang } 6059b54502bSHong Zhang PetscFunctionReturn(0); 6069b54502bSHong Zhang } 6079b54502bSHong Zhang 6089b54502bSHong Zhang #undef __FUNCT__ 6092401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting" 6109b54502bSHong Zhang /*@ 6112401956bSBarry Smith PCFactorSetPivoting - Determines when pivoting is done during LU. 6129b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 6139b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 6149b54502bSHong Zhang 6159b54502bSHong Zhang Collective on PC 6169b54502bSHong Zhang 6179b54502bSHong Zhang Input Parameters: 6189b54502bSHong Zhang + pc - the preconditioner context 6199b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 6209b54502bSHong Zhang 6219b54502bSHong Zhang Options Database Key: 6222401956bSBarry Smith . -pc_factor_pivoting <dtcol> 6239b54502bSHong Zhang 6249b54502bSHong Zhang Level: intermediate 6259b54502bSHong Zhang 6262401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 6279b54502bSHong Zhang @*/ 6282401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 6299b54502bSHong Zhang { 6309b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 6319b54502bSHong Zhang 6329b54502bSHong Zhang PetscFunctionBegin; 6332401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 6349b54502bSHong Zhang if (f) { 6359b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 6369b54502bSHong Zhang } 6379b54502bSHong Zhang PetscFunctionReturn(0); 6389b54502bSHong Zhang } 6399b54502bSHong Zhang 6409b54502bSHong Zhang #undef __FUNCT__ 6412401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks" 6429b54502bSHong Zhang /*@ 6432401956bSBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 6449b54502bSHong Zhang with BAIJ or SBAIJ matrices 6459b54502bSHong Zhang 6469b54502bSHong Zhang Collective on PC 6479b54502bSHong Zhang 6489b54502bSHong Zhang Input Parameters: 6499b54502bSHong Zhang + pc - the preconditioner context 6509b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 6519b54502bSHong Zhang 6529b54502bSHong Zhang Options Database Key: 6532401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> 6549b54502bSHong Zhang 6559b54502bSHong Zhang Level: intermediate 6569b54502bSHong Zhang 6572401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 6589b54502bSHong Zhang @*/ 6592401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 6609b54502bSHong Zhang { 6619b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 6629b54502bSHong Zhang 6639b54502bSHong Zhang PetscFunctionBegin; 6642401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 6659b54502bSHong Zhang if (f) { 6669b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 6679b54502bSHong Zhang } 6689b54502bSHong Zhang PetscFunctionReturn(0); 6699b54502bSHong Zhang } 6709b54502bSHong Zhang 6719b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 6729b54502bSHong Zhang 6739b54502bSHong Zhang /*MC 6749b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 6759b54502bSHong Zhang 6769b54502bSHong Zhang Options Database Keys: 6772401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 678c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 6792401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 68055ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 6812401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 6822401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 6832401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 6849b54502bSHong Zhang stability of factorization. 685f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 686e22d95b2SBarry Smith . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 687f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 688e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 689e22d95b2SBarry Smith diagonal. 6909b54502bSHong Zhang 6919b54502bSHong Zhang Notes: Not all options work for all matrix formats 6929b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 6939b54502bSHong Zhang algorithms 6949b54502bSHong Zhang 6959b54502bSHong Zhang Level: beginner 6969b54502bSHong Zhang 6979b54502bSHong Zhang Concepts: LU factorization, direct solver 6989b54502bSHong Zhang 6999b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 7009b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 7019b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 7029b54502bSHong Zhang 7039b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 704a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 705e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 706e22d95b2SBarry Smith PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 7079b54502bSHong Zhang M*/ 7089b54502bSHong Zhang 7099b54502bSHong Zhang EXTERN_C_BEGIN 7109b54502bSHong Zhang #undef __FUNCT__ 7119b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 712dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 7139b54502bSHong Zhang { 7149b54502bSHong Zhang PetscErrorCode ierr; 7159b54502bSHong Zhang PetscMPIInt size; 7169b54502bSHong Zhang PC_LU *dir; 7179b54502bSHong Zhang 7189b54502bSHong Zhang PetscFunctionBegin; 71938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 7209b54502bSHong Zhang 7219b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 7229b54502bSHong Zhang dir->fact = 0; 7239b54502bSHong Zhang dir->inplace = PETSC_FALSE; 7249b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 7259b54502bSHong Zhang 7269b54502bSHong Zhang dir->info.fill = 5.0; 7279b54502bSHong Zhang dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 7280a29876aSHong Zhang dir->info.shiftnz = 0.0; 7299b54502bSHong Zhang dir->info.zeropivot = 1.e-12; 7309b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 731fbf22428SSatish Balay dir->info.shiftpd = 0.0; /* false */ 7329b54502bSHong Zhang dir->info.shift_fraction = 0.0; 7339b54502bSHong Zhang dir->col = 0; 7349b54502bSHong Zhang dir->row = 0; 7355c9eb25fSBarry Smith 736*43244d56SBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&dir->solvertype);CHKERRQ(ierr); 7377adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 7389b54502bSHong Zhang if (size == 1) { 7399b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 7409b54502bSHong Zhang } else { 7419b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 7429b54502bSHong Zhang } 7439b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 7449b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 7459b54502bSHong Zhang pc->data = (void*)dir; 7469b54502bSHong Zhang 7479b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 7489b54502bSHong Zhang pc->ops->apply = PCApply_LU; 7499b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 7509b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 7519b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 7529b54502bSHong Zhang pc->ops->view = PCView_LU; 7539b54502bSHong Zhang pc->ops->applyrichardson = 0; 754a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU; 7559b54502bSHong Zhang 756c7393fdbSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_LU", 757c7393fdbSBarry Smith PCFactorSetMatSolverPackage_LU);CHKERRQ(ierr); 758afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 759afaefe49SHong Zhang PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 760afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 761afaefe49SHong Zhang PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 762afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 763afaefe49SHong Zhang PCFactorSetShiftPd_LU);CHKERRQ(ierr); 764afaefe49SHong Zhang 76555ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 76655ba2a51SBarry Smith PCFactorSetFill_LU);CHKERRQ(ierr); 7672401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 7682401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 769e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 770e5a9bf91SBarry Smith PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 7712401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 7722401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 7732401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 7742401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 7752401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 7762401956bSBarry Smith PCFactorSetPivoting_LU);CHKERRQ(ierr); 7772401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 7782401956bSBarry Smith PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 7792401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 7802401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 7819b54502bSHong Zhang PetscFunctionReturn(0); 7829b54502bSHong Zhang } 7839b54502bSHong Zhang EXTERN_C_END 784