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 9*075768bcSBarry Smith #include "src/ksp/pc/impls/factor/lu/lu.h" /*I "petscpc.h" I*/ 109b54502bSHong Zhang 119b54502bSHong Zhang EXTERN_C_BEGIN 129b54502bSHong Zhang #undef __FUNCT__ 13c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage_LU" 14c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage_LU(PC pc,const MatSolverPackage stype) 15c62fe420SBarry Smith { 16c62fe420SBarry Smith PetscErrorCode ierr; 17c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 18c62fe420SBarry Smith 19c62fe420SBarry Smith PetscFunctionBegin; 20*075768bcSBarry Smith if (((PC_Factor*)lu)->fact) { 2135bd34faSBarry Smith const MatSolverPackage ltype; 2235bd34faSBarry Smith PetscTruth flg; 23*075768bcSBarry Smith ierr = MatFactorGetSolverPackage(((PC_Factor*)lu)->fact,<ype);CHKERRQ(ierr); 2435bd34faSBarry Smith ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr); 2535bd34faSBarry Smith if (!flg) { 2635bd34faSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used"); 2735bd34faSBarry Smith } 2835bd34faSBarry Smith } 29*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)lu)->solvertype);CHKERRQ(ierr); 30*075768bcSBarry Smith ierr = PetscStrallocpy(stype,&((PC_Factor*)lu)->solvertype);CHKERRQ(ierr); 31c62fe420SBarry Smith PetscFunctionReturn(0); 32c62fe420SBarry Smith } 33c62fe420SBarry Smith EXTERN_C_END 34c62fe420SBarry Smith 35c62fe420SBarry Smith EXTERN_C_BEGIN 36c62fe420SBarry Smith #undef __FUNCT__ 37afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU" 38dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z) 39afaefe49SHong Zhang { 40c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 41afaefe49SHong Zhang 42afaefe49SHong Zhang PetscFunctionBegin; 43*075768bcSBarry Smith ((PC_Factor*)lu)->info.zeropivot = z; 44afaefe49SHong Zhang PetscFunctionReturn(0); 45afaefe49SHong Zhang } 46afaefe49SHong Zhang EXTERN_C_END 47afaefe49SHong Zhang 48afaefe49SHong Zhang EXTERN_C_BEGIN 49afaefe49SHong Zhang #undef __FUNCT__ 50afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU" 51dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift) 52afaefe49SHong Zhang { 53c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 54afaefe49SHong Zhang 55afaefe49SHong Zhang PetscFunctionBegin; 56afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 57*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = 1.e-12; 58afaefe49SHong Zhang } else { 59*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = shift; 60afaefe49SHong Zhang } 61afaefe49SHong Zhang PetscFunctionReturn(0); 62afaefe49SHong Zhang } 63afaefe49SHong Zhang EXTERN_C_END 64afaefe49SHong Zhang 65afaefe49SHong Zhang EXTERN_C_BEGIN 66afaefe49SHong Zhang #undef __FUNCT__ 67afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU" 68dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift) 69afaefe49SHong Zhang { 70c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 71afaefe49SHong Zhang 72afaefe49SHong Zhang PetscFunctionBegin; 73fbf22428SSatish Balay if (shift) { 74*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 1.0; 75fbf22428SSatish Balay } else { 76*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 0.0; 77fbf22428SSatish Balay } 78afaefe49SHong Zhang PetscFunctionReturn(0); 79afaefe49SHong Zhang } 80afaefe49SHong Zhang EXTERN_C_END 81afaefe49SHong Zhang 82afaefe49SHong Zhang EXTERN_C_BEGIN 83afaefe49SHong Zhang #undef __FUNCT__ 842401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 852401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 869b54502bSHong Zhang { 879b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 889b54502bSHong Zhang 899b54502bSHong Zhang PetscFunctionBegin; 909b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 919b54502bSHong Zhang if (z == PETSC_DECIDE) { 929b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 939b54502bSHong Zhang } else { 949b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 959b54502bSHong Zhang } 969b54502bSHong Zhang PetscFunctionReturn(0); 979b54502bSHong Zhang } 989b54502bSHong Zhang EXTERN_C_END 999b54502bSHong Zhang 1009b54502bSHong Zhang EXTERN_C_BEGIN 1019b54502bSHong Zhang #undef __FUNCT__ 1022401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 1032401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag) 1049b54502bSHong Zhang { 105c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 1069b54502bSHong Zhang 1079b54502bSHong Zhang PetscFunctionBegin; 1089b54502bSHong Zhang lu->reuseordering = flag; 1099b54502bSHong Zhang PetscFunctionReturn(0); 1109b54502bSHong Zhang } 1119b54502bSHong Zhang EXTERN_C_END 1129b54502bSHong Zhang 1139b54502bSHong Zhang EXTERN_C_BEGIN 1149b54502bSHong Zhang #undef __FUNCT__ 1152401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 1162401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag) 1179b54502bSHong Zhang { 118c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 1199b54502bSHong Zhang 1209b54502bSHong Zhang PetscFunctionBegin; 1219b54502bSHong Zhang lu->reusefill = flag; 1229b54502bSHong Zhang PetscFunctionReturn(0); 1239b54502bSHong Zhang } 1249b54502bSHong Zhang EXTERN_C_END 1259b54502bSHong Zhang 1269b54502bSHong Zhang #undef __FUNCT__ 1279b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 1289b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 1299b54502bSHong Zhang { 1309b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1319b54502bSHong Zhang PetscErrorCode ierr; 1329b54502bSHong Zhang PetscTruth flg,set; 1339989ab13SBarry Smith char tname[256], solvertype[64]; 1349b54502bSHong Zhang PetscFList ordlist; 1359b54502bSHong Zhang PetscReal tol; 1369b54502bSHong Zhang 1379b54502bSHong Zhang PetscFunctionBegin; 1389b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1399b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 1402401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 1419b54502bSHong Zhang if (flg) { 1422401956bSBarry Smith ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 1439b54502bSHong Zhang } 144*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in LU/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);CHKERRQ(ierr); 1459b54502bSHong Zhang 1469f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 1479b54502bSHong Zhang if (flg) { 148afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1499b54502bSHong Zhang } 150*075768bcSBarry 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); 1519f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 1529b54502bSHong Zhang if (flg) { 153afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 1549b54502bSHong Zhang } 155*075768bcSBarry 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); 1569b54502bSHong Zhang 1572401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 1589b54502bSHong Zhang if (flg) { 1592401956bSBarry Smith ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1609b54502bSHong Zhang } 1612401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 1629b54502bSHong Zhang if (flg) { 1632401956bSBarry Smith ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1649b54502bSHong Zhang } 1659b54502bSHong Zhang 1669b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 167*075768bcSBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);CHKERRQ(ierr); 1689b54502bSHong Zhang if (flg) { 169e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 1709b54502bSHong Zhang } 1719b54502bSHong Zhang 1725c9eb25fSBarry Smith /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */ 173*075768bcSBarry Smith ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific LU solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);CHKERRQ(ierr); 1745c9eb25fSBarry Smith if (flg) { 175c7393fdbSBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr); 1765c9eb25fSBarry Smith } 1775c9eb25fSBarry Smith 1782401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1799b54502bSHong Zhang if (flg) { 1809b54502bSHong Zhang tol = PETSC_DECIDE; 1812401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1822401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1839b54502bSHong Zhang } 1849b54502bSHong Zhang 185*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",((PC_Factor*)lu)->info.dtcol,&((PC_Factor*)lu)->info.dtcol,&flg);CHKERRQ(ierr); 1869b54502bSHong Zhang 187*075768bcSBarry Smith flg = ((PC_Factor*)lu)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 1882401956bSBarry Smith ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 1899b54502bSHong Zhang if (set) { 1902401956bSBarry Smith ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 1919b54502bSHong Zhang } 1929b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1939b54502bSHong Zhang PetscFunctionReturn(0); 1949b54502bSHong Zhang } 1959b54502bSHong Zhang 1969b54502bSHong Zhang #undef __FUNCT__ 1979b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 1989b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 1999b54502bSHong Zhang { 2009b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 2019b54502bSHong Zhang PetscErrorCode ierr; 2029b54502bSHong Zhang PetscTruth iascii,isstring; 2039b54502bSHong Zhang 2049b54502bSHong Zhang PetscFunctionBegin; 2059b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 2069b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 2079b54502bSHong Zhang if (iascii) { 2089b54502bSHong Zhang 2099b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 2109b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 211*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr); 212*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr); 213*075768bcSBarry Smith if (((PC_Factor*)lu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 214*075768bcSBarry Smith if (((PC_Factor*)lu)->fact) { 215f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 216f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 217f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 218f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 219f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 220f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 221*075768bcSBarry Smith ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr); 2229b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 223f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 224f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 225f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2269b54502bSHong Zhang } 2279b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 2289b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 2299b54502bSHong Zhang } else if (isstring) { 230*075768bcSBarry Smith ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 2319b54502bSHong Zhang } else { 2329b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 2339b54502bSHong Zhang } 2349b54502bSHong Zhang PetscFunctionReturn(0); 2359b54502bSHong Zhang } 2369b54502bSHong Zhang 2379b54502bSHong Zhang #undef __FUNCT__ 238a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU" 239a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat) 2409b54502bSHong Zhang { 2419b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2429b54502bSHong Zhang 2439b54502bSHong Zhang PetscFunctionBegin; 244*075768bcSBarry Smith if (!((PC_Factor*)dir)->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 245*075768bcSBarry Smith *mat = ((PC_Factor*)dir)->fact; 2469b54502bSHong Zhang PetscFunctionReturn(0); 2479b54502bSHong Zhang } 2489b54502bSHong Zhang 2499b54502bSHong Zhang #undef __FUNCT__ 2509b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 2519b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 2529b54502bSHong Zhang { 2539b54502bSHong Zhang PetscErrorCode ierr; 2549b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2559b54502bSHong Zhang 2569b54502bSHong Zhang PetscFunctionBegin; 257*075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 2589b54502bSHong Zhang 2599b54502bSHong Zhang if (dir->inplace) { 2609b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2619b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 262*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 26303c60df9SBarry Smith if (dir->row) { 26403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 26503c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 26603c60df9SBarry Smith } 267*075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 268*075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 2699b54502bSHong Zhang } else { 2709b54502bSHong Zhang MatInfo info; 2719b54502bSHong Zhang if (!pc->setupcalled) { 272*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2739b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2749b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2759b54502bSHong Zhang } 27603c60df9SBarry Smith if (dir->row) { 27703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 27803c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 27903c60df9SBarry Smith } 280*075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 281*075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 282*075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2839b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 284*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 2859b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2869b54502bSHong Zhang if (!dir->reuseordering) { 2879b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2889b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 289*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2909b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2919b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2929b54502bSHong Zhang } 29303c60df9SBarry Smith if (dir->row) { 29403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 29503c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 29603c60df9SBarry Smith } 2979b54502bSHong Zhang } 298*075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 299*075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 300*075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 301*075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 3029b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 303*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 3049b54502bSHong Zhang } 305*075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3069b54502bSHong Zhang } 3079b54502bSHong Zhang PetscFunctionReturn(0); 3089b54502bSHong Zhang } 3099b54502bSHong Zhang 3109b54502bSHong Zhang #undef __FUNCT__ 3119b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 3129b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 3139b54502bSHong Zhang { 3149b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3159b54502bSHong Zhang PetscErrorCode ierr; 3169b54502bSHong Zhang 3179b54502bSHong Zhang PetscFunctionBegin; 318*075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 3199b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 3209b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 321*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 322*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3239b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 3249b54502bSHong Zhang PetscFunctionReturn(0); 3259b54502bSHong Zhang } 3269b54502bSHong Zhang 3279b54502bSHong Zhang #undef __FUNCT__ 3289b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 3299b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 3309b54502bSHong Zhang { 3319b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3329b54502bSHong Zhang PetscErrorCode ierr; 3339b54502bSHong Zhang 3349b54502bSHong Zhang PetscFunctionBegin; 3359b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 336*075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 3379b54502bSHong Zhang PetscFunctionReturn(0); 3389b54502bSHong Zhang } 3399b54502bSHong Zhang 3409b54502bSHong Zhang #undef __FUNCT__ 3419b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 3429b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 3439b54502bSHong Zhang { 3449b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3459b54502bSHong Zhang PetscErrorCode ierr; 3469b54502bSHong Zhang 3479b54502bSHong Zhang PetscFunctionBegin; 3489b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 349*075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 3509b54502bSHong Zhang PetscFunctionReturn(0); 3519b54502bSHong Zhang } 3529b54502bSHong Zhang 3539b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3549b54502bSHong Zhang 3559b54502bSHong Zhang EXTERN_C_BEGIN 3569b54502bSHong Zhang #undef __FUNCT__ 35755ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU" 35855ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 3599b54502bSHong Zhang { 360c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 3619b54502bSHong Zhang 3629b54502bSHong Zhang PetscFunctionBegin; 363*075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = fill; 3649b54502bSHong Zhang PetscFunctionReturn(0); 3659b54502bSHong Zhang } 3669b54502bSHong Zhang EXTERN_C_END 3679b54502bSHong Zhang 3689b54502bSHong Zhang EXTERN_C_BEGIN 3699b54502bSHong Zhang #undef __FUNCT__ 3702401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 3712401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 3729b54502bSHong Zhang { 373c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 3749b54502bSHong Zhang 3759b54502bSHong Zhang PetscFunctionBegin; 3769b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3779b54502bSHong Zhang PetscFunctionReturn(0); 3789b54502bSHong Zhang } 3799b54502bSHong Zhang EXTERN_C_END 3809b54502bSHong Zhang 3819b54502bSHong Zhang EXTERN_C_BEGIN 3829b54502bSHong Zhang #undef __FUNCT__ 383e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 384e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,const MatOrderingType ordering) 3859b54502bSHong Zhang { 3869b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3879b54502bSHong Zhang PetscErrorCode ierr; 3889b54502bSHong Zhang 3899b54502bSHong Zhang PetscFunctionBegin; 390*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 391*075768bcSBarry Smith ierr = PetscStrallocpy(ordering,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3929b54502bSHong Zhang PetscFunctionReturn(0); 3939b54502bSHong Zhang } 3949b54502bSHong Zhang EXTERN_C_END 3959b54502bSHong Zhang 3969b54502bSHong Zhang EXTERN_C_BEGIN 3979b54502bSHong Zhang #undef __FUNCT__ 3982401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU" 3992401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 4009b54502bSHong Zhang { 4019b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 4029b54502bSHong Zhang 4039b54502bSHong Zhang PetscFunctionBegin; 404a83599f4SBarry 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); 405*075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = dtcol; 4069b54502bSHong Zhang PetscFunctionReturn(0); 4079b54502bSHong Zhang } 4089b54502bSHong Zhang EXTERN_C_END 4099b54502bSHong Zhang 4109b54502bSHong Zhang EXTERN_C_BEGIN 4119b54502bSHong Zhang #undef __FUNCT__ 4122401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 4132401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 4149b54502bSHong Zhang { 4159b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 4169b54502bSHong Zhang 4179b54502bSHong Zhang PetscFunctionBegin; 418*075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = pivot ? 1.0 : 0.0; 4199b54502bSHong Zhang PetscFunctionReturn(0); 4209b54502bSHong Zhang } 4219b54502bSHong Zhang EXTERN_C_END 4229b54502bSHong Zhang 4239b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 4249b54502bSHong Zhang 4259b54502bSHong Zhang #undef __FUNCT__ 4262401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 4279b54502bSHong Zhang /*@ 4282401956bSBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 4299b54502bSHong Zhang 4309b54502bSHong Zhang Collective on PC 4319b54502bSHong Zhang 4329b54502bSHong Zhang Input Parameters: 4339b54502bSHong Zhang + pc - the preconditioner context 4349b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 4359b54502bSHong Zhang 4369b54502bSHong Zhang Options Database Key: 4372401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal 4389b54502bSHong Zhang 4399b54502bSHong Zhang Level: intermediate 4409b54502bSHong Zhang 4419b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4429b54502bSHong Zhang 44355ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 4449b54502bSHong Zhang @*/ 4452401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 4469b54502bSHong Zhang { 4479b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4489b54502bSHong Zhang 4499b54502bSHong Zhang PetscFunctionBegin; 4509b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4512401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 4529b54502bSHong Zhang if (f) { 4539b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 4549b54502bSHong Zhang } 4559b54502bSHong Zhang PetscFunctionReturn(0); 4569b54502bSHong Zhang } 4579b54502bSHong Zhang 4589b54502bSHong Zhang #undef __FUNCT__ 459c7393fdbSBarry Smith #define __FUNCT__ "PCFactorSetMatSolverPackage" 460c62fe420SBarry Smith /*@ 461c7393fdbSBarry Smith PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization 462c62fe420SBarry Smith 463c62fe420SBarry Smith Collective on PC 464c62fe420SBarry Smith 465c62fe420SBarry Smith Input Parameters: 466c62fe420SBarry Smith + pc - the preconditioner context 46735bd34faSBarry Smith - stype - for example, spooles, superlu, superlu_dist 468c62fe420SBarry Smith 469c62fe420SBarry Smith Options Database Key: 470c7393fdbSBarry Smith . -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps 471c62fe420SBarry Smith 472c62fe420SBarry Smith Level: intermediate 473c62fe420SBarry Smith 474c62fe420SBarry Smith Note: 475c62fe420SBarry Smith By default this will use the PETSc factorization if it exists 476c62fe420SBarry Smith 47735bd34faSBarry Smith Use PCFactorGetMatrix(pc,&mat); MatFactorGetPackage(mat,&package); to see what package is being used 47835bd34faSBarry Smith 47935bd34faSBarry Smith 480c62fe420SBarry Smith .keywords: PC, set, factorization, direct, fill 481c62fe420SBarry Smith 482c7393fdbSBarry Smith .seealso: MatGetFactor(), MatSolverPackage 483c62fe420SBarry Smith 484c62fe420SBarry Smith @*/ 485c7393fdbSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype) 486c62fe420SBarry Smith { 487c7393fdbSBarry Smith PetscErrorCode ierr,(*f)(PC,const MatSolverPackage); 488c62fe420SBarry Smith 489c62fe420SBarry Smith PetscFunctionBegin; 490c62fe420SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 491c7393fdbSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",(void (**)(void))&f);CHKERRQ(ierr); 492c62fe420SBarry Smith if (f) { 493c62fe420SBarry Smith ierr = (*f)(pc,stype);CHKERRQ(ierr); 494c62fe420SBarry Smith } 495c62fe420SBarry Smith PetscFunctionReturn(0); 496c62fe420SBarry Smith } 497c62fe420SBarry Smith 498c62fe420SBarry Smith #undef __FUNCT__ 49955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill" 5009b54502bSHong Zhang /*@ 50155ba2a51SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 5029b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 5039b54502bSHong Zhang 5049b54502bSHong Zhang Collective on PC 5059b54502bSHong Zhang 5069b54502bSHong Zhang Input Parameters: 5079b54502bSHong Zhang + pc - the preconditioner context 5089b54502bSHong Zhang - fill - amount of expected fill 5099b54502bSHong Zhang 5109b54502bSHong Zhang Options Database Key: 51155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 5129b54502bSHong Zhang 5139b54502bSHong Zhang Level: intermediate 5149b54502bSHong Zhang 5159b54502bSHong Zhang Note: 5169b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 5176cf91177SBarry Smith fill to expect. By running with the option -info PETSc will print the 5189b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 5199b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 5209b54502bSHong Zhang 5219b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 5229b54502bSHong Zhang 5239b54502bSHong Zhang @*/ 52455ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 5259b54502bSHong Zhang { 5269b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5279b54502bSHong Zhang 5289b54502bSHong Zhang PetscFunctionBegin; 5299b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5309b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 53155ba2a51SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5329b54502bSHong Zhang if (f) { 5339b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 5349b54502bSHong Zhang } 5359b54502bSHong Zhang PetscFunctionReturn(0); 5369b54502bSHong Zhang } 5379b54502bSHong Zhang 5389b54502bSHong Zhang #undef __FUNCT__ 5392401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace" 5409b54502bSHong Zhang /*@ 5412401956bSBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 5429b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 5439b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 5449b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 5459b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 5469b54502bSHong Zhang space. 5479b54502bSHong Zhang 5489b54502bSHong Zhang Collective on PC 5499b54502bSHong Zhang 5509b54502bSHong Zhang Input Parameters: 5519b54502bSHong Zhang . pc - the preconditioner context 5529b54502bSHong Zhang 5539b54502bSHong Zhang Options Database Key: 5542401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 5559b54502bSHong Zhang 5569b54502bSHong Zhang Notes: 5572401956bSBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 5589b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 5599b54502bSHong Zhang a call to KSPSetOperators(). 5609b54502bSHong Zhang This is because the Krylov space methods require an application of the 5619b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 5629b54502bSHong Zhang been factored in-place, replacing the original matrix. 5639b54502bSHong Zhang 5649b54502bSHong Zhang Level: intermediate 5659b54502bSHong Zhang 5669b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 5679b54502bSHong Zhang 5689b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 5699b54502bSHong Zhang @*/ 5702401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 5719b54502bSHong Zhang { 5729b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 5739b54502bSHong Zhang 5749b54502bSHong Zhang PetscFunctionBegin; 5759b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5762401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 5779b54502bSHong Zhang if (f) { 5789b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 5799b54502bSHong Zhang } 5809b54502bSHong Zhang PetscFunctionReturn(0); 5819b54502bSHong Zhang } 5829b54502bSHong Zhang 5839b54502bSHong Zhang #undef __FUNCT__ 584e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType" 5859b54502bSHong Zhang /*@C 586e5a9bf91SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 5879b54502bSHong Zhang be used in the LU factorization. 5889b54502bSHong Zhang 5899b54502bSHong Zhang Collective on PC 5909b54502bSHong Zhang 5919b54502bSHong Zhang Input Parameters: 5929b54502bSHong Zhang + pc - the preconditioner context 5939b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 5949b54502bSHong Zhang 5959b54502bSHong Zhang Options Database Key: 5962401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 5979b54502bSHong Zhang 5989b54502bSHong Zhang Level: intermediate 5999b54502bSHong Zhang 6009b54502bSHong Zhang Notes: nested dissection is used by default 6019b54502bSHong Zhang 6022b6de112SBarry Smith For Cholesky and ICC and the SBAIJ format reorderings are not available, 6032b6de112SBarry Smith since only the upper triangular part of the matrix is stored. You can use the 6042b6de112SBarry Smith SeqAIJ format in this case to get reorderings. 6052b6de112SBarry Smith 6069b54502bSHong Zhang @*/ 607e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering) 6089b54502bSHong Zhang { 609e36ff184SSatish Balay PetscErrorCode ierr,(*f)(PC,const MatOrderingType); 6109b54502bSHong Zhang 6119b54502bSHong Zhang PetscFunctionBegin; 612e5a9bf91SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 6139b54502bSHong Zhang if (f) { 6149b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 6159b54502bSHong Zhang } 6169b54502bSHong Zhang PetscFunctionReturn(0); 6179b54502bSHong Zhang } 6189b54502bSHong Zhang 6199b54502bSHong Zhang #undef __FUNCT__ 6202401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting" 6219b54502bSHong Zhang /*@ 6222401956bSBarry Smith PCFactorSetPivoting - Determines when pivoting is done during LU. 6239b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 6249b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 6259b54502bSHong Zhang 6269b54502bSHong Zhang Collective on PC 6279b54502bSHong Zhang 6289b54502bSHong Zhang Input Parameters: 6299b54502bSHong Zhang + pc - the preconditioner context 6309b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 6319b54502bSHong Zhang 6329b54502bSHong Zhang Options Database Key: 6332401956bSBarry Smith . -pc_factor_pivoting <dtcol> 6349b54502bSHong Zhang 6359b54502bSHong Zhang Level: intermediate 6369b54502bSHong Zhang 6372401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 6389b54502bSHong Zhang @*/ 6392401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 6409b54502bSHong Zhang { 6419b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 6429b54502bSHong Zhang 6439b54502bSHong Zhang PetscFunctionBegin; 6442401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 6459b54502bSHong Zhang if (f) { 6469b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 6479b54502bSHong Zhang } 6489b54502bSHong Zhang PetscFunctionReturn(0); 6499b54502bSHong Zhang } 6509b54502bSHong Zhang 6519b54502bSHong Zhang #undef __FUNCT__ 6522401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks" 6539b54502bSHong Zhang /*@ 6542401956bSBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 6559b54502bSHong Zhang with BAIJ or SBAIJ matrices 6569b54502bSHong Zhang 6579b54502bSHong Zhang Collective on PC 6589b54502bSHong Zhang 6599b54502bSHong Zhang Input Parameters: 6609b54502bSHong Zhang + pc - the preconditioner context 6619b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 6629b54502bSHong Zhang 6639b54502bSHong Zhang Options Database Key: 6642401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> 6659b54502bSHong Zhang 6669b54502bSHong Zhang Level: intermediate 6679b54502bSHong Zhang 6682401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 6699b54502bSHong Zhang @*/ 6702401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 6719b54502bSHong Zhang { 6729b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 6739b54502bSHong Zhang 6749b54502bSHong Zhang PetscFunctionBegin; 6752401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 6769b54502bSHong Zhang if (f) { 6779b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 6789b54502bSHong Zhang } 6799b54502bSHong Zhang PetscFunctionReturn(0); 6809b54502bSHong Zhang } 6819b54502bSHong Zhang 6829b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 6839b54502bSHong Zhang 6849b54502bSHong Zhang /*MC 6859b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 6869b54502bSHong Zhang 6879b54502bSHong Zhang Options Database Keys: 6882401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 689c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 6902401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 69155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 6922401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 6932401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 6942401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 6959b54502bSHong Zhang stability of factorization. 696f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 697e22d95b2SBarry Smith . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 698f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 699e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 700e22d95b2SBarry Smith diagonal. 7019b54502bSHong Zhang 7029b54502bSHong Zhang Notes: Not all options work for all matrix formats 7039b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 7049b54502bSHong Zhang algorithms 7059b54502bSHong Zhang 7069b54502bSHong Zhang Level: beginner 7079b54502bSHong Zhang 7089b54502bSHong Zhang Concepts: LU factorization, direct solver 7099b54502bSHong Zhang 7109b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 7119b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 7129b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 7139b54502bSHong Zhang 7149b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 715a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 716e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 717e22d95b2SBarry Smith PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 7189b54502bSHong Zhang M*/ 7199b54502bSHong Zhang 7209b54502bSHong Zhang EXTERN_C_BEGIN 7219b54502bSHong Zhang #undef __FUNCT__ 7229b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 723dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 7249b54502bSHong Zhang { 7259b54502bSHong Zhang PetscErrorCode ierr; 7269b54502bSHong Zhang PetscMPIInt size; 7279b54502bSHong Zhang PC_LU *dir; 7289b54502bSHong Zhang 7299b54502bSHong Zhang PetscFunctionBegin; 73038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 7319b54502bSHong Zhang 732*075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 733*075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 7349b54502bSHong Zhang dir->inplace = PETSC_FALSE; 7359b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 7369b54502bSHong Zhang 737*075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 738*075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 739*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = 0.0; 740*075768bcSBarry Smith ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 741*075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 742*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 0.0; /* false */ 7439b54502bSHong Zhang dir->col = 0; 7449b54502bSHong Zhang dir->row = 0; 7455c9eb25fSBarry Smith 746*075768bcSBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 7477adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 7489b54502bSHong Zhang if (size == 1) { 749*075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 7509b54502bSHong Zhang } else { 751*075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 7529b54502bSHong Zhang } 7539b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 7549b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 7559b54502bSHong Zhang pc->data = (void*)dir; 7569b54502bSHong Zhang 7579b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 7589b54502bSHong Zhang pc->ops->apply = PCApply_LU; 7599b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 7609b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 7619b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 7629b54502bSHong Zhang pc->ops->view = PCView_LU; 7639b54502bSHong Zhang pc->ops->applyrichardson = 0; 764a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU; 7659b54502bSHong Zhang 766c7393fdbSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_LU", 767c7393fdbSBarry Smith PCFactorSetMatSolverPackage_LU);CHKERRQ(ierr); 768afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 769afaefe49SHong Zhang PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 770afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 771afaefe49SHong Zhang PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 772afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 773afaefe49SHong Zhang PCFactorSetShiftPd_LU);CHKERRQ(ierr); 774afaefe49SHong Zhang 77555ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 77655ba2a51SBarry Smith PCFactorSetFill_LU);CHKERRQ(ierr); 7772401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 7782401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 779e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 780e5a9bf91SBarry Smith PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 7812401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 7822401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 7832401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 7842401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 7852401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 7862401956bSBarry Smith PCFactorSetPivoting_LU);CHKERRQ(ierr); 7872401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 7882401956bSBarry Smith PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 7892401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 7902401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 7919b54502bSHong Zhang PetscFunctionReturn(0); 7929b54502bSHong Zhang } 7939b54502bSHong Zhang EXTERN_C_END 794