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" 147087cfbeSBarry Smith PetscErrorCode 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" 327087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool 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" 457087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool 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; 61ace3abfcSBarry Smith PetscBool 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; 84ace3abfcSBarry Smith PetscBool iascii; 859b54502bSHong Zhang 869b54502bSHong Zhang PetscFunctionBegin; 872692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 889b54502bSHong Zhang if (iascii) { 89914a5d51SHong Zhang if (lu->inplace) { 90914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 91914a5d51SHong Zhang } else { 92914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 93145b9266SHong Zhang } 94145b9266SHong Zhang 959b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 969b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 979b54502bSHong Zhang } 98914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 999b54502bSHong Zhang PetscFunctionReturn(0); 1009b54502bSHong Zhang } 1019b54502bSHong Zhang 1029b54502bSHong Zhang #undef __FUNCT__ 1039b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 1049b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 1059b54502bSHong Zhang { 1069b54502bSHong Zhang PetscErrorCode ierr; 1079b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1089b54502bSHong Zhang 1099b54502bSHong Zhang PetscFunctionBegin; 110075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1119b54502bSHong Zhang 1129b54502bSHong Zhang if (dir->inplace) { 1139b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1149b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 115075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 11603c60df9SBarry Smith if (dir->row) { 11703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 11803c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 11903c60df9SBarry Smith } 120075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 121075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1229b54502bSHong Zhang } else { 1239b54502bSHong Zhang MatInfo info; 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) { 13003c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 13103c60df9SBarry Smith ierr = PetscLogObjectParent(pc,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; 139075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1409b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1419b54502bSHong Zhang if (!dir->reuseordering) { 1429b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1439b54502bSHong Zhang if (dir->col) {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) { 14903c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 15003c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 15103c60df9SBarry Smith } 1529b54502bSHong Zhang } 153075768bcSBarry 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; 158075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1599b54502bSHong Zhang } 160075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1619b54502bSHong Zhang } 1629b54502bSHong Zhang PetscFunctionReturn(0); 1639b54502bSHong Zhang } 1649b54502bSHong Zhang 1659b54502bSHong Zhang #undef __FUNCT__ 166*574deadeSBarry Smith #define __FUNCT__ "PCReset_LU" 167*574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1689b54502bSHong Zhang { 1699b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1709b54502bSHong Zhang PetscErrorCode ierr; 1719b54502bSHong Zhang 1729b54502bSHong Zhang PetscFunctionBegin; 173075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 1749b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1759b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 176*574deadeSBarry Smith PetscFunctionReturn(0); 177*574deadeSBarry Smith } 178*574deadeSBarry Smith 179*574deadeSBarry Smith #undef __FUNCT__ 180*574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU" 181*574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 182*574deadeSBarry Smith { 183*574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 184*574deadeSBarry Smith PetscErrorCode ierr; 185*574deadeSBarry Smith 186*574deadeSBarry Smith PetscFunctionBegin; 187*574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 188503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 189503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 190c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1919b54502bSHong Zhang PetscFunctionReturn(0); 1929b54502bSHong Zhang } 1939b54502bSHong Zhang 1949b54502bSHong Zhang #undef __FUNCT__ 1959b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 1969b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1979b54502bSHong Zhang { 1989b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1999b54502bSHong Zhang PetscErrorCode ierr; 2009b54502bSHong Zhang 2019b54502bSHong Zhang PetscFunctionBegin; 2029b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 203075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang #undef __FUNCT__ 2089b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2099b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_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 = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 216075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2179b54502bSHong Zhang PetscFunctionReturn(0); 2189b54502bSHong Zhang } 2199b54502bSHong Zhang 2209b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2219b54502bSHong Zhang 2229b54502bSHong Zhang EXTERN_C_BEGIN 2239b54502bSHong Zhang #undef __FUNCT__ 2242401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2257087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 2269b54502bSHong Zhang { 227c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2289b54502bSHong Zhang 2299b54502bSHong Zhang PetscFunctionBegin; 2309b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2319b54502bSHong Zhang PetscFunctionReturn(0); 2329b54502bSHong Zhang } 2339b54502bSHong Zhang EXTERN_C_END 2349b54502bSHong Zhang 2359b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2369b54502bSHong Zhang 2379b54502bSHong Zhang /*MC 2389b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2399b54502bSHong Zhang 2409b54502bSHong Zhang Options Database Keys: 2412401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 242c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 2432401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 24455ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2452401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2462401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2472401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2489b54502bSHong Zhang stability of factorization. 249145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 250145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 251e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 252e22d95b2SBarry Smith diagonal. 2539b54502bSHong Zhang 2549b54502bSHong Zhang Notes: Not all options work for all matrix formats 2559b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2569b54502bSHong Zhang algorithms 2579b54502bSHong Zhang 2589b54502bSHong Zhang Level: beginner 2599b54502bSHong Zhang 2609b54502bSHong Zhang Concepts: LU factorization, direct solver 2619b54502bSHong Zhang 2629b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2639b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2649b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2659b54502bSHong Zhang 2669b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 267a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2688ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 269145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2708ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2719b54502bSHong Zhang M*/ 2729b54502bSHong Zhang 2739b54502bSHong Zhang EXTERN_C_BEGIN 2749b54502bSHong Zhang #undef __FUNCT__ 2759b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 2767087cfbeSBarry Smith PetscErrorCode PCCreate_LU(PC pc) 2779b54502bSHong Zhang { 2789b54502bSHong Zhang PetscErrorCode ierr; 2799b54502bSHong Zhang PetscMPIInt size; 2809b54502bSHong Zhang PC_LU *dir; 2819b54502bSHong Zhang 2829b54502bSHong Zhang PetscFunctionBegin; 28338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 2849b54502bSHong Zhang 285075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 286879e8a4dSBarry Smith ((PC_Factor*)dir)->fact = PETSC_NULL; 287d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2889b54502bSHong Zhang dir->inplace = PETSC_FALSE; 2899b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2909b54502bSHong Zhang 291075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 292075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 293f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 294915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 295075768bcSBarry Smith ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 296075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 2979b54502bSHong Zhang dir->col = 0; 2989b54502bSHong Zhang dir->row = 0; 2995c9eb25fSBarry Smith 300a1f19f5aSHong Zhang ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 3017adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 3029b54502bSHong Zhang if (size == 1) { 3032692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3049b54502bSHong Zhang } else { 3052692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3069b54502bSHong Zhang } 3079b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3089b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3099b54502bSHong Zhang pc->data = (void*)dir; 3109b54502bSHong Zhang 311*574deadeSBarry Smith pc->ops->reset = PCReset_LU; 3129b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3139b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3149b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3159b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3169b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3179b54502bSHong Zhang pc->ops->view = PCView_LU; 3189b54502bSHong Zhang pc->ops->applyrichardson = 0; 31985317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3209b54502bSHong Zhang 321f8260c8fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 322f8260c8fSBarry Smith PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 3237112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3247112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 32585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 32685317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 32785317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 32885317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 329d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 330d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 331d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 332d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 33385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 33485317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3352401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 3362401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 33785317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 33885317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3392401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 3402401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 3412401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 3422401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 3438ff23777SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 3448ff23777SHong Zhang PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 34585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 34685317021SBarry Smith PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 3472401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 3482401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3499b54502bSHong Zhang PetscFunctionReturn(0); 3509b54502bSHong Zhang } 3519b54502bSHong Zhang EXTERN_C_END 352