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 */ 69b54502bSHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 79b54502bSHong Zhang 89b54502bSHong Zhang typedef struct { 99b54502bSHong Zhang Mat fact; /* factored matrix */ 109b54502bSHong Zhang PetscReal actualfill; /* actual fill in factor */ 119b54502bSHong Zhang PetscTruth inplace; /* flag indicating in-place factorization */ 129b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 139b54502bSHong Zhang MatOrderingType ordering; /* matrix ordering */ 149b54502bSHong Zhang PetscTruth reuseordering; /* reuses previous reordering computed */ 159b54502bSHong Zhang PetscTruth reusefill; /* reuse fill from previous LU */ 169b54502bSHong Zhang MatFactorInfo info; 179b54502bSHong Zhang PetscTruth nonzerosalongdiagonal; 189b54502bSHong Zhang PetscReal nonzerosalongdiagonaltol; 199b54502bSHong Zhang } PC_LU; 209b54502bSHong Zhang 219b54502bSHong Zhang 229b54502bSHong Zhang EXTERN_C_BEGIN 239b54502bSHong Zhang #undef __FUNCT__ 249b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal_LU" 259b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 269b54502bSHong Zhang { 279b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 289b54502bSHong Zhang 299b54502bSHong Zhang PetscFunctionBegin; 309b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 319b54502bSHong Zhang if (z == PETSC_DECIDE) { 329b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 339b54502bSHong Zhang } else { 349b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 359b54502bSHong Zhang } 369b54502bSHong Zhang PetscFunctionReturn(0); 379b54502bSHong Zhang } 389b54502bSHong Zhang EXTERN_C_END 399b54502bSHong Zhang 409b54502bSHong Zhang EXTERN_C_BEGIN 419b54502bSHong Zhang #undef __FUNCT__ 429b54502bSHong Zhang #define __FUNCT__ "PCLUSetSetZeroPivot" 439b54502bSHong Zhang PetscErrorCode PCLUSetZeroPivot_LU(PC pc,PetscReal z) 449b54502bSHong Zhang { 459b54502bSHong Zhang PC_LU *lu; 469b54502bSHong Zhang 479b54502bSHong Zhang PetscFunctionBegin; 489b54502bSHong Zhang lu = (PC_LU*)pc->data; 499b54502bSHong Zhang lu->info.zeropivot = z; 509b54502bSHong Zhang PetscFunctionReturn(0); 519b54502bSHong Zhang } 529b54502bSHong Zhang EXTERN_C_END 539b54502bSHong Zhang 549b54502bSHong Zhang EXTERN_C_BEGIN 559b54502bSHong Zhang #undef __FUNCT__ 569b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering_LU" 579b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag) 589b54502bSHong Zhang { 599b54502bSHong Zhang PC_LU *lu; 609b54502bSHong Zhang 619b54502bSHong Zhang PetscFunctionBegin; 629b54502bSHong Zhang lu = (PC_LU*)pc->data; 639b54502bSHong Zhang lu->reuseordering = flag; 649b54502bSHong Zhang PetscFunctionReturn(0); 659b54502bSHong Zhang } 669b54502bSHong Zhang EXTERN_C_END 679b54502bSHong Zhang 689b54502bSHong Zhang EXTERN_C_BEGIN 699b54502bSHong Zhang #undef __FUNCT__ 709b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill_LU" 719b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag) 729b54502bSHong Zhang { 739b54502bSHong Zhang PC_LU *lu; 749b54502bSHong Zhang 759b54502bSHong Zhang PetscFunctionBegin; 769b54502bSHong Zhang lu = (PC_LU*)pc->data; 779b54502bSHong Zhang lu->reusefill = flag; 789b54502bSHong Zhang PetscFunctionReturn(0); 799b54502bSHong Zhang } 809b54502bSHong Zhang EXTERN_C_END 819b54502bSHong Zhang 829b54502bSHong Zhang EXTERN_C_BEGIN 839b54502bSHong Zhang #undef __FUNCT__ 849b54502bSHong Zhang #define __FUNCT__ "PCLUSetShift_LU" 859b54502bSHong Zhang PetscErrorCode PCLUSetShift_LU(PC pc,PetscTruth shift) 869b54502bSHong Zhang { 879b54502bSHong Zhang PC_LU *dir; 889b54502bSHong Zhang 899b54502bSHong Zhang PetscFunctionBegin; 909b54502bSHong Zhang dir = (PC_LU*)pc->data; 91*0a29876aSHong Zhang dir->info.shiftnz = shift; 929b54502bSHong Zhang if (shift) dir->info.shift_fraction = 0.0; 939b54502bSHong Zhang PetscFunctionReturn(0); 949b54502bSHong Zhang } 959b54502bSHong Zhang EXTERN_C_END 969b54502bSHong Zhang 979b54502bSHong Zhang #undef __FUNCT__ 989b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 999b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 1009b54502bSHong Zhang { 1019b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1029b54502bSHong Zhang PetscErrorCode ierr; 1039b54502bSHong Zhang PetscTruth flg,set; 1049b54502bSHong Zhang char tname[256]; 1059b54502bSHong Zhang PetscFList ordlist; 1069b54502bSHong Zhang PetscReal tol; 1079b54502bSHong Zhang 1089b54502bSHong Zhang PetscFunctionBegin; 1099b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1109b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 1119b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);CHKERRQ(ierr); 1129b54502bSHong Zhang if (flg) { 1139b54502bSHong Zhang ierr = PCLUSetUseInPlace(pc);CHKERRQ(ierr); 1149b54502bSHong Zhang } 1159b54502bSHong 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); 1169b54502bSHong Zhang 1179b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",&flg);CHKERRQ(ierr); 1189b54502bSHong Zhang if (flg) { 1199b54502bSHong Zhang ierr = PCLUSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr); 1209b54502bSHong Zhang } 121*0a29876aSHong Zhang ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);CHKERRQ(ierr); 1229b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_shift","Manteuffel shift applied to diagonal","PCLUSetShift",&flg);CHKERRQ(ierr); 1239b54502bSHong Zhang if (flg) { 1249b54502bSHong Zhang ierr = PCLUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 1259b54502bSHong Zhang } 1269b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);CHKERRQ(ierr); 1279b54502bSHong Zhang 1289b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);CHKERRQ(ierr); 1299b54502bSHong Zhang if (flg) { 1309b54502bSHong Zhang ierr = PCLUSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr); 1319b54502bSHong Zhang } 1329b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);CHKERRQ(ierr); 1339b54502bSHong Zhang if (flg) { 1349b54502bSHong Zhang ierr = PCLUSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr); 1359b54502bSHong Zhang } 1369b54502bSHong Zhang 1379b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 1389b54502bSHong Zhang ierr = PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);CHKERRQ(ierr); 1399b54502bSHong Zhang if (flg) { 1409b54502bSHong Zhang ierr = PCLUSetMatOrdering(pc,tname);CHKERRQ(ierr); 1419b54502bSHong Zhang } 1429b54502bSHong Zhang 1439b54502bSHong Zhang ierr = PetscOptionsName("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1449b54502bSHong Zhang if (flg) { 1459b54502bSHong Zhang tol = PETSC_DECIDE; 1469b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCLUReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1479b54502bSHong Zhang ierr = PCLUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1489b54502bSHong Zhang } 1499b54502bSHong Zhang 1509b54502bSHong Zhang ierr = PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);CHKERRQ(ierr); 1519b54502bSHong Zhang 1529b54502bSHong Zhang flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 1539b54502bSHong Zhang ierr = PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 1549b54502bSHong Zhang if (set) { 1559b54502bSHong Zhang ierr = PCLUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 1569b54502bSHong Zhang } 1579b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1589b54502bSHong Zhang PetscFunctionReturn(0); 1599b54502bSHong Zhang } 1609b54502bSHong Zhang 1619b54502bSHong Zhang #undef __FUNCT__ 1629b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 1639b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 1649b54502bSHong Zhang { 1659b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 1669b54502bSHong Zhang PetscErrorCode ierr; 1679b54502bSHong Zhang PetscTruth iascii,isstring; 1689b54502bSHong Zhang 1699b54502bSHong Zhang PetscFunctionBegin; 1709b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1719b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1729b54502bSHong Zhang if (iascii) { 1739b54502bSHong Zhang MatInfo info; 1749b54502bSHong Zhang 1759b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 1769b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 1779b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);CHKERRQ(ierr); 1789b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %g\n",lu->info.zeropivot);CHKERRQ(ierr); 179*0a29876aSHong Zhang if (lu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 1809b54502bSHong Zhang if (lu->fact) { 1819b54502bSHong Zhang ierr = MatGetInfo(lu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1829b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU nonzeros %g\n",info.nz_used);CHKERRQ(ierr); 1839b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);CHKERRQ(ierr); 1849b54502bSHong Zhang ierr = MatView(lu->fact,viewer);CHKERRQ(ierr); 1859b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1869b54502bSHong Zhang } 1879b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 1889b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1899b54502bSHong Zhang } else if (isstring) { 1909b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 1919b54502bSHong Zhang } else { 1929b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 1939b54502bSHong Zhang } 1949b54502bSHong Zhang PetscFunctionReturn(0); 1959b54502bSHong Zhang } 1969b54502bSHong Zhang 1979b54502bSHong Zhang #undef __FUNCT__ 1989b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_LU" 1999b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat) 2009b54502bSHong Zhang { 2019b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2029b54502bSHong Zhang 2039b54502bSHong Zhang PetscFunctionBegin; 2049b54502bSHong Zhang if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 2059b54502bSHong Zhang *mat = dir->fact; 2069b54502bSHong Zhang PetscFunctionReturn(0); 2079b54502bSHong Zhang } 2089b54502bSHong Zhang 2099b54502bSHong Zhang #undef __FUNCT__ 2109b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 2119b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 2129b54502bSHong Zhang { 2139b54502bSHong Zhang PetscErrorCode ierr; 2149b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2159b54502bSHong Zhang 2169b54502bSHong Zhang PetscFunctionBegin; 2179b54502bSHong Zhang if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill; 2189b54502bSHong Zhang 2199b54502bSHong Zhang if (dir->inplace) { 2209b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2219b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2229b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2239b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2249b54502bSHong Zhang ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);CHKERRQ(ierr); 2259b54502bSHong Zhang dir->fact = pc->pmat; 2269b54502bSHong Zhang } else { 2279b54502bSHong Zhang MatInfo info; 2289b54502bSHong Zhang if (!pc->setupcalled) { 2299b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2309b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2319b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2329b54502bSHong Zhang } 2339b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2349b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2359b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2369b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 2379b54502bSHong Zhang PetscLogObjectParent(pc,dir->fact); 2389b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2399b54502bSHong Zhang if (!dir->reuseordering) { 2409b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2419b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2429b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 2439b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 2449b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 2459b54502bSHong Zhang } 2469b54502bSHong Zhang if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);} 2479b54502bSHong Zhang } 2489b54502bSHong Zhang ierr = MatDestroy(dir->fact);CHKERRQ(ierr); 2499b54502bSHong Zhang ierr = MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);CHKERRQ(ierr); 2509b54502bSHong Zhang ierr = MatGetInfo(dir->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2519b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 2529b54502bSHong Zhang PetscLogObjectParent(pc,dir->fact); 2539b54502bSHong Zhang } 2549b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&dir->info,&dir->fact);CHKERRQ(ierr); 2559b54502bSHong Zhang } 2569b54502bSHong Zhang PetscFunctionReturn(0); 2579b54502bSHong Zhang } 2589b54502bSHong Zhang 2599b54502bSHong Zhang #undef __FUNCT__ 2609b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 2619b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 2629b54502bSHong Zhang { 2639b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2649b54502bSHong Zhang PetscErrorCode ierr; 2659b54502bSHong Zhang 2669b54502bSHong Zhang PetscFunctionBegin; 2679b54502bSHong Zhang if (!dir->inplace && dir->fact) {ierr = MatDestroy(dir->fact);CHKERRQ(ierr);} 2689b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2699b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 2709b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 2719b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 2729b54502bSHong Zhang PetscFunctionReturn(0); 2739b54502bSHong Zhang } 2749b54502bSHong Zhang 2759b54502bSHong Zhang #undef __FUNCT__ 2769b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 2779b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 2789b54502bSHong Zhang { 2799b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2809b54502bSHong Zhang PetscErrorCode ierr; 2819b54502bSHong Zhang 2829b54502bSHong Zhang PetscFunctionBegin; 2839b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 2849b54502bSHong Zhang else {ierr = MatSolve(dir->fact,x,y);CHKERRQ(ierr);} 2859b54502bSHong Zhang PetscFunctionReturn(0); 2869b54502bSHong Zhang } 2879b54502bSHong Zhang 2889b54502bSHong Zhang #undef __FUNCT__ 2899b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2909b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2919b54502bSHong Zhang { 2929b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2939b54502bSHong Zhang PetscErrorCode ierr; 2949b54502bSHong Zhang 2959b54502bSHong Zhang PetscFunctionBegin; 2969b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 2979b54502bSHong Zhang else {ierr = MatSolveTranspose(dir->fact,x,y);CHKERRQ(ierr);} 2989b54502bSHong Zhang PetscFunctionReturn(0); 2999b54502bSHong Zhang } 3009b54502bSHong Zhang 3019b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3029b54502bSHong Zhang 3039b54502bSHong Zhang EXTERN_C_BEGIN 3049b54502bSHong Zhang #undef __FUNCT__ 3059b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill_LU" 3069b54502bSHong Zhang PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill) 3079b54502bSHong Zhang { 3089b54502bSHong Zhang PC_LU *dir; 3099b54502bSHong Zhang 3109b54502bSHong Zhang PetscFunctionBegin; 3119b54502bSHong Zhang dir = (PC_LU*)pc->data; 3129b54502bSHong Zhang dir->info.fill = fill; 3139b54502bSHong Zhang PetscFunctionReturn(0); 3149b54502bSHong Zhang } 3159b54502bSHong Zhang EXTERN_C_END 3169b54502bSHong Zhang 3179b54502bSHong Zhang EXTERN_C_BEGIN 3189b54502bSHong Zhang #undef __FUNCT__ 3199b54502bSHong Zhang #define __FUNCT__ "PCLUSetDamping_LU" 3209b54502bSHong Zhang PetscErrorCode PCLUSetDamping_LU(PC pc,PetscReal damping) 3219b54502bSHong Zhang { 3229b54502bSHong Zhang PC_LU *dir; 3239b54502bSHong Zhang 3249b54502bSHong Zhang PetscFunctionBegin; 3259b54502bSHong Zhang dir = (PC_LU*)pc->data; 3269b54502bSHong Zhang if (damping == (PetscReal) PETSC_DECIDE) { 327*0a29876aSHong Zhang dir->info.shiftnz = 1.e-12; 3289b54502bSHong Zhang } else { 329*0a29876aSHong Zhang dir->info.shiftnz = damping; 3309b54502bSHong Zhang } 3319b54502bSHong Zhang PetscFunctionReturn(0); 3329b54502bSHong Zhang } 3339b54502bSHong Zhang EXTERN_C_END 3349b54502bSHong Zhang 3359b54502bSHong Zhang EXTERN_C_BEGIN 3369b54502bSHong Zhang #undef __FUNCT__ 3379b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace_LU" 3389b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace_LU(PC pc) 3399b54502bSHong Zhang { 3409b54502bSHong Zhang PC_LU *dir; 3419b54502bSHong Zhang 3429b54502bSHong Zhang PetscFunctionBegin; 3439b54502bSHong Zhang dir = (PC_LU*)pc->data; 3449b54502bSHong Zhang dir->inplace = PETSC_TRUE; 3459b54502bSHong Zhang PetscFunctionReturn(0); 3469b54502bSHong Zhang } 3479b54502bSHong Zhang EXTERN_C_END 3489b54502bSHong Zhang 3499b54502bSHong Zhang EXTERN_C_BEGIN 3509b54502bSHong Zhang #undef __FUNCT__ 3519b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering_LU" 3529b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering) 3539b54502bSHong Zhang { 3549b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3559b54502bSHong Zhang PetscErrorCode ierr; 3569b54502bSHong Zhang 3579b54502bSHong Zhang PetscFunctionBegin; 3589b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 3599b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 3609b54502bSHong Zhang PetscFunctionReturn(0); 3619b54502bSHong Zhang } 3629b54502bSHong Zhang EXTERN_C_END 3639b54502bSHong Zhang 3649b54502bSHong Zhang EXTERN_C_BEGIN 3659b54502bSHong Zhang #undef __FUNCT__ 3669b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting_LU" 3679b54502bSHong Zhang PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol) 3689b54502bSHong Zhang { 3699b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3709b54502bSHong Zhang 3719b54502bSHong Zhang PetscFunctionBegin; 3729b54502bSHong 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); 3739b54502bSHong Zhang dir->info.dtcol = dtcol; 3749b54502bSHong Zhang PetscFunctionReturn(0); 3759b54502bSHong Zhang } 3769b54502bSHong Zhang EXTERN_C_END 3779b54502bSHong Zhang 3789b54502bSHong Zhang EXTERN_C_BEGIN 3799b54502bSHong Zhang #undef __FUNCT__ 3809b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks_LU" 3819b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot) 3829b54502bSHong Zhang { 3839b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 3849b54502bSHong Zhang 3859b54502bSHong Zhang PetscFunctionBegin; 3869b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 3879b54502bSHong Zhang PetscFunctionReturn(0); 3889b54502bSHong Zhang } 3899b54502bSHong Zhang EXTERN_C_END 3909b54502bSHong Zhang 3919b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 3929b54502bSHong Zhang 3939b54502bSHong Zhang #undef __FUNCT__ 3949b54502bSHong Zhang #define __FUNCT__ "PCLUReorderForNonzeroDiagonal" 3959b54502bSHong Zhang /*@ 3969b54502bSHong Zhang PCLUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 3979b54502bSHong Zhang 3989b54502bSHong Zhang Collective on PC 3999b54502bSHong Zhang 4009b54502bSHong Zhang Input Parameters: 4019b54502bSHong Zhang + pc - the preconditioner context 4029b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 4039b54502bSHong Zhang 4049b54502bSHong Zhang Options Database Key: 4059b54502bSHong Zhang . -pc_lu_nonzeros_along_diagonal 4069b54502bSHong Zhang 4079b54502bSHong Zhang Level: intermediate 4089b54502bSHong Zhang 4099b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4109b54502bSHong Zhang 4119b54502bSHong Zhang .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot(), MatReorderForNonzeroDiagonal() 4129b54502bSHong Zhang @*/ 4139b54502bSHong Zhang PetscErrorCode PCLUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 4149b54502bSHong Zhang { 4159b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4169b54502bSHong Zhang 4179b54502bSHong Zhang PetscFunctionBegin; 4189b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4199b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 4209b54502bSHong Zhang if (f) { 4219b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 4229b54502bSHong Zhang } 4239b54502bSHong Zhang PetscFunctionReturn(0); 4249b54502bSHong Zhang } 4259b54502bSHong Zhang 4269b54502bSHong Zhang #undef __FUNCT__ 4279b54502bSHong Zhang #define __FUNCT__ "PCLUSetZeroPivot" 4289b54502bSHong Zhang /*@ 4299b54502bSHong Zhang PCLUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 4309b54502bSHong Zhang 4319b54502bSHong Zhang Collective on PC 4329b54502bSHong Zhang 4339b54502bSHong Zhang Input Parameters: 4349b54502bSHong Zhang + pc - the preconditioner context 4359b54502bSHong Zhang - zero - all pivots smaller than this will be considered zero 4369b54502bSHong Zhang 4379b54502bSHong Zhang Options Database Key: 4389b54502bSHong Zhang . -pc_lu_zeropivot <zero> - Sets the zero pivot size 4399b54502bSHong Zhang 4409b54502bSHong Zhang Level: intermediate 4419b54502bSHong Zhang 4429b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4439b54502bSHong Zhang 4449b54502bSHong Zhang .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot() 4459b54502bSHong Zhang @*/ 4469b54502bSHong Zhang PetscErrorCode PCLUSetZeroPivot(PC pc,PetscReal zero) 4479b54502bSHong Zhang { 4489b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4499b54502bSHong Zhang 4509b54502bSHong Zhang PetscFunctionBegin; 4519b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4529b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 4539b54502bSHong Zhang if (f) { 4549b54502bSHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 4559b54502bSHong Zhang } 4569b54502bSHong Zhang PetscFunctionReturn(0); 4579b54502bSHong Zhang } 4589b54502bSHong Zhang 4599b54502bSHong Zhang #undef __FUNCT__ 4609b54502bSHong Zhang #define __FUNCT__ "PCLUSetShift" 4619b54502bSHong Zhang /*@ 4629b54502bSHong Zhang PCLUSetShift - specify whether to use Manteuffel shifting of LU. 4639b54502bSHong Zhang If an LU factorisation breaks down because of nonpositive pivots, 4649b54502bSHong Zhang adding sufficient identity to the diagonal will remedy this. 4659b54502bSHong Zhang Setting this causes a bisection method to find the minimum shift that 4669b54502bSHong Zhang will lead to a well-defined LU. 4679b54502bSHong Zhang 4689b54502bSHong Zhang Input parameters: 4699b54502bSHong Zhang + pc - the preconditioner context 4709b54502bSHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 4719b54502bSHong Zhang 4729b54502bSHong Zhang Options Database Key: 4739b54502bSHong Zhang . -pc_lu_shift - Activate PCLUSetShift() 4749b54502bSHong Zhang 4759b54502bSHong Zhang Level: intermediate 4769b54502bSHong Zhang 4779b54502bSHong Zhang .keywords: PC, indefinite, factorization 4789b54502bSHong Zhang 4799b54502bSHong Zhang .seealso: PCLUSetDamping(), PCILUSetShift() 4809b54502bSHong Zhang @*/ 4819b54502bSHong Zhang PetscErrorCode PCLUSetShift(PC pc,PetscTruth shifting) 4829b54502bSHong Zhang { 4839b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4849b54502bSHong Zhang 4859b54502bSHong Zhang PetscFunctionBegin; 4869b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4879b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 4889b54502bSHong Zhang if (f) { 4899b54502bSHong Zhang ierr = (*f)(pc,shifting);CHKERRQ(ierr); 4909b54502bSHong Zhang } 4919b54502bSHong Zhang PetscFunctionReturn(0); 4929b54502bSHong Zhang } 4939b54502bSHong Zhang 4949b54502bSHong Zhang 4959b54502bSHong Zhang #undef __FUNCT__ 4969b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseOrdering" 4979b54502bSHong Zhang /*@ 4989b54502bSHong Zhang PCLUSetReuseOrdering - When similar matrices are factored, this 4999b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 5009b54502bSHong Zhang following factors; applies to both fill and drop tolerance LUs. 5019b54502bSHong Zhang 5029b54502bSHong Zhang Collective on PC 5039b54502bSHong Zhang 5049b54502bSHong Zhang Input Parameters: 5059b54502bSHong Zhang + pc - the preconditioner context 5069b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 5079b54502bSHong Zhang 5089b54502bSHong Zhang Options Database Key: 5099b54502bSHong Zhang . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 5109b54502bSHong Zhang 5119b54502bSHong Zhang Level: intermediate 5129b54502bSHong Zhang 5139b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 5149b54502bSHong Zhang 5159b54502bSHong Zhang .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill() 5169b54502bSHong Zhang @*/ 5179b54502bSHong Zhang PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag) 5189b54502bSHong Zhang { 5199b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 5209b54502bSHong Zhang 5219b54502bSHong Zhang PetscFunctionBegin; 5229b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5239b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 5249b54502bSHong Zhang if (f) { 5259b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 5269b54502bSHong Zhang } 5279b54502bSHong Zhang PetscFunctionReturn(0); 5289b54502bSHong Zhang } 5299b54502bSHong Zhang 5309b54502bSHong Zhang #undef __FUNCT__ 5319b54502bSHong Zhang #define __FUNCT__ "PCLUSetReuseFill" 5329b54502bSHong Zhang /*@ 5339b54502bSHong Zhang PCLUSetReuseFill - When matrices with same nonzero structure are LU factored, 5349b54502bSHong Zhang this causes later ones to use the fill computed in the initial factorization. 5359b54502bSHong Zhang 5369b54502bSHong Zhang Collective on PC 5379b54502bSHong Zhang 5389b54502bSHong Zhang Input Parameters: 5399b54502bSHong Zhang + pc - the preconditioner context 5409b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 5419b54502bSHong Zhang 5429b54502bSHong Zhang Options Database Key: 5439b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 5449b54502bSHong Zhang 5459b54502bSHong Zhang Level: intermediate 5469b54502bSHong Zhang 5479b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, LU 5489b54502bSHong Zhang 5499b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill() 5509b54502bSHong Zhang @*/ 5519b54502bSHong Zhang PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag) 5529b54502bSHong Zhang { 5539b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 5549b54502bSHong Zhang 5559b54502bSHong Zhang PetscFunctionBegin; 5569b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5579b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5589b54502bSHong Zhang if (f) { 5599b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 5609b54502bSHong Zhang } 5619b54502bSHong Zhang PetscFunctionReturn(0); 5629b54502bSHong Zhang } 5639b54502bSHong Zhang 5649b54502bSHong Zhang #undef __FUNCT__ 5659b54502bSHong Zhang #define __FUNCT__ "PCLUSetFill" 5669b54502bSHong Zhang /*@ 5679b54502bSHong Zhang PCLUSetFill - Indicate the amount of fill you expect in the factored matrix, 5689b54502bSHong Zhang fill = number nonzeros in factor/number nonzeros in original matrix. 5699b54502bSHong Zhang 5709b54502bSHong Zhang Collective on PC 5719b54502bSHong Zhang 5729b54502bSHong Zhang Input Parameters: 5739b54502bSHong Zhang + pc - the preconditioner context 5749b54502bSHong Zhang - fill - amount of expected fill 5759b54502bSHong Zhang 5769b54502bSHong Zhang Options Database Key: 5779b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 5789b54502bSHong Zhang 5799b54502bSHong Zhang Level: intermediate 5809b54502bSHong Zhang 5819b54502bSHong Zhang Note: 5829b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 5839b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 5849b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 5859b54502bSHong Zhang future runs. Default PETSc uses a value of 5.0 5869b54502bSHong Zhang 5879b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 5889b54502bSHong Zhang 5899b54502bSHong Zhang .seealso: PCILUSetFill() 5909b54502bSHong Zhang @*/ 5919b54502bSHong Zhang PetscErrorCode PCLUSetFill(PC pc,PetscReal fill) 5929b54502bSHong Zhang { 5939b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5949b54502bSHong Zhang 5959b54502bSHong Zhang PetscFunctionBegin; 5969b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5979b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 5989b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5999b54502bSHong Zhang if (f) { 6009b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 6019b54502bSHong Zhang } 6029b54502bSHong Zhang PetscFunctionReturn(0); 6039b54502bSHong Zhang } 6049b54502bSHong Zhang 6059b54502bSHong Zhang #undef __FUNCT__ 6069b54502bSHong Zhang #define __FUNCT__ "PCLUSetDamping" 6079b54502bSHong Zhang /*@ 6089b54502bSHong Zhang PCLUSetDamping - adds this quantity to the diagonal of the matrix during the 6099b54502bSHong Zhang LU numerical factorization 6109b54502bSHong Zhang 6119b54502bSHong Zhang Collective on PC 6129b54502bSHong Zhang 6139b54502bSHong Zhang Input Parameters: 6149b54502bSHong Zhang + pc - the preconditioner context 6159b54502bSHong Zhang - damping - amount of damping (use PETSC_DECIDE for default of 1.e-12) 6169b54502bSHong Zhang 6179b54502bSHong Zhang Options Database Key: 6189b54502bSHong Zhang . -pc_lu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 6199b54502bSHong Zhang 6209b54502bSHong Zhang Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 6219b54502bSHong Zhang pivot, then the damping is doubled until this is alleviated. 6229b54502bSHong Zhang 6239b54502bSHong Zhang Level: intermediate 6249b54502bSHong Zhang 6259b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 6269b54502bSHong Zhang .seealso: PCILUSetFill(), PCILUSetDamp() 6279b54502bSHong Zhang @*/ 6289b54502bSHong Zhang PetscErrorCode PCLUSetDamping(PC pc,PetscReal damping) 6299b54502bSHong Zhang { 6309b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 6319b54502bSHong Zhang 6329b54502bSHong Zhang PetscFunctionBegin; 6339b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6349b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 6359b54502bSHong Zhang if (f) { 6369b54502bSHong Zhang ierr = (*f)(pc,damping);CHKERRQ(ierr); 6379b54502bSHong Zhang } 6389b54502bSHong Zhang PetscFunctionReturn(0); 6399b54502bSHong Zhang } 6409b54502bSHong Zhang 6419b54502bSHong Zhang #undef __FUNCT__ 6429b54502bSHong Zhang #define __FUNCT__ "PCLUSetUseInPlace" 6439b54502bSHong Zhang /*@ 6449b54502bSHong Zhang PCLUSetUseInPlace - Tells the system to do an in-place factorization. 6459b54502bSHong Zhang For dense matrices, this enables the solution of much larger problems. 6469b54502bSHong Zhang For sparse matrices the factorization cannot be done truly in-place 6479b54502bSHong Zhang so this does not save memory during the factorization, but after the matrix 6489b54502bSHong Zhang is factored, the original unfactored matrix is freed, thus recovering that 6499b54502bSHong Zhang space. 6509b54502bSHong Zhang 6519b54502bSHong Zhang Collective on PC 6529b54502bSHong Zhang 6539b54502bSHong Zhang Input Parameters: 6549b54502bSHong Zhang . pc - the preconditioner context 6559b54502bSHong Zhang 6569b54502bSHong Zhang Options Database Key: 6579b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 6589b54502bSHong Zhang 6599b54502bSHong Zhang Notes: 6609b54502bSHong Zhang PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 6619b54502bSHong Zhang a different matrix is provided for the multiply and the preconditioner in 6629b54502bSHong Zhang a call to KSPSetOperators(). 6639b54502bSHong Zhang This is because the Krylov space methods require an application of the 6649b54502bSHong Zhang matrix multiplication, which is not possible here because the matrix has 6659b54502bSHong Zhang been factored in-place, replacing the original matrix. 6669b54502bSHong Zhang 6679b54502bSHong Zhang Level: intermediate 6689b54502bSHong Zhang 6699b54502bSHong Zhang .keywords: PC, set, factorization, direct, inplace, in-place, LU 6709b54502bSHong Zhang 6719b54502bSHong Zhang .seealso: PCILUSetUseInPlace() 6729b54502bSHong Zhang @*/ 6739b54502bSHong Zhang PetscErrorCode PCLUSetUseInPlace(PC pc) 6749b54502bSHong Zhang { 6759b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 6769b54502bSHong Zhang 6779b54502bSHong Zhang PetscFunctionBegin; 6789b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6799b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 6809b54502bSHong Zhang if (f) { 6819b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 6829b54502bSHong Zhang } 6839b54502bSHong Zhang PetscFunctionReturn(0); 6849b54502bSHong Zhang } 6859b54502bSHong Zhang 6869b54502bSHong Zhang #undef __FUNCT__ 6879b54502bSHong Zhang #define __FUNCT__ "PCLUSetMatOrdering" 6889b54502bSHong Zhang /*@C 6899b54502bSHong Zhang PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 6909b54502bSHong Zhang be used in the LU factorization. 6919b54502bSHong Zhang 6929b54502bSHong Zhang Collective on PC 6939b54502bSHong Zhang 6949b54502bSHong Zhang Input Parameters: 6959b54502bSHong Zhang + pc - the preconditioner context 6969b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 6979b54502bSHong Zhang 6989b54502bSHong Zhang Options Database Key: 6999b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 7009b54502bSHong Zhang 7019b54502bSHong Zhang Level: intermediate 7029b54502bSHong Zhang 7039b54502bSHong Zhang Notes: nested dissection is used by default 7049b54502bSHong Zhang 7059b54502bSHong Zhang .seealso: PCILUSetMatOrdering() 7069b54502bSHong Zhang @*/ 7079b54502bSHong Zhang PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering) 7089b54502bSHong Zhang { 7099b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 7109b54502bSHong Zhang 7119b54502bSHong Zhang PetscFunctionBegin; 7129b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 7139b54502bSHong Zhang if (f) { 7149b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 7159b54502bSHong Zhang } 7169b54502bSHong Zhang PetscFunctionReturn(0); 7179b54502bSHong Zhang } 7189b54502bSHong Zhang 7199b54502bSHong Zhang #undef __FUNCT__ 7209b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivoting" 7219b54502bSHong Zhang /*@ 7229b54502bSHong Zhang PCLUSetPivoting - Determines when pivoting is done during LU. 7239b54502bSHong Zhang For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 7249b54502bSHong Zhang it is never done. For the Matlab and SuperLU factorization this is used. 7259b54502bSHong Zhang 7269b54502bSHong Zhang Collective on PC 7279b54502bSHong Zhang 7289b54502bSHong Zhang Input Parameters: 7299b54502bSHong Zhang + pc - the preconditioner context 7309b54502bSHong Zhang - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 7319b54502bSHong Zhang 7329b54502bSHong Zhang Options Database Key: 7339b54502bSHong Zhang . -pc_lu_pivoting <dtcol> 7349b54502bSHong Zhang 7359b54502bSHong Zhang Level: intermediate 7369b54502bSHong Zhang 7379b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks() 7389b54502bSHong Zhang @*/ 7399b54502bSHong Zhang PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol) 7409b54502bSHong Zhang { 7419b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 7429b54502bSHong Zhang 7439b54502bSHong Zhang PetscFunctionBegin; 7449b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);CHKERRQ(ierr); 7459b54502bSHong Zhang if (f) { 7469b54502bSHong Zhang ierr = (*f)(pc,dtcol);CHKERRQ(ierr); 7479b54502bSHong Zhang } 7489b54502bSHong Zhang PetscFunctionReturn(0); 7499b54502bSHong Zhang } 7509b54502bSHong Zhang 7519b54502bSHong Zhang #undef __FUNCT__ 7529b54502bSHong Zhang #define __FUNCT__ "PCLUSetPivotInBlocks" 7539b54502bSHong Zhang /*@ 7549b54502bSHong Zhang PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block 7559b54502bSHong Zhang with BAIJ or SBAIJ matrices 7569b54502bSHong Zhang 7579b54502bSHong Zhang Collective on PC 7589b54502bSHong Zhang 7599b54502bSHong Zhang Input Parameters: 7609b54502bSHong Zhang + pc - the preconditioner context 7619b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 7629b54502bSHong Zhang 7639b54502bSHong Zhang Options Database Key: 7649b54502bSHong Zhang . -pc_lu_pivot_in_blocks <true,false> 7659b54502bSHong Zhang 7669b54502bSHong Zhang Level: intermediate 7679b54502bSHong Zhang 7689b54502bSHong Zhang .seealso: PCILUSetMatOrdering(), PCLUSetPivoting() 7699b54502bSHong Zhang @*/ 7709b54502bSHong Zhang PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot) 7719b54502bSHong Zhang { 7729b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 7739b54502bSHong Zhang 7749b54502bSHong Zhang PetscFunctionBegin; 7759b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 7769b54502bSHong Zhang if (f) { 7779b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 7789b54502bSHong Zhang } 7799b54502bSHong Zhang PetscFunctionReturn(0); 7809b54502bSHong Zhang } 7819b54502bSHong Zhang 7829b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 7839b54502bSHong Zhang 7849b54502bSHong Zhang /*MC 7859b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 7869b54502bSHong Zhang 7879b54502bSHong Zhang Options Database Keys: 7889b54502bSHong Zhang + -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering() 7899b54502bSHong Zhang . -pc_lu_reuse_fill - Activates PCLUSetReuseFill() 7909b54502bSHong Zhang . -pc_lu_fill <fill> - Sets fill amount 7919b54502bSHong Zhang . -pc_lu_damping <damping> - Sets damping amount 7929b54502bSHong Zhang . -pc_lu_in_place - Activates in-place factorization 7939b54502bSHong Zhang . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 7949b54502bSHong Zhang . -pc_lu_pivoting <dtcol> 7959b54502bSHong Zhang - -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 7969b54502bSHong Zhang stability of factorization. 7979b54502bSHong Zhang 7989b54502bSHong Zhang Notes: Not all options work for all matrix formats 7999b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 8009b54502bSHong Zhang algorithms 8019b54502bSHong Zhang 8029b54502bSHong Zhang Level: beginner 8039b54502bSHong Zhang 8049b54502bSHong Zhang Concepts: LU factorization, direct solver 8059b54502bSHong Zhang 8069b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 8079b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 8089b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 8099b54502bSHong Zhang 8109b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 8119b54502bSHong Zhang PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(), 8129b54502bSHong Zhang PCLUSetFill(), PCLUSetDamping(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCLUSetPivoting(), 8139b54502bSHong Zhang PCLUSetPivotingInBlocks() 8149b54502bSHong Zhang M*/ 8159b54502bSHong Zhang 8169b54502bSHong Zhang EXTERN_C_BEGIN 8179b54502bSHong Zhang #undef __FUNCT__ 8189b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 8199b54502bSHong Zhang PetscErrorCode PCCreate_LU(PC pc) 8209b54502bSHong Zhang { 8219b54502bSHong Zhang PetscErrorCode ierr; 8229b54502bSHong Zhang PetscMPIInt size; 8239b54502bSHong Zhang PC_LU *dir; 8249b54502bSHong Zhang 8259b54502bSHong Zhang PetscFunctionBegin; 8269b54502bSHong Zhang ierr = PetscNew(PC_LU,&dir);CHKERRQ(ierr); 8279b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_LU)); 8289b54502bSHong Zhang 8299b54502bSHong Zhang ierr = MatFactorInfoInitialize(&dir->info);CHKERRQ(ierr); 8309b54502bSHong Zhang dir->fact = 0; 8319b54502bSHong Zhang dir->inplace = PETSC_FALSE; 8329b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 8339b54502bSHong Zhang 8349b54502bSHong Zhang dir->info.fill = 5.0; 8359b54502bSHong Zhang dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 836*0a29876aSHong Zhang dir->info.shiftnz = 0.0; 8379b54502bSHong Zhang dir->info.zeropivot = 1.e-12; 8389b54502bSHong Zhang dir->info.pivotinblocks = 1.0; 839*0a29876aSHong Zhang dir->info.shiftpd = PETSC_FALSE; 8409b54502bSHong Zhang dir->info.shift_fraction = 0.0; 8419b54502bSHong Zhang dir->col = 0; 8429b54502bSHong Zhang dir->row = 0; 8439b54502bSHong Zhang ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 8449b54502bSHong Zhang if (size == 1) { 8459b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_ND,&dir->ordering);CHKERRQ(ierr); 8469b54502bSHong Zhang } else { 8479b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);CHKERRQ(ierr); 8489b54502bSHong Zhang } 8499b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 8509b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 8519b54502bSHong Zhang pc->data = (void*)dir; 8529b54502bSHong Zhang 8539b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 8549b54502bSHong Zhang pc->ops->apply = PCApply_LU; 8559b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 8569b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 8579b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 8589b54502bSHong Zhang pc->ops->view = PCView_LU; 8599b54502bSHong Zhang pc->ops->applyrichardson = 0; 8609b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU; 8619b54502bSHong Zhang 8629b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU", 8639b54502bSHong Zhang PCLUSetFill_LU);CHKERRQ(ierr); 8649b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU", 8659b54502bSHong Zhang PCLUSetDamping_LU);CHKERRQ(ierr); 8669b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetShift_C","PCLUSetShift_LU", 8679b54502bSHong Zhang PCLUSetShift_LU);CHKERRQ(ierr); 8689b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU", 8699b54502bSHong Zhang PCLUSetUseInPlace_LU);CHKERRQ(ierr); 8709b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU", 8719b54502bSHong Zhang PCLUSetMatOrdering_LU);CHKERRQ(ierr); 8729b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU", 8739b54502bSHong Zhang PCLUSetReuseOrdering_LU);CHKERRQ(ierr); 8749b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU", 8759b54502bSHong Zhang PCLUSetReuseFill_LU);CHKERRQ(ierr); 8769b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU", 8779b54502bSHong Zhang PCLUSetPivoting_LU);CHKERRQ(ierr); 8789b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU", 8799b54502bSHong Zhang PCLUSetPivotInBlocks_LU);CHKERRQ(ierr); 8809b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetZeroPivot_C","PCLUSetZeroPivot_LU", 8819b54502bSHong Zhang PCLUSetZeroPivot_LU);CHKERRQ(ierr); 8829b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUReorderForNonzeroDiagonal_C","PCLUReorderForNonzeroDiagonal_LU", 8839b54502bSHong Zhang PCLUReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 8849b54502bSHong Zhang PetscFunctionReturn(0); 8859b54502bSHong Zhang } 8869b54502bSHong Zhang EXTERN_C_END 887