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 97c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/lu/lu.h" /*I "petscpc.h" I*/ 109b54502bSHong Zhang 119b54502bSHong Zhang EXTERN_C_BEGIN 129b54502bSHong Zhang #undef __FUNCT__ 132401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 159b54502bSHong Zhang { 169b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 179b54502bSHong Zhang 189b54502bSHong Zhang PetscFunctionBegin; 199b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 209b54502bSHong Zhang if (z == PETSC_DECIDE) { 219b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 229b54502bSHong Zhang } else { 239b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 249b54502bSHong Zhang } 259b54502bSHong Zhang PetscFunctionReturn(0); 269b54502bSHong Zhang } 279b54502bSHong Zhang EXTERN_C_END 289b54502bSHong Zhang 299b54502bSHong Zhang EXTERN_C_BEGIN 309b54502bSHong Zhang #undef __FUNCT__ 312401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 322401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag) 339b54502bSHong Zhang { 34c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 359b54502bSHong Zhang 369b54502bSHong Zhang PetscFunctionBegin; 379b54502bSHong Zhang lu->reuseordering = flag; 389b54502bSHong Zhang PetscFunctionReturn(0); 399b54502bSHong Zhang } 409b54502bSHong Zhang EXTERN_C_END 419b54502bSHong Zhang 429b54502bSHong Zhang EXTERN_C_BEGIN 439b54502bSHong Zhang #undef __FUNCT__ 442401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 452401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag) 469b54502bSHong Zhang { 47c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 489b54502bSHong Zhang 499b54502bSHong Zhang PetscFunctionBegin; 509b54502bSHong Zhang lu->reusefill = flag; 519b54502bSHong Zhang PetscFunctionReturn(0); 529b54502bSHong Zhang } 539b54502bSHong Zhang EXTERN_C_END 549b54502bSHong Zhang 559b54502bSHong Zhang #undef __FUNCT__ 569b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 579b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 589b54502bSHong Zhang { 599b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 609b54502bSHong Zhang PetscErrorCode ierr; 618ff23777SHong Zhang PetscTruth flg = PETSC_FALSE; 629b54502bSHong Zhang PetscReal tol; 639b54502bSHong Zhang 649b54502bSHong Zhang PetscFunctionBegin; 659b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 668ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 675c9eb25fSBarry Smith 682401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 699b54502bSHong Zhang if (flg) { 709b54502bSHong Zhang tol = PETSC_DECIDE; 712401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 722401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 739b54502bSHong Zhang } 749b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 759b54502bSHong Zhang PetscFunctionReturn(0); 769b54502bSHong Zhang } 779b54502bSHong Zhang 789b54502bSHong Zhang #undef __FUNCT__ 799b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 809b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 819b54502bSHong Zhang { 829b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 839b54502bSHong Zhang PetscErrorCode ierr; 849b54502bSHong Zhang PetscTruth iascii,isstring; 859b54502bSHong Zhang 869b54502bSHong Zhang PetscFunctionBegin; 879b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 889b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 899b54502bSHong Zhang if (iascii) { 909b54502bSHong Zhang 919b54502bSHong Zhang if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 929b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 93075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr); 94075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr); 95*145b9266SHong Zhang 96*145b9266SHong Zhang if (((PC_Factor*)lu)->info.shifttype==MAT_SHIFT_POSITIVE_DEFINITE) { 97*145b9266SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr); 98*145b9266SHong Zhang } 99*145b9266SHong Zhang if (((PC_Factor*)lu)->info.shifttype==MAT_SHIFT_NONZERO) { 100*145b9266SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr); 101*145b9266SHong Zhang } 102*145b9266SHong Zhang if (((PC_Factor*)lu)->info.shifttype==MAT_SHIFT_INBLOCKS) { 103*145b9266SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr); 104*145b9266SHong Zhang } 105*145b9266SHong Zhang 106075768bcSBarry Smith if (((PC_Factor*)lu)->fact) { 107f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 108f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 109f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 110f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 111f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 112f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 113075768bcSBarry Smith ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr); 1149b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 115f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 116f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 117f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1189b54502bSHong Zhang } 1199b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 1209b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1219b54502bSHong Zhang } else if (isstring) { 122075768bcSBarry Smith ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 1239b54502bSHong Zhang } else { 1249b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 1259b54502bSHong Zhang } 1269b54502bSHong Zhang PetscFunctionReturn(0); 1279b54502bSHong Zhang } 1289b54502bSHong Zhang 1299b54502bSHong Zhang #undef __FUNCT__ 1309b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 1319b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 1329b54502bSHong Zhang { 1339b54502bSHong Zhang PetscErrorCode ierr; 1349b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1359b54502bSHong Zhang 1369b54502bSHong Zhang PetscFunctionBegin; 137075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1389b54502bSHong Zhang 1399b54502bSHong Zhang if (dir->inplace) { 1409b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1419b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 142075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 14303c60df9SBarry Smith if (dir->row) { 14403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 14503c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 14603c60df9SBarry Smith } 147075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 148075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1499b54502bSHong Zhang } else { 1509b54502bSHong Zhang MatInfo info; 1519b54502bSHong Zhang if (!pc->setupcalled) { 152075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1539b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1549b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1559b54502bSHong Zhang } 15603c60df9SBarry Smith if (dir->row) { 15703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 15803c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 15903c60df9SBarry Smith } 160075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 161075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 162075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1639b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 164075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1659b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1669b54502bSHong Zhang if (!dir->reuseordering) { 1679b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1689b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 169075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1709b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1719b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1729b54502bSHong Zhang } 17303c60df9SBarry Smith if (dir->row) { 17403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 17503c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 17603c60df9SBarry Smith } 1779b54502bSHong Zhang } 178075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 179075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 180075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 181075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1829b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 183075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1849b54502bSHong Zhang } 185075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1869b54502bSHong Zhang } 1879b54502bSHong Zhang PetscFunctionReturn(0); 1889b54502bSHong Zhang } 1899b54502bSHong Zhang 1909b54502bSHong Zhang #undef __FUNCT__ 1919b54502bSHong Zhang #define __FUNCT__ "PCDestroy_LU" 1929b54502bSHong Zhang static PetscErrorCode PCDestroy_LU(PC pc) 1939b54502bSHong Zhang { 1949b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1959b54502bSHong Zhang PetscErrorCode ierr; 1969b54502bSHong Zhang 1979b54502bSHong Zhang PetscFunctionBegin; 198075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 1999b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 2009b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 201075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 202075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 2039b54502bSHong Zhang ierr = PetscFree(dir);CHKERRQ(ierr); 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang #undef __FUNCT__ 2089b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 2099b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 2109b54502bSHong Zhang { 2119b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2129b54502bSHong Zhang PetscErrorCode ierr; 2139b54502bSHong Zhang 2149b54502bSHong Zhang PetscFunctionBegin; 2159b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 216075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2179b54502bSHong Zhang PetscFunctionReturn(0); 2189b54502bSHong Zhang } 2199b54502bSHong Zhang 2209b54502bSHong Zhang #undef __FUNCT__ 2219b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2229b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2239b54502bSHong Zhang { 2249b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2259b54502bSHong Zhang PetscErrorCode ierr; 2269b54502bSHong Zhang 2279b54502bSHong Zhang PetscFunctionBegin; 2289b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 229075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2309b54502bSHong Zhang PetscFunctionReturn(0); 2319b54502bSHong Zhang } 2329b54502bSHong Zhang 2339b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2349b54502bSHong Zhang 2359b54502bSHong Zhang EXTERN_C_BEGIN 2369b54502bSHong Zhang #undef __FUNCT__ 2372401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2382401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 2399b54502bSHong Zhang { 240c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2419b54502bSHong Zhang 2429b54502bSHong Zhang PetscFunctionBegin; 2439b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2449b54502bSHong Zhang PetscFunctionReturn(0); 2459b54502bSHong Zhang } 2469b54502bSHong Zhang EXTERN_C_END 2479b54502bSHong Zhang 2489b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2499b54502bSHong Zhang 2509b54502bSHong Zhang /*MC 2519b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2529b54502bSHong Zhang 2539b54502bSHong Zhang Options Database Keys: 2542401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 255c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 2562401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 25755ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2582401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2592401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2602401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2619b54502bSHong Zhang stability of factorization. 262*145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 263*145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 264e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 265e22d95b2SBarry Smith diagonal. 2669b54502bSHong Zhang 2679b54502bSHong Zhang Notes: Not all options work for all matrix formats 2689b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2699b54502bSHong Zhang algorithms 2709b54502bSHong Zhang 2719b54502bSHong Zhang Level: beginner 2729b54502bSHong Zhang 2739b54502bSHong Zhang Concepts: LU factorization, direct solver 2749b54502bSHong Zhang 2759b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2769b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2779b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2789b54502bSHong Zhang 2799b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 280a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2818ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 282*145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2838ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2849b54502bSHong Zhang M*/ 2859b54502bSHong Zhang 2869b54502bSHong Zhang EXTERN_C_BEGIN 2879b54502bSHong Zhang #undef __FUNCT__ 2889b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 289dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 2909b54502bSHong Zhang { 2919b54502bSHong Zhang PetscErrorCode ierr; 2929b54502bSHong Zhang PetscMPIInt size; 2939b54502bSHong Zhang PC_LU *dir; 2949b54502bSHong Zhang 2959b54502bSHong Zhang PetscFunctionBegin; 29638f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 2979b54502bSHong Zhang 298075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 299075768bcSBarry Smith ((PC_Factor*)dir)->fact = 0; 3009b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3019b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 3029b54502bSHong Zhang 303075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 304075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 305915743fcSHong Zhang ((PC_Factor*)dir)->info.shifttype = MAT_SHIFT_NONE; 306915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 307075768bcSBarry Smith ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 308075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3099b54502bSHong Zhang dir->col = 0; 3109b54502bSHong Zhang dir->row = 0; 3115c9eb25fSBarry Smith 312075768bcSBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 3137adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 3149b54502bSHong Zhang if (size == 1) { 315075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3169b54502bSHong Zhang } else { 317075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3189b54502bSHong Zhang } 3199b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3209b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3219b54502bSHong Zhang pc->data = (void*)dir; 3229b54502bSHong Zhang 3239b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3249b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3259b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3269b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3279b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3289b54502bSHong Zhang pc->ops->view = PCView_LU; 3299b54502bSHong Zhang pc->ops->applyrichardson = 0; 33085317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3319b54502bSHong Zhang 3327112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3337112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 33485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 33585317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 33685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 33785317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 338d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 339d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 340d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 341d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 34285317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 34385317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3442401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 3452401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 34685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 34785317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3482401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 3492401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 3502401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 3512401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 3528ff23777SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 3538ff23777SHong Zhang PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 35485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 35585317021SBarry Smith PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 3562401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 3572401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3589b54502bSHong Zhang PetscFunctionReturn(0); 3599b54502bSHong Zhang } 3609b54502bSHong Zhang EXTERN_C_END 361