1dba47a55SKris Buschelman 29b54502bSHong Zhang /* 39b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation 49b54502bSHong Zhang Note: this need not be consided a preconditioner since it supplies 59b54502bSHong Zhang a direct solver. 69b54502bSHong Zhang */ 7ee45ca4aSHong Zhang 8c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/ 99b54502bSHong Zhang 10680c5173SHong Zhang 119b54502bSHong Zhang #undef __FUNCT__ 122401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 137087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 149b54502bSHong Zhang { 159b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 169b54502bSHong Zhang 179b54502bSHong Zhang PetscFunctionBegin; 189b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 192fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 202fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 219b54502bSHong Zhang PetscFunctionReturn(0); 229b54502bSHong Zhang } 239b54502bSHong Zhang 249b54502bSHong Zhang #undef __FUNCT__ 252401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 267087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag) 279b54502bSHong Zhang { 28c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 299b54502bSHong Zhang 309b54502bSHong Zhang PetscFunctionBegin; 319b54502bSHong Zhang lu->reuseordering = flag; 329b54502bSHong Zhang PetscFunctionReturn(0); 339b54502bSHong Zhang } 349b54502bSHong Zhang 359b54502bSHong Zhang #undef __FUNCT__ 362401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 377087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag) 389b54502bSHong Zhang { 39c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 409b54502bSHong Zhang 419b54502bSHong Zhang PetscFunctionBegin; 429b54502bSHong Zhang lu->reusefill = flag; 439b54502bSHong Zhang PetscFunctionReturn(0); 449b54502bSHong Zhang } 459b54502bSHong Zhang 469b54502bSHong Zhang #undef __FUNCT__ 479b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 484416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc) 499b54502bSHong Zhang { 509b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 519b54502bSHong Zhang PetscErrorCode ierr; 52ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 539b54502bSHong Zhang PetscReal tol; 549b54502bSHong Zhang 559b54502bSHong Zhang PetscFunctionBegin; 56e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr); 57e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 585c9eb25fSBarry Smith 592401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 609b54502bSHong Zhang if (flg) { 619b54502bSHong Zhang tol = PETSC_DECIDE; 622401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 632401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 649b54502bSHong Zhang } 659b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 669b54502bSHong Zhang PetscFunctionReturn(0); 679b54502bSHong Zhang } 689b54502bSHong Zhang 699b54502bSHong Zhang #undef __FUNCT__ 709b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 719b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 729b54502bSHong Zhang { 739b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 749b54502bSHong Zhang PetscErrorCode ierr; 75ace3abfcSBarry Smith PetscBool iascii; 769b54502bSHong Zhang 779b54502bSHong Zhang PetscFunctionBegin; 78251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 799b54502bSHong Zhang if (iascii) { 80914a5d51SHong Zhang if (lu->inplace) { 81914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 82914a5d51SHong Zhang } else { 83914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 84145b9266SHong Zhang } 85145b9266SHong Zhang 869b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 879b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 889b54502bSHong Zhang } 89914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 909b54502bSHong Zhang PetscFunctionReturn(0); 919b54502bSHong Zhang } 929b54502bSHong Zhang 939b54502bSHong Zhang #undef __FUNCT__ 949b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 959b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 969b54502bSHong Zhang { 979b54502bSHong Zhang PetscErrorCode ierr; 989b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 9900c67f3bSHong Zhang const MatSolverPackage stype; 10000e125f8SBarry Smith MatFactorError err; 1019b54502bSHong Zhang PetscFunctionBegin; 102075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1039b54502bSHong Zhang 10484d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 1059b54502bSHong Zhang if (dir->inplace) { 106fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 107fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 108075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 10903c60df9SBarry Smith if (dir->row) { 1103bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 11203c60df9SBarry Smith } 113075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 11400e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 11500e125f8SBarry Smith if (err) { /* Factor() fails */ 11600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1176baea169SHong Zhang PetscFunctionReturn(0); 1186baea169SHong Zhang } 1196baea169SHong Zhang 120075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1219b54502bSHong Zhang } else { 1229b54502bSHong Zhang MatInfo info; 12300e125f8SBarry Smith 1249b54502bSHong Zhang if (!pc->setupcalled) { 125075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1269b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1279b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1289b54502bSHong Zhang } 12903c60df9SBarry Smith if (dir->row) { 1303bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1313bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 13203c60df9SBarry Smith } 133d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 134075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 135a1f19f5aSHong Zhang } 136075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 137075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1389b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1393bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1409b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1419b54502bSHong Zhang if (!dir->reuseordering) { 142fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 143fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 144075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1459b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1469b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1479b54502bSHong Zhang } 14803c60df9SBarry Smith if (dir->row) { 1493bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1503bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 15103c60df9SBarry Smith } 1529b54502bSHong Zhang } 1536bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 154075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 155075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 156075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1579b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 1583bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 15904545d6dSBarry Smith } else { 160*b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 161*b8b68cfdSBarry Smith if (err == PC_FACTOR_NUMERIC_ZEROPIVOT) { 162*b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 163*b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 16404545d6dSBarry Smith } 1659b54502bSHong Zhang } 16600e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 16700e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 16800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1698c1cd74cSHong Zhang PetscFunctionReturn(0); 1708c1cd74cSHong Zhang } 1718c1cd74cSHong Zhang 172075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 17300e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 17400e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 17500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1768c1cd74cSHong Zhang } 177680c5173SHong Zhang 1789b54502bSHong Zhang } 17900c67f3bSHong Zhang 18000c67f3bSHong Zhang ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr); 18100c67f3bSHong Zhang if (!stype) { 18200e125f8SBarry Smith const MatSolverPackage solverpackage; 18300e125f8SBarry Smith ierr = MatFactorGetSolverPackage(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 18400e125f8SBarry Smith ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr); 18500c67f3bSHong Zhang } 1869b54502bSHong Zhang PetscFunctionReturn(0); 1879b54502bSHong Zhang } 1889b54502bSHong Zhang 1899b54502bSHong Zhang #undef __FUNCT__ 190574deadeSBarry Smith #define __FUNCT__ "PCReset_LU" 191574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1929b54502bSHong Zhang { 1939b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1949b54502bSHong Zhang PetscErrorCode ierr; 1959b54502bSHong Zhang 1969b54502bSHong Zhang PetscFunctionBegin; 1976bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 198fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 199fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 200574deadeSBarry Smith PetscFunctionReturn(0); 201574deadeSBarry Smith } 202574deadeSBarry Smith 203574deadeSBarry Smith #undef __FUNCT__ 204574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU" 205574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 206574deadeSBarry Smith { 207574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 208574deadeSBarry Smith PetscErrorCode ierr; 209574deadeSBarry Smith 210574deadeSBarry Smith PetscFunctionBegin; 211574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 212503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 213503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 214c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2159b54502bSHong Zhang PetscFunctionReturn(0); 2169b54502bSHong Zhang } 2179b54502bSHong Zhang 2189b54502bSHong Zhang #undef __FUNCT__ 2199b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 2209b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 2219b54502bSHong Zhang { 2229b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2239b54502bSHong Zhang PetscErrorCode ierr; 2249b54502bSHong Zhang 2259b54502bSHong Zhang PetscFunctionBegin; 2262fa5cd67SKarl Rupp if (dir->inplace) { 2272fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2282fa5cd67SKarl Rupp } else { 2292fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2302fa5cd67SKarl Rupp } 2319b54502bSHong Zhang PetscFunctionReturn(0); 2329b54502bSHong Zhang } 2339b54502bSHong Zhang 2349b54502bSHong Zhang #undef __FUNCT__ 2359b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2369b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2379b54502bSHong Zhang { 2389b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2399b54502bSHong Zhang PetscErrorCode ierr; 2409b54502bSHong Zhang 2419b54502bSHong Zhang PetscFunctionBegin; 2422fa5cd67SKarl Rupp if (dir->inplace) { 2432fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2442fa5cd67SKarl Rupp } else { 2452fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2462fa5cd67SKarl Rupp } 2479b54502bSHong Zhang PetscFunctionReturn(0); 2489b54502bSHong Zhang } 2499b54502bSHong Zhang 2509b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2519b54502bSHong Zhang 2529b54502bSHong Zhang #undef __FUNCT__ 2532401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2548e37d05fSBarry Smith PetscErrorCode PCFactorSetUseInPlace_LU(PC pc,PetscBool flg) 2559b54502bSHong Zhang { 256c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2579b54502bSHong Zhang 2589b54502bSHong Zhang PetscFunctionBegin; 2598e37d05fSBarry Smith dir->inplace = flg; 2608e37d05fSBarry Smith PetscFunctionReturn(0); 2618e37d05fSBarry Smith } 2628e37d05fSBarry Smith 2638e37d05fSBarry Smith #undef __FUNCT__ 2648e37d05fSBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_LU" 2658e37d05fSBarry Smith PetscErrorCode PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg) 2668e37d05fSBarry Smith { 2678e37d05fSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2688e37d05fSBarry Smith 2698e37d05fSBarry Smith PetscFunctionBegin; 2708e37d05fSBarry Smith *flg = dir->inplace; 2719b54502bSHong Zhang PetscFunctionReturn(0); 2729b54502bSHong Zhang } 2739b54502bSHong Zhang 2749b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2759b54502bSHong Zhang 2769b54502bSHong Zhang /*MC 2779b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2789b54502bSHong Zhang 2799b54502bSHong Zhang Options Database Keys: 2802401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 281d45987f3SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2822401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 28355ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2842401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2852401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2862401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2879b54502bSHong Zhang stability of factorization. 288145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 289145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 290e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 291e22d95b2SBarry Smith diagonal. 2929b54502bSHong Zhang 2939b54502bSHong Zhang Notes: Not all options work for all matrix formats 2949b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2959b54502bSHong Zhang algorithms 2969b54502bSHong Zhang 2979b54502bSHong Zhang Level: beginner 2989b54502bSHong Zhang 2999b54502bSHong Zhang Concepts: LU factorization, direct solver 3009b54502bSHong Zhang 3019b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 3029b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 3039b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 3049b54502bSHong Zhang 3059b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 306a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 3078ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 308145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 3098ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 3109b54502bSHong Zhang M*/ 3119b54502bSHong Zhang 3129b54502bSHong Zhang #undef __FUNCT__ 3139b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 3148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 3159b54502bSHong Zhang { 3169b54502bSHong Zhang PetscErrorCode ierr; 3179b54502bSHong Zhang PetscMPIInt size; 3189b54502bSHong Zhang PC_LU *dir; 3199b54502bSHong Zhang 3209b54502bSHong Zhang PetscFunctionBegin; 321b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 3229b54502bSHong Zhang 323075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 3242fa5cd67SKarl Rupp 3250298fd71SBarry Smith ((PC_Factor*)dir)->fact = NULL; 326d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 3279b54502bSHong Zhang dir->inplace = PETSC_FALSE; 3289b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 3299b54502bSHong Zhang 330075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 331075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 332f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 333915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3340ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 335075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3369b54502bSHong Zhang dir->col = 0; 3379b54502bSHong Zhang dir->row = 0; 3385c9eb25fSBarry Smith 339ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 3409b54502bSHong Zhang if (size == 1) { 34119fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3429b54502bSHong Zhang } else { 34319fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3449b54502bSHong Zhang } 3459b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3469b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3479b54502bSHong Zhang pc->data = (void*)dir; 3489b54502bSHong Zhang 349574deadeSBarry Smith pc->ops->reset = PCReset_LU; 3509b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3519b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3529b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3539b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3549b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3559b54502bSHong Zhang pc->ops->view = PCView_LU; 3569b54502bSHong Zhang pc->ops->applyrichardson = 0; 35785317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3589b54502bSHong Zhang 359bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 360bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 361bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 362bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 363bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 364bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 365bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 366bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 3678e37d05fSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);CHKERRQ(ierr); 368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 369bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 370bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr); 371bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 372bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 373bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3749b54502bSHong Zhang PetscFunctionReturn(0); 3759b54502bSHong Zhang } 376