19b54502bSHong Zhang /* 29b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 39b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 49b54502bSHong Zhang a direct solver. 59b54502bSHong Zhang */ 6ee45ca4aSHong Zhang 79b54502bSHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 8ee45ca4aSHong Zhang #include "src/ksp/pc/impls/factor/lu/lu.h" 99b54502bSHong Zhang 109b54502bSHong Zhang EXTERN_C_BEGIN 119b54502bSHong Zhang #undef __FUNCT__ 12afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU" 13afaefe49SHong Zhang PetscErrorCode PCFactorSetZeroPivot_LU(PC pc,PetscReal z) 14afaefe49SHong Zhang { 15afaefe49SHong Zhang PC_LU *lu; 16afaefe49SHong Zhang 17afaefe49SHong Zhang PetscFunctionBegin; 18afaefe49SHong Zhang lu = (PC_LU*)pc->data; 19afaefe49SHong Zhang lu->info.zeropivot = z; 20afaefe49SHong Zhang PetscFunctionReturn(0); 21afaefe49SHong Zhang } 22afaefe49SHong Zhang EXTERN_C_END 23afaefe49SHong Zhang 24afaefe49SHong Zhang EXTERN_C_BEGIN 25afaefe49SHong Zhang #undef __FUNCT__ 26afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU" 27afaefe49SHong Zhang PetscErrorCode PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift) 28afaefe49SHong Zhang { 29afaefe49SHong Zhang PC_LU *dir; 30afaefe49SHong Zhang 31afaefe49SHong Zhang PetscFunctionBegin; 32afaefe49SHong Zhang dir = (PC_LU*)pc->data; 33afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 34afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 35afaefe49SHong Zhang } else { 36afaefe49SHong Zhang dir->info.shiftnz = shift; 37afaefe49SHong Zhang } 38afaefe49SHong Zhang PetscFunctionReturn(0); 39afaefe49SHong Zhang } 40afaefe49SHong Zhang EXTERN_C_END 41afaefe49SHong Zhang 42afaefe49SHong Zhang EXTERN_C_BEGIN 43afaefe49SHong Zhang #undef __FUNCT__ 44afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU" 45afaefe49SHong Zhang PetscErrorCode PCFactorSetShiftPd_LU(PC pc,PetscTruth shift) 46afaefe49SHong Zhang { 47afaefe49SHong Zhang PC_LU *dir; 48afaefe49SHong Zhang 49afaefe49SHong Zhang PetscFunctionBegin; 50afaefe49SHong Zhang dir = (PC_LU*)pc->data; 51afaefe49SHong Zhang dir->info.shiftpd = shift; 52afaefe49SHong Zhang if (shift) dir->info.shift_fraction = 0.0; 53afaefe49SHong Zhang PetscFunctionReturn(0); 54afaefe49SHong Zhang } 55afaefe49SHong Zhang EXTERN_C_END 56afaefe49SHong Zhang 57afaefe49SHong Zhang EXTERN_C_BEGIN 58afaefe49SHong Zhang #undef __FUNCT__ 599b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU" 609b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 619b54502bSHong Zhang { 629b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 639b54502bSHong Zhang 649b54502bSHong Zhang PetscFunctionBegin; 659b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 669b54502bSHong Zhang if (z == PETSC_DECIDE) { 679b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 689b54502bSHong Zhang } else { 699b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 709b54502bSHong Zhang } 719b54502bSHong Zhang PetscFunctionReturn(0); 729b54502bSHong Zhang } 739b54502bSHong Zhang EXTERN_C_END 749b54502bSHong Zhang 759b54502bSHong Zhang EXTERN_C_BEGIN 769b54502bSHong Zhang #undef __FUNCT__ 779b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering_LU" 789b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag) 799b54502bSHong Zhang { 809b54502bSHong Zhang PC_LU *lu; 819b54502bSHong Zhang 829b54502bSHong Zhang PetscFunctionBegin; 839b54502bSHong Zhang lu = (PC_LU*)pc->data; 849b54502bSHong Zhang lu->reuseordering = flag; 859b54502bSHong Zhang PetscFunctionReturn(0); 869b54502bSHong Zhang } 879b54502bSHong Zhang EXTERN_C_END 889b54502bSHong Zhang 899b54502bSHong Zhang EXTERN_C_BEGIN 909b54502bSHong Zhang #undef __FUNCT__ 919b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill_LU" 929b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag) 939b54502bSHong Zhang { 949b54502bSHong Zhang PC_LU *lu; 959b54502bSHong Zhang 969b54502bSHong Zhang PetscFunctionBegin; 979b54502bSHong Zhang lu = (PC_LU*)pc->data; 989b54502bSHong Zhang lu->reusefill = flag; 999b54502bSHong Zhang PetscFunctionReturn(0); 1009b54502bSHong Zhang } 1019b54502bSHong Zhang EXTERN_C_END 1029b54502bSHong Zhang 1039b54502bSHong Zhang #undef __FUNCT__ 1049b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 1059b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 1069b54502bSHong Zhang { 1079b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1089b54502bSHong Zhang PetscErrorCode ierr; 1099b54502bSHong Zhang PetscTruth flg,set; 1109b54502bSHong Zhang char tname[256]; 1119b54502bSHong Zhang PetscFList ordlist; 1129b54502bSHong Zhang PetscReal tol; 1139b54502bSHong Zhang 1149b54502bSHong Zhang PetscFunctionBegin; 1159b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1169b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 1179b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr); 1189b54502bSHong Zhang if (flg) { 1199b54502bSHong Zhang ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr); 1209b54502bSHong Zhang } 1219b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);CHKERRQ(ierr); 1229b54502bSHong Zhang 123ee45ca4aSHong Zhang ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 1249b54502bSHong Zhang if (flg) { 125afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1269b54502bSHong Zhang } 1270a29876aSHong Zhang ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 128ee45ca4aSHong Zhang ierr = PetscOptionsName("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 1299b54502bSHong Zhang if (flg) { 130afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 1319b54502bSHong Zhang } 132ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 1339b54502bSHong Zhang 1349b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr); 1359b54502bSHong Zhang if (flg) { 1369b54502bSHong Zhang ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1379b54502bSHong Zhang } 1389b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr); 1399b54502bSHong Zhang if (flg) { 1409b54502bSHong Zhang ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1419b54502bSHong Zhang } 1429b54502bSHong Zhang 1439b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 1449b54502bSHong Zhang ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 1459b54502bSHong Zhang if (flg) { 1469b54502bSHong Zhang ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr); 1479b54502bSHong Zhang } 1489b54502bSHong Zhang 1499b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1509b54502bSHong Zhang if (flg) { 1519b54502bSHong Zhang tol = PETSC_DECIDE; 1529b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1539b54502bSHong Zhang ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1549b54502bSHong Zhang } 1559b54502bSHong Zhang 1569b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 1579b54502bSHong Zhang 1589b54502bSHong Zhang flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 1599b54502bSHong Zhang ierr = PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 1609b54502bSHong Zhang if (set) { 1619b54502bSHong Zhang ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1649b54502bSHong Zhang PetscFunctionReturn(0); 1659b54502bSHong Zhang } 1669b54502bSHong Zhang 1679b54502bSHong Zhang #undef __FUNCT__ 1689b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 1699b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 1709b54502bSHong Zhang { 1719b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1729b54502bSHong Zhang PetscErrorCode ierr; 1739b54502bSHong Zhang PetscTruth iascii,isstring; 1749b54502bSHong Zhang 1759b54502bSHong Zhang PetscFunctionBegin; 1769b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1779b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1789b54502bSHong Zhang if (iascii) { 1799b54502bSHong Zhang MatInfo info; 1809b54502bSHong Zhang 1819b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 1829b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 1839b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 1849b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %g\n",lu->info.zeropivot);CHKERRQ(ierr); 1850a29876aSHong Zhang if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 1869b54502bSHong Zhang if (lu->fact) { 1879b54502bSHong Zhang ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1889b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU nonzeros %g\n",info.nz_used);CHKERRQ(ierr); 1899b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr); 1909b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 1919b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1929b54502bSHong Zhang } 1939b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 1949b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1959b54502bSHong Zhang } else if (isstring) { 1969b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 1979b54502bSHong Zhang } else { 1989b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 1999b54502bSHong Zhang } 2009b54502bSHong Zhang PetscFunctionReturn(0); 2019b54502bSHong Zhang } 2029b54502bSHong Zhang 2039b54502bSHong Zhang #undef __FUNCT__ 2049b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_LU" 2059b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat) 2069b54502bSHong Zhang { 2079b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2089b54502bSHong Zhang 2099b54502bSHong Zhang PetscFunctionBegin; 2109b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 2119b54502bSHong Zhang *mat = dir->fact; 2129b54502bSHong Zhang PetscFunctionReturn(0); 2139b54502bSHong Zhang } 2149b54502bSHong Zhang 2159b54502bSHong Zhang #undef __FUNCT__ 2169b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 2179b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 2189b54502bSHong Zhang { 2199b54502bSHong Zhang PetscErrorCode ierr; 2209b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2219b54502bSHong Zhang 2229b54502bSHong Zhang PetscFunctionBegin; 2239b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 2249b54502bSHong Zhang 2259b54502bSHong Zhang if (dir->inplace) { 2269b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2279b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2289b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2299b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2309b54502bSHong Zhang ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 2319b54502bSHong Zhang dir->fact = pc->pmat; 2329b54502bSHong Zhang } else { 2339b54502bSHong Zhang MatInfo info; 2349b54502bSHong Zhang if (!pc->setupcalled) { 2359b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2369b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2379b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2389b54502bSHong Zhang } 2399b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2409b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2419b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2429b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 243*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2449b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2459b54502bSHong Zhang if (!dir->reuseordering) { 2469b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2479b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2489b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2499b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2509b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2519b54502bSHong Zhang } 2529b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2539b54502bSHong Zhang } 2549b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 2559b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2569b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2579b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 258*52e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2599b54502bSHong Zhang } 2609b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 2619b54502bSHong Zhang } 2629b54502bSHong Zhang PetscFunctionReturn(0); 2639b54502bSHong Zhang } 2649b54502bSHong Zhang 2659b54502bSHong Zhang #undef __FUNCT__ 2669b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 2679b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 2689b54502bSHong Zhang { 2699b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2709b54502bSHong Zhang PetscErrorCode ierr; 2719b54502bSHong Zhang 2729b54502bSHong Zhang PetscFunctionBegin; 2739b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 2749b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2759b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2769b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 2779b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 2789b54502bSHong Zhang PetscFunctionReturn(0); 2799b54502bSHong Zhang } 2809b54502bSHong Zhang 2819b54502bSHong Zhang #undef __FUNCT__ 2829b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 2839b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 2849b54502bSHong Zhang { 2859b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2869b54502bSHong Zhang PetscErrorCode ierr; 2879b54502bSHong Zhang 2889b54502bSHong Zhang PetscFunctionBegin; 2899b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 2909b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 2919b54502bSHong Zhang PetscFunctionReturn(0); 2929b54502bSHong Zhang } 2939b54502bSHong Zhang 2949b54502bSHong Zhang #undef __FUNCT__ 2959b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2969b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2979b54502bSHong Zhang { 2989b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2999b54502bSHong Zhang PetscErrorCode ierr; 3009b54502bSHong Zhang 3019b54502bSHong Zhang PetscFunctionBegin; 3029b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 3039b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 3049b54502bSHong Zhang PetscFunctionReturn(0); 3059b54502bSHong Zhang } 3069b54502bSHong Zhang 3079b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3089b54502bSHong Zhang 3099b54502bSHong Zhang EXTERN_C_BEGIN 3109b54502bSHong Zhang #undef __FUNCT__ 3119b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill_LU" 3129b54502bSHong Zhang PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill) 3139b54502bSHong Zhang { 3149b54502bSHong Zhang PC_LU *dir; 3159b54502bSHong Zhang 3169b54502bSHong Zhang PetscFunctionBegin; 3179b54502bSHong Zhang dir = (PC_LU*)pc->data; 3189b54502bSHong Zhang dir->info.fill = fill; 3199b54502bSHong Zhang PetscFunctionReturn(0); 3209b54502bSHong Zhang } 3219b54502bSHong Zhang EXTERN_C_END 3229b54502bSHong Zhang 3239b54502bSHong Zhang EXTERN_C_BEGIN 3249b54502bSHong Zhang #undef __FUNCT__ 3259b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace_LU" 3269b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace_LU(PC pc) 3279b54502bSHong Zhang { 3289b54502bSHong Zhang PC_LU *dir; 3299b54502bSHong Zhang 3309b54502bSHong Zhang PetscFunctionBegin; 3319b54502bSHong Zhang dir = (PC_LU*)pc->data; 3329b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3339b54502bSHong Zhang PetscFunctionReturn(0); 3349b54502bSHong Zhang } 3359b54502bSHong Zhang EXTERN_C_END 3369b54502bSHong Zhang 3379b54502bSHong Zhang EXTERN_C_BEGIN 3389b54502bSHong Zhang #undef __FUNCT__ 3399b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering_LU" 3409b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering) 3419b54502bSHong Zhang { 3429b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3439b54502bSHong Zhang PetscErrorCode ierr; 3449b54502bSHong Zhang 3459b54502bSHong Zhang PetscFunctionBegin; 3469b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3479b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 3489b54502bSHong Zhang PetscFunctionReturn(0); 3499b54502bSHong Zhang } 3509b54502bSHong Zhang EXTERN_C_END 3519b54502bSHong Zhang 3529b54502bSHong Zhang EXTERN_C_BEGIN 3539b54502bSHong Zhang #undef __FUNCT__ 3549b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting_LU" 3559b54502bSHong Zhang PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol) 3569b54502bSHong Zhang { 3579b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3589b54502bSHong Zhang 3599b54502bSHong Zhang PetscFunctionBegin; 3609b54502bSHong Zhang if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",dtcol); 3619b54502bSHong Zhang dir->info.dtcol = dtcol; 3629b54502bSHong Zhang PetscFunctionReturn(0); 3639b54502bSHong Zhang } 3649b54502bSHong Zhang EXTERN_C_END 3659b54502bSHong Zhang 3669b54502bSHong Zhang EXTERN_C_BEGIN 3679b54502bSHong Zhang #undef __FUNCT__ 3689b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks_LU" 3699b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 3709b54502bSHong Zhang { 3719b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3729b54502bSHong Zhang 3739b54502bSHong Zhang PetscFunctionBegin; 3749b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 3759b54502bSHong Zhang PetscFunctionReturn(0); 3769b54502bSHong Zhang } 3779b54502bSHong Zhang EXTERN_C_END 3789b54502bSHong Zhang 3799b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3809b54502bSHong Zhang 3819b54502bSHong Zhang #undef __FUNCT__ 3829b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal" 3839b54502bSHong Zhang /*@ 3849b54502bSHong Zhang PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 3859b54502bSHong Zhang 3869b54502bSHong Zhang Collective on PC 3879b54502bSHong Zhang 3889b54502bSHong Zhang Input Parameters: 3899b54502bSHong Zhang + pc - the preconditioner context 3909b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 3919b54502bSHong Zhang 3929b54502bSHong Zhang Options Database Key: 3939b54502bSHong Zhang . -pc_lu_nonzeros_along_diagonal 3949b54502bSHong Zhang 3959b54502bSHong Zhang Level: intermediate 3969b54502bSHong Zhang 3979b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 3989b54502bSHong Zhang 399ee45ca4aSHong Zhang .seealso: PCLUSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 4009b54502bSHong Zhang @*/ 4019b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 4029b54502bSHong Zhang { 4039b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4049b54502bSHong Zhang 4059b54502bSHong Zhang PetscFunctionBegin; 4069b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4079b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 4089b54502bSHong Zhang if (f) { 4099b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 4109b54502bSHong Zhang } 4119b54502bSHong Zhang PetscFunctionReturn(0); 4129b54502bSHong Zhang } 4139b54502bSHong Zhang 4149b54502bSHong Zhang #undef __FUNCT__ 4159b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering" 4169b54502bSHong Zhang /*@ 4179b54502bSHong Zhang PCLUSetReuseOrdering - When similar matrices are factored, this 4189b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 4199b54502bSHong Zhang following factors; applies to both fill and drop tolerance LUs. 4209b54502bSHong Zhang 4219b54502bSHong Zhang Collective on PC 4229b54502bSHong Zhang 4239b54502bSHong Zhang Input Parameters: 4249b54502bSHong Zhang + pc - the preconditioner context 4259b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 4269b54502bSHong Zhang 4279b54502bSHong Zhang Options Database Key: 4289b54502bSHong Zhang . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 4299b54502bSHong Zhang 4309b54502bSHong Zhang Level: intermediate 4319b54502bSHong Zhang 4329b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 4339b54502bSHong Zhang 4349b54502bSHong Zhang .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill() 4359b54502bSHong Zhang @*/ 4369b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag) 4379b54502bSHong Zhang { 4389b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4399b54502bSHong Zhang 4409b54502bSHong Zhang PetscFunctionBegin; 4419b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4429b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 4439b54502bSHong Zhang if (f) { 4449b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 4459b54502bSHong Zhang } 4469b54502bSHong Zhang PetscFunctionReturn(0); 4479b54502bSHong Zhang } 4489b54502bSHong Zhang 4499b54502bSHong Zhang #undef __FUNCT__ 4509b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill" 4519b54502bSHong Zhang /*@ 4529b54502bSHong Zhang PCLUSetReuseFill - When matrices with same nonzero structure are LU factored, 4539b54502bSHong Zhang this causes later ones to use the fill computed in the initial factorization. 4549b54502bSHong Zhang 4559b54502bSHong Zhang Collective on PC 4569b54502bSHong Zhang 4579b54502bSHong Zhang Input Parameters: 4589b54502bSHong Zhang + pc - the preconditioner context 4599b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 4609b54502bSHong Zhang 4619b54502bSHong Zhang Options Database Key: 4629b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 4639b54502bSHong Zhang 4649b54502bSHong Zhang Level: intermediate 4659b54502bSHong Zhang 4669b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 4679b54502bSHong Zhang 4689b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill() 4699b54502bSHong Zhang @*/ 4709b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag) 4719b54502bSHong Zhang { 4729b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4739b54502bSHong Zhang 4749b54502bSHong Zhang PetscFunctionBegin; 4759b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4769b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 4779b54502bSHong Zhang if (f) { 4789b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 4799b54502bSHong Zhang } 4809b54502bSHong Zhang PetscFunctionReturn(0); 4819b54502bSHong Zhang } 4829b54502bSHong Zhang 4839b54502bSHong Zhang #undef __FUNCT__ 4849b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill" 4859b54502bSHong Zhang /*@ 4869b54502bSHong Zhang PCLUSetFill - Indicate the amount of fill you expect in the factored matrix, 4879b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 4889b54502bSHong Zhang 4899b54502bSHong Zhang Collective on PC 4909b54502bSHong Zhang 4919b54502bSHong Zhang Input Parameters: 4929b54502bSHong Zhang + pc - the preconditioner context 4939b54502bSHong Zhang - fill - amount of expected fill 4949b54502bSHong Zhang 4959b54502bSHong Zhang Options Database Key: 4969b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 4979b54502bSHong Zhang 4989b54502bSHong Zhang Level: intermediate 4999b54502bSHong Zhang 5009b54502bSHong Zhang Note: 5019b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 5029b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 5039b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 5049b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 5059b54502bSHong Zhang 5069b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 5079b54502bSHong Zhang 5089b54502bSHong Zhang .seealso: PCILUSetFill() 5099b54502bSHong Zhang @*/ 5109b54502bSHong Zhang PetscErrorCode PCLUSetFill(PC pc,PetscReal fill) 5119b54502bSHong Zhang { 5129b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5139b54502bSHong Zhang 5149b54502bSHong Zhang PetscFunctionBegin; 5159b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5169b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 5179b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5189b54502bSHong Zhang if (f) { 5199b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 5209b54502bSHong Zhang } 5219b54502bSHong Zhang PetscFunctionReturn(0); 5229b54502bSHong Zhang } 5239b54502bSHong Zhang 5249b54502bSHong Zhang #undef __FUNCT__ 5259b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace" 5269b54502bSHong Zhang /*@ 5279b54502bSHong Zhang PCLUSetUseInPlace - Tells the system to do an in-place factorization. 5289b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 5299b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 5309b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 5319b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 5329b54502bSHong Zhang space. 5339b54502bSHong Zhang 5349b54502bSHong Zhang Collective on PC 5359b54502bSHong Zhang 5369b54502bSHong Zhang Input Parameters: 5379b54502bSHong Zhang . pc - the preconditioner context 5389b54502bSHong Zhang 5399b54502bSHong Zhang Options Database Key: 5409b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 5419b54502bSHong Zhang 5429b54502bSHong Zhang Notes: 5439b54502bSHong Zhang PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 5449b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 5459b54502bSHong Zhang a call to KSPSetOperators(). 5469b54502bSHong Zhang This is because the Krylov space methods require an application of the 5479b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 5489b54502bSHong Zhang been factored in-place, replacing the original matrix. 5499b54502bSHong Zhang 5509b54502bSHong Zhang Level: intermediate 5519b54502bSHong Zhang 5529b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 5539b54502bSHong Zhang 5549b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 5559b54502bSHong Zhang @*/ 5569b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace(PC pc) 5579b54502bSHong Zhang { 5589b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 5599b54502bSHong Zhang 5609b54502bSHong Zhang PetscFunctionBegin; 5619b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5629b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 5639b54502bSHong Zhang if (f) { 5649b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 5659b54502bSHong Zhang } 5669b54502bSHong Zhang PetscFunctionReturn(0); 5679b54502bSHong Zhang } 5689b54502bSHong Zhang 5699b54502bSHong Zhang #undef __FUNCT__ 5709b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering" 5719b54502bSHong Zhang /*@C 5729b54502bSHong Zhang PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 5739b54502bSHong Zhang be used in the LU factorization. 5749b54502bSHong Zhang 5759b54502bSHong Zhang Collective on PC 5769b54502bSHong Zhang 5779b54502bSHong Zhang Input Parameters: 5789b54502bSHong Zhang + pc - the preconditioner context 5799b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 5809b54502bSHong Zhang 5819b54502bSHong Zhang Options Database Key: 5829b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 5839b54502bSHong Zhang 5849b54502bSHong Zhang Level: intermediate 5859b54502bSHong Zhang 5869b54502bSHong Zhang Notes: nested dissection is used by default 5879b54502bSHong Zhang 5889b54502bSHong Zhang .seealso: PCILUSetMatOrdering() 5899b54502bSHong Zhang @*/ 5909b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering) 5919b54502bSHong Zhang { 5929b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 5939b54502bSHong Zhang 5949b54502bSHong Zhang PetscFunctionBegin; 5959b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 5969b54502bSHong Zhang if (f) { 5979b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 5989b54502bSHong Zhang } 5999b54502bSHong Zhang PetscFunctionReturn(0); 6009b54502bSHong Zhang } 6019b54502bSHong Zhang 6029b54502bSHong Zhang #undef __FUNCT__ 6039b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting" 6049b54502bSHong Zhang /*@ 6059b54502bSHong Zhang PCLUSetPivoting - Determines when pivoting is done during LU. 6069b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 6079b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 6089b54502bSHong Zhang 6099b54502bSHong Zhang Collective on PC 6109b54502bSHong Zhang 6119b54502bSHong Zhang Input Parameters: 6129b54502bSHong Zhang + pc - the preconditioner context 6139b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 6149b54502bSHong Zhang 6159b54502bSHong Zhang Options Database Key: 6169b54502bSHong Zhang . -pc_lu_pivoting <dtcol> 6179b54502bSHong Zhang 6189b54502bSHong Zhang Level: intermediate 6199b54502bSHong Zhang 6209b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks() 6219b54502bSHong Zhang @*/ 6229b54502bSHong Zhang PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol) 6239b54502bSHong Zhang { 6249b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 6259b54502bSHong Zhang 6269b54502bSHong Zhang PetscFunctionBegin; 6279b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 6289b54502bSHong Zhang if (f) { 6299b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 6309b54502bSHong Zhang } 6319b54502bSHong Zhang PetscFunctionReturn(0); 6329b54502bSHong Zhang } 6339b54502bSHong Zhang 6349b54502bSHong Zhang #undef __FUNCT__ 6359b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks" 6369b54502bSHong Zhang /*@ 6379b54502bSHong Zhang PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block 6389b54502bSHong Zhang with BAIJ or SBAIJ matrices 6399b54502bSHong Zhang 6409b54502bSHong Zhang Collective on PC 6419b54502bSHong Zhang 6429b54502bSHong Zhang Input Parameters: 6439b54502bSHong Zhang + pc - the preconditioner context 6449b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 6459b54502bSHong Zhang 6469b54502bSHong Zhang Options Database Key: 6479b54502bSHong Zhang . -pc_lu_pivot_in_blocks <true,false> 6489b54502bSHong Zhang 6499b54502bSHong Zhang Level: intermediate 6509b54502bSHong Zhang 6519b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivoting() 6529b54502bSHong Zhang @*/ 6539b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot) 6549b54502bSHong Zhang { 6559b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 6569b54502bSHong Zhang 6579b54502bSHong Zhang PetscFunctionBegin; 6589b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 6599b54502bSHong Zhang if (f) { 6609b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 6619b54502bSHong Zhang } 6629b54502bSHong Zhang PetscFunctionReturn(0); 6639b54502bSHong Zhang } 6649b54502bSHong Zhang 6659b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 6669b54502bSHong Zhang 6679b54502bSHong Zhang /*MC 6689b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 6699b54502bSHong Zhang 6709b54502bSHong Zhang Options Database Keys: 6719b54502bSHong Zhang + -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 6729b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 6739b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 6749b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 6759b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 6769b54502bSHong Zhang - -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 6779b54502bSHong Zhang stability of factorization. 6789b54502bSHong Zhang 6799b54502bSHong Zhang Notes: Not all options work for all matrix formats 6809b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 6819b54502bSHong Zhang algorithms 6829b54502bSHong Zhang 6839b54502bSHong Zhang Level: beginner 6849b54502bSHong Zhang 6859b54502bSHong Zhang Concepts: LU factorization, direct solver 6869b54502bSHong Zhang 6879b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 6889b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 6899b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 6909b54502bSHong Zhang 6919b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 6929b54502bSHong Zhang PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(), 693ee45ca4aSHong Zhang PCLUSetFill(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCFactorSetPivoting(), 6949b54502bSHong Zhang PCLUSetPivotingInBlocks() 6959b54502bSHong Zhang M*/ 6969b54502bSHong Zhang 6979b54502bSHong Zhang EXTERN_C_BEGIN 6989b54502bSHong Zhang #undef __FUNCT__ 6999b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 7009b54502bSHong Zhang PetscErrorCode PCCreate_LU(PC pc) 7019b54502bSHong Zhang { 7029b54502bSHong Zhang PetscErrorCode ierr; 7039b54502bSHong Zhang PetscMPIInt size; 7049b54502bSHong Zhang PC_LU *dir; 7059b54502bSHong Zhang 7069b54502bSHong Zhang PetscFunctionBegin; 7079b54502bSHong Zhang ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr); 708*52e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_LU));CHKERRQ(ierr); 7099b54502bSHong Zhang 7109b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 7119b54502bSHong Zhang dir->fact = 0; 7129b54502bSHong Zhang dir->inplace = PETSC_FALSE; 7139b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 7149b54502bSHong Zhang 7159b54502bSHong Zhang dir->info.fill = 5.0; 7169b54502bSHong Zhang dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 7170a29876aSHong Zhang dir->info.shiftnz = 0.0; 7189b54502bSHong Zhang dir->info.zeropivot = 1.e-12; 7199b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 7200a29876aSHong Zhang dir->info.shiftpd = PETSC_FALSE; 7219b54502bSHong Zhang dir->info.shift_fraction = 0.0; 7229b54502bSHong Zhang dir->col = 0; 7239b54502bSHong Zhang dir->row = 0; 7249b54502bSHong Zhang ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 7259b54502bSHong Zhang if (size == 1) { 7269b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 7279b54502bSHong Zhang } else { 7289b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 7299b54502bSHong Zhang } 7309b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 7319b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 7329b54502bSHong Zhang pc->data = (void*)dir; 7339b54502bSHong Zhang 7349b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 7359b54502bSHong Zhang pc->ops->apply = PCApply_LU; 7369b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 7379b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 7389b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 7399b54502bSHong Zhang pc->ops->view = PCView_LU; 7409b54502bSHong Zhang pc->ops->applyrichardson = 0; 7419b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU; 7429b54502bSHong Zhang 743afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 744afaefe49SHong Zhang PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 745afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 746afaefe49SHong Zhang PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 747afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 748afaefe49SHong Zhang PCFactorSetShiftPd_LU);CHKERRQ(ierr); 749afaefe49SHong Zhang 7509b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU", 7519b54502bSHong Zhang PCLUSetFill_LU);CHKERRQ(ierr); 7529b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU", 7539b54502bSHong Zhang PCLUSetUseInPlace_LU);CHKERRQ(ierr); 7549b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU", 7559b54502bSHong Zhang PCLUSetMatOrdering_LU);CHKERRQ(ierr); 7569b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU", 7579b54502bSHong Zhang PCLUSetReuseOrdering_LU);CHKERRQ(ierr); 7589b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU", 7599b54502bSHong Zhang PCLUSetReuseFill_LU);CHKERRQ(ierr); 7609b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU", 7619b54502bSHong Zhang PCLUSetPivoting_LU);CHKERRQ(ierr); 7629b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU", 7639b54502bSHong Zhang PCLUSetPivotInBlocks_LU);CHKERRQ(ierr); 7649b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU", 7659b54502bSHong Zhang PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 7669b54502bSHong Zhang PetscFunctionReturn(0); 7679b54502bSHong Zhang } 7689b54502bSHong Zhang EXTERN_C_END 769