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__ 14afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_LU" 15dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_LU(PC pc,PetscReal z) 16afaefe49SHong Zhang { 17afaefe49SHong Zhang PC_LU *lu; 18afaefe49SHong Zhang 19afaefe49SHong Zhang PetscFunctionBegin; 20afaefe49SHong Zhang lu = (PC_LU*)pc->data; 21afaefe49SHong Zhang lu->info.zeropivot = z; 22afaefe49SHong Zhang PetscFunctionReturn(0); 23afaefe49SHong Zhang } 24afaefe49SHong Zhang EXTERN_C_END 25afaefe49SHong Zhang 26afaefe49SHong Zhang EXTERN_C_BEGIN 27afaefe49SHong Zhang #undef __FUNCT__ 28afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_LU" 29dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_LU(PC pc,PetscReal shift) 30afaefe49SHong Zhang { 31afaefe49SHong Zhang PC_LU *dir; 32afaefe49SHong Zhang 33afaefe49SHong Zhang PetscFunctionBegin; 34afaefe49SHong Zhang dir = (PC_LU*)pc->data; 35afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 36afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 37afaefe49SHong Zhang } else { 38afaefe49SHong Zhang dir->info.shiftnz = shift; 39afaefe49SHong Zhang } 40afaefe49SHong Zhang PetscFunctionReturn(0); 41afaefe49SHong Zhang } 42afaefe49SHong Zhang EXTERN_C_END 43afaefe49SHong Zhang 44afaefe49SHong Zhang EXTERN_C_BEGIN 45afaefe49SHong Zhang #undef __FUNCT__ 46afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_LU" 47dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_LU(PC pc,PetscTruth shift) 48afaefe49SHong Zhang { 49afaefe49SHong Zhang PC_LU *dir; 50afaefe49SHong Zhang 51afaefe49SHong Zhang PetscFunctionBegin; 52afaefe49SHong Zhang dir = (PC_LU*)pc->data; 53fbf22428SSatish Balay if (shift) { 54fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 55fbf22428SSatish Balay dir->info.shiftpd = 1.0; 56fbf22428SSatish Balay } else { 57fbf22428SSatish Balay dir->info.shiftpd = 0.0; 58fbf22428SSatish Balay } 59afaefe49SHong Zhang PetscFunctionReturn(0); 60afaefe49SHong Zhang } 61afaefe49SHong Zhang EXTERN_C_END 62afaefe49SHong Zhang 63afaefe49SHong Zhang EXTERN_C_BEGIN 64afaefe49SHong Zhang #undef __FUNCT__ 652401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 662401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 679b54502bSHong Zhang { 689b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 699b54502bSHong Zhang 709b54502bSHong Zhang PetscFunctionBegin; 719b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 729b54502bSHong Zhang if (z == PETSC_DECIDE) { 739b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 749b54502bSHong Zhang } else { 759b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 769b54502bSHong Zhang } 779b54502bSHong Zhang PetscFunctionReturn(0); 789b54502bSHong Zhang } 799b54502bSHong Zhang EXTERN_C_END 809b54502bSHong Zhang 819b54502bSHong Zhang EXTERN_C_BEGIN 829b54502bSHong Zhang #undef __FUNCT__ 832401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 842401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag) 859b54502bSHong Zhang { 869b54502bSHong Zhang PC_LU *lu; 879b54502bSHong Zhang 889b54502bSHong Zhang PetscFunctionBegin; 899b54502bSHong Zhang lu = (PC_LU*)pc->data; 909b54502bSHong Zhang lu->reuseordering = flag; 919b54502bSHong Zhang PetscFunctionReturn(0); 929b54502bSHong Zhang } 939b54502bSHong Zhang EXTERN_C_END 949b54502bSHong Zhang 959b54502bSHong Zhang EXTERN_C_BEGIN 969b54502bSHong Zhang #undef __FUNCT__ 972401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 982401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag) 999b54502bSHong Zhang { 1009b54502bSHong Zhang PC_LU *lu; 1019b54502bSHong Zhang 1029b54502bSHong Zhang PetscFunctionBegin; 1039b54502bSHong Zhang lu = (PC_LU*)pc->data; 1049b54502bSHong Zhang lu->reusefill = flag; 1059b54502bSHong Zhang PetscFunctionReturn(0); 1069b54502bSHong Zhang } 1079b54502bSHong Zhang EXTERN_C_END 1089b54502bSHong Zhang 1099b54502bSHong Zhang #undef __FUNCT__ 1109b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 1119b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 1129b54502bSHong Zhang { 1139b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1149b54502bSHong Zhang PetscErrorCode ierr; 1159b54502bSHong Zhang PetscTruth flg,set; 1169b54502bSHong Zhang char tname[256]; 1179b54502bSHong Zhang PetscFList ordlist; 1189b54502bSHong Zhang PetscReal tol; 1199b54502bSHong Zhang 1209b54502bSHong Zhang PetscFunctionBegin; 1219b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1229b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 1232401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_in_place","Form LU in the same memory as the matrix","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 1249b54502bSHong Zhang if (flg) { 1252401956bSBarry Smith ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr); 1269b54502bSHong Zhang } 12755ba2a51SBarry 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); 1289b54502bSHong Zhang 1299f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 1309b54502bSHong Zhang if (flg) { 131afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1329b54502bSHong Zhang } 1339f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 1349f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 1359b54502bSHong Zhang if (flg) { 136afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 1379b54502bSHong Zhang } 138ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 1399b54502bSHong Zhang 1402401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 1419b54502bSHong Zhang if (flg) { 1422401956bSBarry Smith ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1439b54502bSHong Zhang } 1442401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 1459b54502bSHong Zhang if (flg) { 1462401956bSBarry Smith ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1479b54502bSHong Zhang } 1489b54502bSHong Zhang 1499b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 150e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in LU","PCFactorSetMatOrderingType",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 1519b54502bSHong Zhang if (flg) { 152e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 1539b54502bSHong Zhang } 1549b54502bSHong Zhang 1552401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1569b54502bSHong Zhang if (flg) { 1579b54502bSHong Zhang tol = PETSC_DECIDE; 1582401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1592401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1609b54502bSHong Zhang } 1619b54502bSHong Zhang 1622401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_pivoting","Pivoting tolerance (used only for some factorization)","PCFactorSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 1639b54502bSHong Zhang 1649b54502bSHong Zhang flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 1652401956bSBarry Smith ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 1669b54502bSHong Zhang if (set) { 1672401956bSBarry Smith ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 1689b54502bSHong Zhang } 1699b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1709b54502bSHong Zhang PetscFunctionReturn(0); 1719b54502bSHong Zhang } 1729b54502bSHong Zhang 1739b54502bSHong Zhang #undef __FUNCT__ 1749b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 1759b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 1769b54502bSHong Zhang { 1779b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1789b54502bSHong Zhang PetscErrorCode ierr; 1799b54502bSHong Zhang PetscTruth iascii,isstring; 1809b54502bSHong Zhang 1819b54502bSHong Zhang PetscFunctionBegin; 1829b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1839b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1849b54502bSHong Zhang if (iascii) { 1859b54502bSHong Zhang 1869b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 1879b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 1889b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 189a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",lu->info.zeropivot);CHKERRQ(ierr); 1900a29876aSHong Zhang if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 1919b54502bSHong Zhang if (lu->fact) { 192f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 193f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 194f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 195f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 196f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 197f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 1989b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 1999b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 200f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 201f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 202f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2039b54502bSHong Zhang } 2049b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 2059b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 2069b54502bSHong Zhang } else if (isstring) { 2079b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 2089b54502bSHong Zhang } else { 2099b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 2109b54502bSHong Zhang } 2119b54502bSHong Zhang PetscFunctionReturn(0); 2129b54502bSHong Zhang } 2139b54502bSHong Zhang 2149b54502bSHong Zhang #undef __FUNCT__ 215a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_LU" 216a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_LU(PC pc,Mat *mat) 2179b54502bSHong Zhang { 2189b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2199b54502bSHong Zhang 2209b54502bSHong Zhang PetscFunctionBegin; 2219b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 2229b54502bSHong Zhang *mat = dir->fact; 2239b54502bSHong Zhang PetscFunctionReturn(0); 2249b54502bSHong Zhang } 2259b54502bSHong Zhang 2269b54502bSHong Zhang #undef __FUNCT__ 2279b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 2289b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 2299b54502bSHong Zhang { 2309b54502bSHong Zhang PetscErrorCode ierr; 2319b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2329b54502bSHong Zhang 2339b54502bSHong Zhang PetscFunctionBegin; 2349b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 2359b54502bSHong Zhang 2369b54502bSHong Zhang if (dir->inplace) { 2379b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2389b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2399b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 24003c60df9SBarry Smith if (dir->row) { 24103c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 24203c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 24303c60df9SBarry Smith } 2449b54502bSHong Zhang ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 2459b54502bSHong Zhang dir->fact = pc->pmat; 2469b54502bSHong Zhang } else { 2479b54502bSHong Zhang MatInfo info; 2489b54502bSHong Zhang if (!pc->setupcalled) { 2499b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2509b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2519b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2529b54502bSHong Zhang } 25303c60df9SBarry Smith if (dir->row) { 25403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 25503c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 25603c60df9SBarry Smith } 2579b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2589b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2599b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 26052e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2619b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2629b54502bSHong Zhang if (!dir->reuseordering) { 2639b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2649b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2659b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2669b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2679b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2689b54502bSHong Zhang } 26903c60df9SBarry Smith if (dir->row) { 27003c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 27103c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 27203c60df9SBarry Smith } 2739b54502bSHong Zhang } 2749b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 2759b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2769b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2779b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 27852e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,dir->fact);CHKERRQ(ierr); 2799b54502bSHong Zhang } 2809b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 2819b54502bSHong Zhang } 2829b54502bSHong Zhang PetscFunctionReturn(0); 2839b54502bSHong Zhang } 2849b54502bSHong Zhang 2859b54502bSHong Zhang #undef __FUNCT__ 2869b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 2879b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 2889b54502bSHong Zhang { 2899b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2909b54502bSHong Zhang PetscErrorCode ierr; 2919b54502bSHong Zhang 2929b54502bSHong Zhang PetscFunctionBegin; 2939b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 2949b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2959b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2969b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 2979b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 2989b54502bSHong Zhang PetscFunctionReturn(0); 2999b54502bSHong Zhang } 3009b54502bSHong Zhang 3019b54502bSHong Zhang #undef __FUNCT__ 3029b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 3039b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 3049b54502bSHong Zhang { 3059b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3069b54502bSHong Zhang PetscErrorCode ierr; 3079b54502bSHong Zhang 3089b54502bSHong Zhang PetscFunctionBegin; 3099b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 3109b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 3119b54502bSHong Zhang PetscFunctionReturn(0); 3129b54502bSHong Zhang } 3139b54502bSHong Zhang 3149b54502bSHong Zhang #undef __FUNCT__ 3159b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 3169b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 3179b54502bSHong Zhang { 3189b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3199b54502bSHong Zhang PetscErrorCode ierr; 3209b54502bSHong Zhang 3219b54502bSHong Zhang PetscFunctionBegin; 3229b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 3239b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 3249b54502bSHong Zhang PetscFunctionReturn(0); 3259b54502bSHong Zhang } 3269b54502bSHong Zhang 3279b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3289b54502bSHong Zhang 3299b54502bSHong Zhang EXTERN_C_BEGIN 3309b54502bSHong Zhang #undef __FUNCT__ 33155ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_LU" 33255ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_LU(PC pc,PetscReal fill) 3339b54502bSHong Zhang { 3349b54502bSHong Zhang PC_LU *dir; 3359b54502bSHong Zhang 3369b54502bSHong Zhang PetscFunctionBegin; 3379b54502bSHong Zhang dir = (PC_LU*)pc->data; 3389b54502bSHong Zhang dir->info.fill = fill; 3399b54502bSHong Zhang PetscFunctionReturn(0); 3409b54502bSHong Zhang } 3419b54502bSHong Zhang EXTERN_C_END 3429b54502bSHong Zhang 3439b54502bSHong Zhang EXTERN_C_BEGIN 3449b54502bSHong Zhang #undef __FUNCT__ 3452401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 3462401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 3479b54502bSHong Zhang { 3489b54502bSHong Zhang PC_LU *dir; 3499b54502bSHong Zhang 3509b54502bSHong Zhang PetscFunctionBegin; 3519b54502bSHong Zhang dir = (PC_LU*)pc->data; 3529b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3539b54502bSHong Zhang PetscFunctionReturn(0); 3549b54502bSHong Zhang } 3559b54502bSHong Zhang EXTERN_C_END 3569b54502bSHong Zhang 3579b54502bSHong Zhang EXTERN_C_BEGIN 3589b54502bSHong Zhang #undef __FUNCT__ 359e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_LU" 360e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_LU(PC pc,MatOrderingType ordering) 3619b54502bSHong Zhang { 3629b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3639b54502bSHong Zhang PetscErrorCode ierr; 3649b54502bSHong Zhang 3659b54502bSHong Zhang PetscFunctionBegin; 3669b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3679b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 3689b54502bSHong Zhang PetscFunctionReturn(0); 3699b54502bSHong Zhang } 3709b54502bSHong Zhang EXTERN_C_END 3719b54502bSHong Zhang 3729b54502bSHong Zhang EXTERN_C_BEGIN 3739b54502bSHong Zhang #undef __FUNCT__ 3742401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting_LU" 3752401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting_LU(PC pc,PetscReal dtcol) 3769b54502bSHong Zhang { 3779b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3789b54502bSHong Zhang 3799b54502bSHong Zhang PetscFunctionBegin; 380a83599f4SBarry 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); 3819b54502bSHong Zhang dir->info.dtcol = dtcol; 3829b54502bSHong Zhang PetscFunctionReturn(0); 3839b54502bSHong Zhang } 3849b54502bSHong Zhang EXTERN_C_END 3859b54502bSHong Zhang 3869b54502bSHong Zhang EXTERN_C_BEGIN 3879b54502bSHong Zhang #undef __FUNCT__ 3882401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_LU" 3892401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 3909b54502bSHong Zhang { 3919b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3929b54502bSHong Zhang 3939b54502bSHong Zhang PetscFunctionBegin; 3949b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 3959b54502bSHong Zhang PetscFunctionReturn(0); 3969b54502bSHong Zhang } 3979b54502bSHong Zhang EXTERN_C_END 3989b54502bSHong Zhang 3999b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 4009b54502bSHong Zhang 4019b54502bSHong Zhang #undef __FUNCT__ 4022401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal" 4039b54502bSHong Zhang /*@ 4042401956bSBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 4059b54502bSHong Zhang 4069b54502bSHong Zhang Collective on PC 4079b54502bSHong Zhang 4089b54502bSHong Zhang Input Parameters: 4099b54502bSHong Zhang + pc - the preconditioner context 4109b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 4119b54502bSHong Zhang 4129b54502bSHong Zhang Options Database Key: 4132401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal 4149b54502bSHong Zhang 4159b54502bSHong Zhang Level: intermediate 4169b54502bSHong Zhang 4179b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4189b54502bSHong Zhang 41955ba2a51SBarry Smith .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal() 4209b54502bSHong Zhang @*/ 4212401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 4229b54502bSHong Zhang { 4239b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4249b54502bSHong Zhang 4259b54502bSHong Zhang PetscFunctionBegin; 4269b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4272401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 4289b54502bSHong Zhang if (f) { 4299b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 4309b54502bSHong Zhang } 4319b54502bSHong Zhang PetscFunctionReturn(0); 4329b54502bSHong Zhang } 4339b54502bSHong Zhang 4349b54502bSHong Zhang #undef __FUNCT__ 43555ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill" 4369b54502bSHong Zhang /*@ 43755ba2a51SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 4389b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 4399b54502bSHong Zhang 4409b54502bSHong Zhang Collective on PC 4419b54502bSHong Zhang 4429b54502bSHong Zhang Input Parameters: 4439b54502bSHong Zhang + pc - the preconditioner context 4449b54502bSHong Zhang - fill - amount of expected fill 4459b54502bSHong Zhang 4469b54502bSHong Zhang Options Database Key: 44755ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 4489b54502bSHong Zhang 4499b54502bSHong Zhang Level: intermediate 4509b54502bSHong Zhang 4519b54502bSHong Zhang Note: 4529b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 4536cf91177SBarry Smith fill to expect. By running with the option -info PETSc will print the 4549b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 4559b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 4569b54502bSHong Zhang 4579b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4589b54502bSHong Zhang 4599b54502bSHong Zhang @*/ 46055ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill(PC pc,PetscReal fill) 4619b54502bSHong Zhang { 4629b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4639b54502bSHong Zhang 4649b54502bSHong Zhang PetscFunctionBegin; 4659b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4669b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 46755ba2a51SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 4689b54502bSHong Zhang if (f) { 4699b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 4709b54502bSHong Zhang } 4719b54502bSHong Zhang PetscFunctionReturn(0); 4729b54502bSHong Zhang } 4739b54502bSHong Zhang 4749b54502bSHong Zhang #undef __FUNCT__ 4752401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace" 4769b54502bSHong Zhang /*@ 4772401956bSBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 4789b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 4799b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 4809b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 4819b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 4829b54502bSHong Zhang space. 4839b54502bSHong Zhang 4849b54502bSHong Zhang Collective on PC 4859b54502bSHong Zhang 4869b54502bSHong Zhang Input Parameters: 4879b54502bSHong Zhang . pc - the preconditioner context 4889b54502bSHong Zhang 4899b54502bSHong Zhang Options Database Key: 4902401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 4919b54502bSHong Zhang 4929b54502bSHong Zhang Notes: 4932401956bSBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 4949b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 4959b54502bSHong Zhang a call to KSPSetOperators(). 4969b54502bSHong Zhang This is because the Krylov space methods require an application of the 4979b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 4989b54502bSHong Zhang been factored in-place, replacing the original matrix. 4999b54502bSHong Zhang 5009b54502bSHong Zhang Level: intermediate 5019b54502bSHong Zhang 5029b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 5039b54502bSHong Zhang 5049b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 5059b54502bSHong Zhang @*/ 5062401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace(PC pc) 5079b54502bSHong Zhang { 5089b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 5099b54502bSHong Zhang 5109b54502bSHong Zhang PetscFunctionBegin; 5119b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5122401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 5139b54502bSHong Zhang if (f) { 5149b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 5159b54502bSHong Zhang } 5169b54502bSHong Zhang PetscFunctionReturn(0); 5179b54502bSHong Zhang } 5189b54502bSHong Zhang 5199b54502bSHong Zhang #undef __FUNCT__ 520e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType" 5219b54502bSHong Zhang /*@C 522e5a9bf91SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 5239b54502bSHong Zhang be used in the LU factorization. 5249b54502bSHong Zhang 5259b54502bSHong Zhang Collective on PC 5269b54502bSHong Zhang 5279b54502bSHong Zhang Input Parameters: 5289b54502bSHong Zhang + pc - the preconditioner context 5299b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 5309b54502bSHong Zhang 5319b54502bSHong Zhang Options Database Key: 5322401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 5339b54502bSHong Zhang 5349b54502bSHong Zhang Level: intermediate 5359b54502bSHong Zhang 5369b54502bSHong Zhang Notes: nested dissection is used by default 5379b54502bSHong Zhang 538*2b6de112SBarry Smith For Cholesky and ICC and the SBAIJ format reorderings are not available, 539*2b6de112SBarry Smith since only the upper triangular part of the matrix is stored. You can use the 540*2b6de112SBarry Smith SeqAIJ format in this case to get reorderings. 541*2b6de112SBarry Smith 5429b54502bSHong Zhang @*/ 543e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 5449b54502bSHong Zhang { 5459b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 5469b54502bSHong Zhang 5479b54502bSHong Zhang PetscFunctionBegin; 548e5a9bf91SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);CHKERRQ(ierr); 5499b54502bSHong Zhang if (f) { 5509b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 5519b54502bSHong Zhang } 5529b54502bSHong Zhang PetscFunctionReturn(0); 5539b54502bSHong Zhang } 5549b54502bSHong Zhang 5559b54502bSHong Zhang #undef __FUNCT__ 5562401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivoting" 5579b54502bSHong Zhang /*@ 5582401956bSBarry Smith PCFactorSetPivoting - Determines when pivoting is done during LU. 5599b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 5609b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 5619b54502bSHong Zhang 5629b54502bSHong Zhang Collective on PC 5639b54502bSHong Zhang 5649b54502bSHong Zhang Input Parameters: 5659b54502bSHong Zhang + pc - the preconditioner context 5669b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 5679b54502bSHong Zhang 5689b54502bSHong Zhang Options Database Key: 5692401956bSBarry Smith . -pc_factor_pivoting <dtcol> 5709b54502bSHong Zhang 5719b54502bSHong Zhang Level: intermediate 5729b54502bSHong Zhang 5732401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks() 5749b54502bSHong Zhang @*/ 5752401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivoting(PC pc,PetscReal dtcol) 5769b54502bSHong Zhang { 5779b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5789b54502bSHong Zhang 5799b54502bSHong Zhang PetscFunctionBegin; 5802401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 5819b54502bSHong Zhang if (f) { 5829b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 5839b54502bSHong Zhang } 5849b54502bSHong Zhang PetscFunctionReturn(0); 5859b54502bSHong Zhang } 5869b54502bSHong Zhang 5879b54502bSHong Zhang #undef __FUNCT__ 5882401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks" 5899b54502bSHong Zhang /*@ 5902401956bSBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 5919b54502bSHong Zhang with BAIJ or SBAIJ matrices 5929b54502bSHong Zhang 5939b54502bSHong Zhang Collective on PC 5949b54502bSHong Zhang 5959b54502bSHong Zhang Input Parameters: 5969b54502bSHong Zhang + pc - the preconditioner context 5979b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 5989b54502bSHong Zhang 5999b54502bSHong Zhang Options Database Key: 6002401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> 6019b54502bSHong Zhang 6029b54502bSHong Zhang Level: intermediate 6039b54502bSHong Zhang 6042401956bSBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting() 6059b54502bSHong Zhang @*/ 6062401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot) 6079b54502bSHong Zhang { 6089b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 6099b54502bSHong Zhang 6109b54502bSHong Zhang PetscFunctionBegin; 6112401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 6129b54502bSHong Zhang if (f) { 6139b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 6149b54502bSHong Zhang } 6159b54502bSHong Zhang PetscFunctionReturn(0); 6169b54502bSHong Zhang } 6179b54502bSHong Zhang 6189b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 6199b54502bSHong Zhang 6209b54502bSHong Zhang /*MC 6219b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 6229b54502bSHong Zhang 6239b54502bSHong Zhang Options Database Keys: 6242401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 6252401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 62655ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 6272401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 6282401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 6292401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 6309b54502bSHong Zhang stability of factorization. 631f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 632e22d95b2SBarry Smith . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 633f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 634e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 635e22d95b2SBarry Smith diagonal. 6369b54502bSHong Zhang 6379b54502bSHong Zhang Notes: Not all options work for all matrix formats 6389b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 6399b54502bSHong Zhang algorithms 6409b54502bSHong Zhang 6419b54502bSHong Zhang Level: beginner 6429b54502bSHong Zhang 6439b54502bSHong Zhang Concepts: LU factorization, direct solver 6449b54502bSHong Zhang 6459b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 6469b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 6479b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 6489b54502bSHong Zhang 6499b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 650a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 651e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetPivoting(), 652e22d95b2SBarry Smith PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), PCFactorReorderForNonzeroDiagonal() 6539b54502bSHong Zhang M*/ 6549b54502bSHong Zhang 6559b54502bSHong Zhang EXTERN_C_BEGIN 6569b54502bSHong Zhang #undef __FUNCT__ 6579b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 658dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 6599b54502bSHong Zhang { 6609b54502bSHong Zhang PetscErrorCode ierr; 6619b54502bSHong Zhang PetscMPIInt size; 6629b54502bSHong Zhang PC_LU *dir; 6639b54502bSHong Zhang 6649b54502bSHong Zhang PetscFunctionBegin; 66538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 6669b54502bSHong Zhang 6679b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 6689b54502bSHong Zhang dir->fact = 0; 6699b54502bSHong Zhang dir->inplace = PETSC_FALSE; 6709b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 6719b54502bSHong Zhang 6729b54502bSHong Zhang dir->info.fill = 5.0; 6739b54502bSHong Zhang dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 6740a29876aSHong Zhang dir->info.shiftnz = 0.0; 6759b54502bSHong Zhang dir->info.zeropivot = 1.e-12; 6769b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 677fbf22428SSatish Balay dir->info.shiftpd = 0.0; /* false */ 6789b54502bSHong Zhang dir->info.shift_fraction = 0.0; 6799b54502bSHong Zhang dir->col = 0; 6809b54502bSHong Zhang dir->row = 0; 6817adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 6829b54502bSHong Zhang if (size == 1) { 6839b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 6849b54502bSHong Zhang } else { 6859b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 6869b54502bSHong Zhang } 6879b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 6889b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 6899b54502bSHong Zhang pc->data = (void*)dir; 6909b54502bSHong Zhang 6919b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 6929b54502bSHong Zhang pc->ops->apply = PCApply_LU; 6939b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 6949b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 6959b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 6969b54502bSHong Zhang pc->ops->view = PCView_LU; 6979b54502bSHong Zhang pc->ops->applyrichardson = 0; 698a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_LU; 6999b54502bSHong Zhang 700afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_LU", 701afaefe49SHong Zhang PCFactorSetZeroPivot_LU);CHKERRQ(ierr); 702afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_LU", 703afaefe49SHong Zhang PCFactorSetShiftNonzero_LU);CHKERRQ(ierr); 704afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_LU", 705afaefe49SHong Zhang PCFactorSetShiftPd_LU);CHKERRQ(ierr); 706afaefe49SHong Zhang 70755ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_LU", 70855ba2a51SBarry Smith PCFactorSetFill_LU);CHKERRQ(ierr); 7092401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 7102401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 711e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_LU", 712e5a9bf91SBarry Smith PCFactorSetMatOrderingType_LU);CHKERRQ(ierr); 7132401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 7142401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 7152401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 7162401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 7172401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivoting_C","PCFactorSetPivoting_LU", 7182401956bSBarry Smith PCFactorSetPivoting_LU);CHKERRQ(ierr); 7192401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_LU", 7202401956bSBarry Smith PCFactorSetPivotInBlocks_LU);CHKERRQ(ierr); 7212401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 7222401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 7239b54502bSHong Zhang PetscFunctionReturn(0); 7249b54502bSHong Zhang } 7259b54502bSHong Zhang EXTERN_C_END 726