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 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; 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 EXTERN_C_END 249b54502bSHong Zhang 259b54502bSHong Zhang EXTERN_C_BEGIN 269b54502bSHong Zhang #undef __FUNCT__ 272401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 287087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag) 299b54502bSHong Zhang { 30c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 319b54502bSHong Zhang 329b54502bSHong Zhang PetscFunctionBegin; 339b54502bSHong Zhang lu->reuseordering = flag; 349b54502bSHong Zhang PetscFunctionReturn(0); 359b54502bSHong Zhang } 369b54502bSHong Zhang EXTERN_C_END 379b54502bSHong Zhang 389b54502bSHong Zhang EXTERN_C_BEGIN 399b54502bSHong Zhang #undef __FUNCT__ 402401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 417087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag) 429b54502bSHong Zhang { 43c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 449b54502bSHong Zhang 459b54502bSHong Zhang PetscFunctionBegin; 469b54502bSHong Zhang lu->reusefill = flag; 479b54502bSHong Zhang PetscFunctionReturn(0); 489b54502bSHong Zhang } 499b54502bSHong Zhang EXTERN_C_END 509b54502bSHong Zhang 519b54502bSHong Zhang #undef __FUNCT__ 529b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 539b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 549b54502bSHong Zhang { 559b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 569b54502bSHong Zhang PetscErrorCode ierr; 57ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 589b54502bSHong Zhang PetscReal tol; 599b54502bSHong Zhang 609b54502bSHong Zhang PetscFunctionBegin; 619b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 628ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 635c9eb25fSBarry Smith 642401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 659b54502bSHong Zhang if (flg) { 669b54502bSHong Zhang tol = PETSC_DECIDE; 672401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 682401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 699b54502bSHong Zhang } 709b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 719b54502bSHong Zhang PetscFunctionReturn(0); 729b54502bSHong Zhang } 739b54502bSHong Zhang 749b54502bSHong Zhang #undef __FUNCT__ 759b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 769b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 779b54502bSHong Zhang { 789b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 799b54502bSHong Zhang PetscErrorCode ierr; 80ace3abfcSBarry Smith PetscBool iascii; 819b54502bSHong Zhang 829b54502bSHong Zhang PetscFunctionBegin; 83251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 849b54502bSHong Zhang if (iascii) { 85914a5d51SHong Zhang if (lu->inplace) { 86914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 87914a5d51SHong Zhang } else { 88914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 89145b9266SHong Zhang } 90145b9266SHong Zhang 919b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 929b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 939b54502bSHong Zhang } 94914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 959b54502bSHong Zhang PetscFunctionReturn(0); 969b54502bSHong Zhang } 979b54502bSHong Zhang 989b54502bSHong Zhang #undef __FUNCT__ 999b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 1009b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 1019b54502bSHong Zhang { 1029b54502bSHong Zhang PetscErrorCode ierr; 1039b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1049b54502bSHong Zhang 1059b54502bSHong Zhang PetscFunctionBegin; 106075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1079b54502bSHong Zhang 1089b54502bSHong Zhang if (dir->inplace) { 109fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 110fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 111075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 11203c60df9SBarry Smith if (dir->row) { 11303c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 11403c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 11503c60df9SBarry Smith } 116075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1172fa5cd67SKarl Rupp 118075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1199b54502bSHong Zhang } else { 1209b54502bSHong Zhang MatInfo info; 1219b54502bSHong Zhang if (!pc->setupcalled) { 122075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1239b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1249b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1259b54502bSHong Zhang } 12603c60df9SBarry Smith if (dir->row) { 12703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 12803c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 12903c60df9SBarry Smith } 130d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 131075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 132a1f19f5aSHong Zhang } 133075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 134075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1359b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 136075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1379b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1389b54502bSHong Zhang if (!dir->reuseordering) { 139fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 140fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 141075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1429b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1439b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1449b54502bSHong Zhang } 14503c60df9SBarry Smith if (dir->row) { 14603c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 14703c60df9SBarry Smith ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 14803c60df9SBarry Smith } 1499b54502bSHong Zhang } 1506bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 151075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 152075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 153075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1549b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 155075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1569b54502bSHong Zhang } 157b4813acdSShri Abhyankar if ((!pc->setupcalled) || (pc->setupcalled && !pc->reuse)) { 158075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1599b54502bSHong Zhang } 160b4813acdSShri Abhyankar } 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; 1726bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 173fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 174fcfd50ebSBarry Smith 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; 2012fa5cd67SKarl Rupp if (dir->inplace) { 2022fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 2032fa5cd67SKarl Rupp } else { 2042fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2052fa5cd67SKarl Rupp } 2069b54502bSHong Zhang PetscFunctionReturn(0); 2079b54502bSHong Zhang } 2089b54502bSHong Zhang 2099b54502bSHong Zhang #undef __FUNCT__ 2109b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2119b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2129b54502bSHong Zhang { 2139b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2149b54502bSHong Zhang PetscErrorCode ierr; 2159b54502bSHong Zhang 2169b54502bSHong Zhang PetscFunctionBegin; 2172fa5cd67SKarl Rupp if (dir->inplace) { 2182fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2192fa5cd67SKarl Rupp } else { 2202fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2212fa5cd67SKarl Rupp } 2229b54502bSHong Zhang PetscFunctionReturn(0); 2239b54502bSHong Zhang } 2249b54502bSHong Zhang 2259b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2269b54502bSHong Zhang 2279b54502bSHong Zhang EXTERN_C_BEGIN 2289b54502bSHong Zhang #undef __FUNCT__ 2292401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2307087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 2319b54502bSHong Zhang { 232c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2339b54502bSHong Zhang 2349b54502bSHong Zhang PetscFunctionBegin; 2359b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2369b54502bSHong Zhang PetscFunctionReturn(0); 2379b54502bSHong Zhang } 2389b54502bSHong Zhang EXTERN_C_END 2399b54502bSHong Zhang 2409b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2419b54502bSHong Zhang 2429b54502bSHong Zhang /*MC 2439b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2449b54502bSHong Zhang 2459b54502bSHong Zhang Options Database Keys: 2462401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 247d45987f3SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2482401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 24955ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2502401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2512401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2522401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2539b54502bSHong Zhang stability of factorization. 254145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 255145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 256e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 257e22d95b2SBarry Smith diagonal. 2589b54502bSHong Zhang 2599b54502bSHong Zhang Notes: Not all options work for all matrix formats 2609b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2619b54502bSHong Zhang algorithms 2629b54502bSHong Zhang 2639b54502bSHong Zhang Level: beginner 2649b54502bSHong Zhang 2659b54502bSHong Zhang Concepts: LU factorization, direct solver 2669b54502bSHong Zhang 2679b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2689b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2699b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2709b54502bSHong Zhang 2719b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 272a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2738ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 274145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2758ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2769b54502bSHong Zhang M*/ 2779b54502bSHong Zhang 2789b54502bSHong Zhang EXTERN_C_BEGIN 2799b54502bSHong Zhang #undef __FUNCT__ 2809b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 2817087cfbeSBarry Smith PetscErrorCode PCCreate_LU(PC pc) 2829b54502bSHong Zhang { 2839b54502bSHong Zhang PetscErrorCode ierr; 2849b54502bSHong Zhang PetscMPIInt size; 2859b54502bSHong Zhang PC_LU *dir; 2869b54502bSHong Zhang 2879b54502bSHong Zhang PetscFunctionBegin; 28838f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 2899b54502bSHong Zhang 290075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 2912fa5cd67SKarl Rupp 2920298fd71SBarry Smith ((PC_Factor*)dir)->fact = NULL; 293d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2949b54502bSHong Zhang dir->inplace = PETSC_FALSE; 2959b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2969b54502bSHong Zhang 297075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 298075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 299f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 300915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 3010ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 302075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 3039b54502bSHong Zhang dir->col = 0; 3049b54502bSHong Zhang dir->row = 0; 3055c9eb25fSBarry Smith 306a1f19f5aSHong Zhang ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 307*ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 3089b54502bSHong Zhang if (size == 1) { 30919fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3109b54502bSHong Zhang } else { 31119fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3129b54502bSHong Zhang } 3139b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3149b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3159b54502bSHong Zhang pc->data = (void*)dir; 3169b54502bSHong Zhang 317574deadeSBarry Smith pc->ops->reset = PCReset_LU; 3189b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3199b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3209b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3219b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3229b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3239b54502bSHong Zhang pc->ops->view = PCView_LU; 3249b54502bSHong Zhang pc->ops->applyrichardson = 0; 32585317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3269b54502bSHong Zhang 327f8260c8fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 328f8260c8fSBarry Smith PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 3297112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 3307112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 33185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 33285317021SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 33385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 33485317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 335d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 336d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 337d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 338d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 33985317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 34085317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 3412401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 3422401956bSBarry Smith PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 34385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 34485317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 3452401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 3462401956bSBarry Smith PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 3472401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 3482401956bSBarry Smith PCFactorSetReuseFill_LU);CHKERRQ(ierr); 3498ff23777SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 3508ff23777SHong Zhang PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 35185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 35285317021SBarry Smith PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 3532401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 3542401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3559b54502bSHong Zhang PetscFunctionReturn(0); 3569b54502bSHong Zhang } 3579b54502bSHong Zhang EXTERN_C_END 358