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 #undef __FUNCT__ 112401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 127087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 139b54502bSHong Zhang { 149b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 159b54502bSHong Zhang 169b54502bSHong Zhang PetscFunctionBegin; 179b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 182fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 192fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 209b54502bSHong Zhang PetscFunctionReturn(0); 219b54502bSHong Zhang } 229b54502bSHong Zhang 239b54502bSHong Zhang #undef __FUNCT__ 242401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 257087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag) 269b54502bSHong Zhang { 27c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 289b54502bSHong Zhang 299b54502bSHong Zhang PetscFunctionBegin; 309b54502bSHong Zhang lu->reuseordering = flag; 319b54502bSHong Zhang PetscFunctionReturn(0); 329b54502bSHong Zhang } 339b54502bSHong Zhang 349b54502bSHong Zhang #undef __FUNCT__ 352401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_LU" 367087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag) 379b54502bSHong Zhang { 38c62fe420SBarry Smith PC_LU *lu = (PC_LU*)pc->data; 399b54502bSHong Zhang 409b54502bSHong Zhang PetscFunctionBegin; 419b54502bSHong Zhang lu->reusefill = flag; 429b54502bSHong Zhang PetscFunctionReturn(0); 439b54502bSHong Zhang } 449b54502bSHong Zhang 459b54502bSHong Zhang #undef __FUNCT__ 469b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_LU" 479b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_LU(PC pc) 489b54502bSHong Zhang { 499b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 509b54502bSHong Zhang PetscErrorCode ierr; 51ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 529b54502bSHong Zhang PetscReal tol; 539b54502bSHong Zhang 549b54502bSHong Zhang PetscFunctionBegin; 559b54502bSHong Zhang ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 568ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 575c9eb25fSBarry Smith 582401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 599b54502bSHong Zhang if (flg) { 609b54502bSHong Zhang tol = PETSC_DECIDE; 612401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 622401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 639b54502bSHong Zhang } 649b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 659b54502bSHong Zhang PetscFunctionReturn(0); 669b54502bSHong Zhang } 679b54502bSHong Zhang 689b54502bSHong Zhang #undef __FUNCT__ 699b54502bSHong Zhang #define __FUNCT__ "PCView_LU" 709b54502bSHong Zhang static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 719b54502bSHong Zhang { 729b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 739b54502bSHong Zhang PetscErrorCode ierr; 74ace3abfcSBarry Smith PetscBool iascii; 759b54502bSHong Zhang 769b54502bSHong Zhang PetscFunctionBegin; 77251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 789b54502bSHong Zhang if (iascii) { 79914a5d51SHong Zhang if (lu->inplace) { 80914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 81914a5d51SHong Zhang } else { 82914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 83145b9266SHong Zhang } 84145b9266SHong Zhang 859b54502bSHong Zhang if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 869b54502bSHong Zhang if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 879b54502bSHong Zhang } 88914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 899b54502bSHong Zhang PetscFunctionReturn(0); 909b54502bSHong Zhang } 919b54502bSHong Zhang 929b54502bSHong Zhang #undef __FUNCT__ 939b54502bSHong Zhang #define __FUNCT__ "PCSetUp_LU" 949b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 959b54502bSHong Zhang { 969b54502bSHong Zhang PetscErrorCode ierr; 979b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 989b54502bSHong Zhang 999b54502bSHong Zhang PetscFunctionBegin; 100075768bcSBarry Smith if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 1019b54502bSHong Zhang 1029b54502bSHong Zhang if (dir->inplace) { 103fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 104fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 105075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 10603c60df9SBarry Smith if (dir->row) { 107*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 108*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 10903c60df9SBarry Smith } 110075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1112fa5cd67SKarl Rupp 112075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 1139b54502bSHong Zhang } else { 1149b54502bSHong Zhang MatInfo info; 1159b54502bSHong Zhang if (!pc->setupcalled) { 116075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1179b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1189b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1199b54502bSHong Zhang } 12003c60df9SBarry Smith if (dir->row) { 121*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 122*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 12303c60df9SBarry Smith } 124d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 125075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 126a1f19f5aSHong Zhang } 127075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 128075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1299b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 130*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1319b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1329b54502bSHong Zhang if (!dir->reuseordering) { 133fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 134fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 135075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1369b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1379b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1389b54502bSHong Zhang } 13903c60df9SBarry Smith if (dir->row) { 140*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 141*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 14203c60df9SBarry Smith } 1439b54502bSHong Zhang } 1446bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 145075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 146075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 147075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1489b54502bSHong Zhang dir->actualfill = info.fill_ratio_needed; 149*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 1509b54502bSHong Zhang } 15100c46b4dSJed Brown if ((!pc->setupcalled) || pc->setupcalled) { 152075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 1539b54502bSHong Zhang } 154b4813acdSShri Abhyankar } 1559b54502bSHong Zhang PetscFunctionReturn(0); 1569b54502bSHong Zhang } 1579b54502bSHong Zhang 1589b54502bSHong Zhang #undef __FUNCT__ 159574deadeSBarry Smith #define __FUNCT__ "PCReset_LU" 160574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1619b54502bSHong Zhang { 1629b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1639b54502bSHong Zhang PetscErrorCode ierr; 1649b54502bSHong Zhang 1659b54502bSHong Zhang PetscFunctionBegin; 1666bf464f9SBarry Smith if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 167fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 168fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 169574deadeSBarry Smith PetscFunctionReturn(0); 170574deadeSBarry Smith } 171574deadeSBarry Smith 172574deadeSBarry Smith #undef __FUNCT__ 173574deadeSBarry Smith #define __FUNCT__ "PCDestroy_LU" 174574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 175574deadeSBarry Smith { 176574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 177574deadeSBarry Smith PetscErrorCode ierr; 178574deadeSBarry Smith 179574deadeSBarry Smith PetscFunctionBegin; 180574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 181503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 182503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 183c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1849b54502bSHong Zhang PetscFunctionReturn(0); 1859b54502bSHong Zhang } 1869b54502bSHong Zhang 1879b54502bSHong Zhang #undef __FUNCT__ 1889b54502bSHong Zhang #define __FUNCT__ "PCApply_LU" 1899b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1909b54502bSHong Zhang { 1919b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1929b54502bSHong Zhang PetscErrorCode ierr; 1939b54502bSHong Zhang 1949b54502bSHong Zhang PetscFunctionBegin; 1952fa5cd67SKarl Rupp if (dir->inplace) { 1962fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1972fa5cd67SKarl Rupp } else { 1982fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1992fa5cd67SKarl Rupp } 2009b54502bSHong Zhang PetscFunctionReturn(0); 2019b54502bSHong Zhang } 2029b54502bSHong Zhang 2039b54502bSHong Zhang #undef __FUNCT__ 2049b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_LU" 2059b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2069b54502bSHong Zhang { 2079b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2089b54502bSHong Zhang PetscErrorCode ierr; 2099b54502bSHong Zhang 2109b54502bSHong Zhang PetscFunctionBegin; 2112fa5cd67SKarl Rupp if (dir->inplace) { 2122fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2132fa5cd67SKarl Rupp } else { 2142fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2152fa5cd67SKarl Rupp } 2169b54502bSHong Zhang PetscFunctionReturn(0); 2179b54502bSHong Zhang } 2189b54502bSHong Zhang 2199b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2209b54502bSHong Zhang 2219b54502bSHong Zhang #undef __FUNCT__ 2222401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_LU" 2237087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 2249b54502bSHong Zhang { 225c62fe420SBarry Smith PC_LU *dir = (PC_LU*)pc->data; 2269b54502bSHong Zhang 2279b54502bSHong Zhang PetscFunctionBegin; 2289b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2299b54502bSHong Zhang PetscFunctionReturn(0); 2309b54502bSHong Zhang } 2319b54502bSHong Zhang 2329b54502bSHong Zhang /* ------------------------------------------------------------------------ */ 2339b54502bSHong Zhang 2349b54502bSHong Zhang /*MC 2359b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2369b54502bSHong Zhang 2379b54502bSHong Zhang Options Database Keys: 2382401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 239d45987f3SHong Zhang . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 2402401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 24155ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2422401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2432401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2442401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2459b54502bSHong Zhang stability of factorization. 246145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 247145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 248e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 249e22d95b2SBarry Smith diagonal. 2509b54502bSHong Zhang 2519b54502bSHong Zhang Notes: Not all options work for all matrix formats 2529b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2539b54502bSHong Zhang algorithms 2549b54502bSHong Zhang 2559b54502bSHong Zhang Level: beginner 2569b54502bSHong Zhang 2579b54502bSHong Zhang Concepts: LU factorization, direct solver 2589b54502bSHong Zhang 2599b54502bSHong Zhang Notes: Usually this will compute an "exact" solution in one iteration and does 2609b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2619b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2629b54502bSHong Zhang 2639b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 264a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2658ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 266145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2678ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2689b54502bSHong Zhang M*/ 2699b54502bSHong Zhang 2709b54502bSHong Zhang #undef __FUNCT__ 2719b54502bSHong Zhang #define __FUNCT__ "PCCreate_LU" 2728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 2739b54502bSHong Zhang { 2749b54502bSHong Zhang PetscErrorCode ierr; 2759b54502bSHong Zhang PetscMPIInt size; 2769b54502bSHong Zhang PC_LU *dir; 2779b54502bSHong Zhang 2789b54502bSHong Zhang PetscFunctionBegin; 27938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 2809b54502bSHong Zhang 281075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 2822fa5cd67SKarl Rupp 2830298fd71SBarry Smith ((PC_Factor*)dir)->fact = NULL; 284d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2859b54502bSHong Zhang dir->inplace = PETSC_FALSE; 2869b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2879b54502bSHong Zhang 288075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 289075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 290f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 291915743fcSHong Zhang ((PC_Factor*)dir)->info.shiftamount = 0.0; 2920ed735ceSHong Zhang ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 293075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 2949b54502bSHong Zhang dir->col = 0; 2959b54502bSHong Zhang dir->row = 0; 2965c9eb25fSBarry Smith 297a1f19f5aSHong Zhang ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 298ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 2999b54502bSHong Zhang if (size == 1) { 30019fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3019b54502bSHong Zhang } else { 30219fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 3039b54502bSHong Zhang } 3049b54502bSHong Zhang dir->reusefill = PETSC_FALSE; 3059b54502bSHong Zhang dir->reuseordering = PETSC_FALSE; 3069b54502bSHong Zhang pc->data = (void*)dir; 3079b54502bSHong Zhang 308574deadeSBarry Smith pc->ops->reset = PCReset_LU; 3099b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 3109b54502bSHong Zhang pc->ops->apply = PCApply_LU; 3119b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 3129b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 3139b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 3149b54502bSHong Zhang pc->ops->view = PCView_LU; 3159b54502bSHong Zhang pc->ops->applyrichardson = 0; 31685317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3179b54502bSHong Zhang 318bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 319bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 320bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 321bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 322bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 323bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 324bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 325bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 326bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 327bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 328bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr); 329bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 330bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 331bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 3329b54502bSHong Zhang PetscFunctionReturn(0); 3339b54502bSHong Zhang } 334