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 8*c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/ 99b54502bSHong Zhang 109b54502bSHong Zhang EXTERN_C_BEGIN 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; 199b54502bSHong Zhang if (z == PETSC_DECIDE) { 209b54502bSHong Zhang lu->nonzerosalongdiagonaltol = 1.e-10; 219b54502bSHong Zhang } else { 229b54502bSHong Zhang lu->nonzerosalongdiagonaltol = z; 239b54502bSHong Zhang } 249b54502bSHong Zhang PetscFunctionReturn(0); 259b54502bSHong Zhang } 269b54502bSHong Zhang EXTERN_C_END 279b54502bSHong Zhang 289b54502bSHong Zhang EXTERN_C_BEGIN 299b54502bSHong Zhang #undef __FUNCT__ 302401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 317087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag) 329b54502bSHong Zhang { 33c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 349b54502bSHong Zhang 359b54502bSHong Zhang PetscFunctionBegin; 369b54502bSHong Zhang lu->reuseordering = flag; 379b54502bSHong Zhang PetscFunctionReturn(0); 389b54502bSHong Zhang } 399b54502bSHong Zhang EXTERN_C_END 409b54502bSHong Zhang 419b54502bSHong Zhang EXTERN_C_BEGIN 429b54502bSHong Zhang #undef __FUNCT__ 432401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 447087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag) 459b54502bSHong Zhang { 46c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 479b54502bSHong Zhang 489b54502bSHong Zhang PetscFunctionBegin; 499b54502bSHong Zhang lu->reusefill = flag; 509b54502bSHong Zhang PetscFunctionReturn(0); 519b54502bSHong Zhang } 529b54502bSHong Zhang EXTERN_C_END 539b54502bSHong Zhang 549b54502bSHong Zhang #undef __FUNCT__ 559b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 569b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 579b54502bSHong Zhang { 589b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 599b54502bSHong Zhang PetscErrorCode ierr; 60ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 619b54502bSHong Zhang PetscReal tol; 629b54502bSHong Zhang 639b54502bSHong Zhang PetscFunctionBegin; 649b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 658ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 665c9eb25fSBarry Smith 672401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 689b54502bSHong Zhang if (flg) { 699b54502bSHong Zhang tol = PETSC_DECIDE; 702401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 712401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 729b54502bSHong Zhang } 739b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 749b54502bSHong Zhang PetscFunctionReturn(0); 759b54502bSHong Zhang } 769b54502bSHong Zhang 779b54502bSHong Zhang #undef __FUNCT__ 789b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 799b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 809b54502bSHong Zhang { 819b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 829b54502bSHong Zhang PetscErrorCode ierr; 83ace3abfcSBarry Smith PetscBool iascii; 849b54502bSHong Zhang 859b54502bSHong Zhang PetscFunctionBegin; 862692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 879b54502bSHong Zhang if (iascii) { 88914a5d51SHong Zhang if (lu->inplace) { 89914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 90914a5d51SHong Zhang } else { 91914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 92145b9266SHong Zhang } 93145b9266SHong Zhang 949b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 959b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 969b54502bSHong Zhang } 97914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 989b54502bSHong Zhang PetscFunctionReturn(0); 999b54502bSHong Zhang } 1009b54502bSHong Zhang 1019b54502bSHong Zhang #undef __FUNCT__ 1029b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 1039b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 1049b54502bSHong Zhang { 1059b54502bSHong Zhang PetscErrorCode ierr; 1069b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1079b54502bSHong Zhang 1089b54502bSHong Zhang PetscFunctionBegin; 109075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1109b54502bSHong Zhang 1119b54502bSHong Zhang if (dir->inplace) { 1129b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1139b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 114075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 11503c60df9SBarry Smith if (dir->row) { 11603c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 11703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 11803c60df9SBarry Smith } 119075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 120075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1219b54502bSHong Zhang } else { 1229b54502bSHong Zhang MatInfo info; 1239b54502bSHong Zhang if (!pc->setupcalled) { 124075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1259b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1269b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1279b54502bSHong Zhang } 12803c60df9SBarry Smith if (dir->row) { 12903c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 13003c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 13103c60df9SBarry Smith } 132d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact){ 133075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 134a1f19f5aSHong Zhang } 135075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 136075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1379b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 138075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1399b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1409b54502bSHong Zhang if (!dir->reuseordering) { 1419b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1429b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 143075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1449b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1459b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1469b54502bSHong Zhang } 14703c60df9SBarry Smith if (dir->row) { 14803c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 14903c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 15003c60df9SBarry Smith } 1519b54502bSHong Zhang } 152075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 153075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 154075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 155075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1569b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 157075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1589b54502bSHong Zhang } 159075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1609b54502bSHong Zhang } 1619b54502bSHong Zhang PetscFunctionReturn(0); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang 1649b54502bSHong Zhang #undef __FUNCT__ 165574deadeSBarry Smith #define __FUNCT__ "PCReset_LU" 166574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1679b54502bSHong Zhang { 1689b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1699b54502bSHong Zhang PetscErrorCode ierr; 1709b54502bSHong Zhang 1719b54502bSHong Zhang PetscFunctionBegin; 172075768bcSBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 1739b54502bSHong Zhang if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 1749b54502bSHong Zhang if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 175574deadeSBarry Smith PetscFunctionReturn(0); 176574deadeSBarry Smith } 177574deadeSBarry Smith 178574deadeSBarry Smith #undef __FUNCT__ 179574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU" 180574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 181574deadeSBarry Smith { 182574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 183574deadeSBarry Smith PetscErrorCode ierr; 184574deadeSBarry Smith 185574deadeSBarry Smith PetscFunctionBegin; 186574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 187503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 188503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 189c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1909b54502bSHong Zhang PetscFunctionReturn(0); 1919b54502bSHong Zhang } 1929b54502bSHong Zhang 1939b54502bSHong Zhang #undef __FUNCT__ 1949b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 1959b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1969b54502bSHong Zhang { 1979b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1989b54502bSHong Zhang PetscErrorCode ierr; 1999b54502bSHong Zhang 2009b54502bSHong Zhang PetscFunctionBegin; 2019b54502bSHong Zhang if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 202075768bcSBarry Smith else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2039b54502bSHong Zhang PetscFunctionReturn(0); 2049b54502bSHong Zhang } 2059b54502bSHong Zhang 2069b54502bSHong Zhang #undef __FUNCT__ 2079b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2089b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2099b54502bSHong Zhang { 2109b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2119b54502bSHong Zhang PetscErrorCode ierr; 2129b54502bSHong Zhang 2139b54502bSHong Zhang PetscFunctionBegin; 2149b54502bSHong Zhang if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 215075768bcSBarry Smith else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 2169b54502bSHong Zhang PetscFunctionReturn(0); 2179b54502bSHong Zhang } 2189b54502bSHong Zhang 2199b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2209b54502bSHong Zhang 2219b54502bSHong Zhang EXTERN_C_BEGIN 2229b54502bSHong Zhang #undef __FUNCT__ 2232401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2247087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 2259b54502bSHong Zhang { 226c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2279b54502bSHong Zhang 2289b54502bSHong Zhang PetscFunctionBegin; 2299b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2309b54502bSHong Zhang PetscFunctionReturn(0); 2319b54502bSHong Zhang } 2329b54502bSHong Zhang EXTERN_C_END 2339b54502bSHong Zhang 2349b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2359b54502bSHong Zhang 2369b54502bSHong Zhang /*MC 2379b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2389b54502bSHong Zhang 2399b54502bSHong Zhang Options Database Keys: 2402401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 241c7393fdbSBarry Smith . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 2422401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 24355ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2442401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2452401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2462401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2479b54502bSHong Zhang stability of factorization. 248145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 249145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 250e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 251e22d95b2SBarry Smith diagonal. 2529b54502bSHong Zhang 2539b54502bSHong Zhang Notes: Not all options work for all matrix formats 2549b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2559b54502bSHong Zhang algorithms 2569b54502bSHong Zhang 2579b54502bSHong Zhang Level: beginner 2589b54502bSHong Zhang 2599b54502bSHong Zhang Concepts: LU factorization, direct solver 2609b54502bSHong Zhang 2619b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2629b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2639b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2649b54502bSHong Zhang 2659b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 266a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2678ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 268145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2698ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2709b54502bSHong Zhang M*/ 2719b54502bSHong Zhang 2729b54502bSHong Zhang EXTERN_C_BEGIN 2739b54502bSHong Zhang #undef __FUNCT__ 2749b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 2757087cfbeSBarry Smith PetscErrorCode PCCreate_LU(PC pc) 2769b54502bSHong Zhang { 2779b54502bSHong Zhang PetscErrorCode ierr; 2789b54502bSHong Zhang PetscMPIInt size; 2799b54502bSHong Zhang PC_LU *dir; 2809b54502bSHong Zhang 2819b54502bSHong Zhang PetscFunctionBegin; 28238f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 2839b54502bSHong Zhang 284075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 285879e8a4dSBarry Smith ((PC_Factor*)dir)->fact = PETSC_NULL; 286d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2879b54502bSHong Zhang dir->inplace = PETSC_FALSE; 2889b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2899b54502bSHong Zhang 290075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 291075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 292f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 293915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 294075768bcSBarry Smith ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 295075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 2969b54502bSHong Zhang dir->col = 0; 2979b54502bSHong Zhang dir->row = 0; 2985c9eb25fSBarry Smith 299a1f19f5aSHong Zhang ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 3007adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 3019b54502bSHong Zhang if (size == 1) { 3022692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3039b54502bSHong Zhang } else { 3042692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3059b54502bSHong Zhang } 3069b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3079b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3089b54502bSHong Zhang pc->data = (void*)dir; 3099b54502bSHong Zhang 310574deadeSBarry Smith pc->ops->reset = PCReset_LU; 3119b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3129b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3139b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3149b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3159b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3169b54502bSHong Zhang pc->ops->view = PCView_LU; 3179b54502bSHong Zhang pc->ops->applyrichardson = 0; 31885317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3199b54502bSHong Zhang 320f8260c8fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 321f8260c8fSBarry Smith PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 3227112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3237112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 32485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 32585317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 32685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 32785317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 328d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 329d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 330d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 331d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 33285317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 33385317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3342401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 3352401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 33685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 33785317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3382401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 3392401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 3402401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 3412401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 3428ff23777SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 3438ff23777SHong Zhang PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 34485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 34585317021SBarry Smith PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 3462401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 3472401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3489b54502bSHong Zhang PetscFunctionReturn(0); 3499b54502bSHong Zhang } 3509b54502bSHong Zhang EXTERN_C_END 351