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 117087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 129b54502bSHong Zhang { 139b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 149b54502bSHong Zhang 159b54502bSHong Zhang PetscFunctionBegin; 169b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE; 172fa5cd67SKarl Rupp if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 182fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z; 199b54502bSHong Zhang PetscFunctionReturn(0); 209b54502bSHong Zhang } 219b54502bSHong Zhang 224416b707SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc) 239b54502bSHong Zhang { 249b54502bSHong Zhang PC_LU *lu = (PC_LU*)pc->data; 259b54502bSHong Zhang PetscErrorCode ierr; 26ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 279b54502bSHong Zhang PetscReal tol; 289b54502bSHong Zhang 299b54502bSHong Zhang PetscFunctionBegin; 30e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr); 31e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 325c9eb25fSBarry Smith 332401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 349b54502bSHong Zhang if (flg) { 359b54502bSHong Zhang tol = PETSC_DECIDE; 360a545947SLisandro Dalcin ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL);CHKERRQ(ierr); 372401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 389b54502bSHong Zhang } 399b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 409b54502bSHong Zhang PetscFunctionReturn(0); 419b54502bSHong Zhang } 429b54502bSHong Zhang 439b54502bSHong Zhang static PetscErrorCode PCSetUp_LU(PC pc) 449b54502bSHong Zhang { 459b54502bSHong Zhang PetscErrorCode ierr; 469b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 47ea799195SBarry Smith MatSolverType stype; 4800e125f8SBarry Smith MatFactorError err; 493d1c1ea0SBarry Smith 509b54502bSHong Zhang PetscFunctionBegin; 51c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 523d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill; 539b54502bSHong Zhang 5484d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 553d1c1ea0SBarry Smith if (dir->hdr.inplace) { 569d76b4d0SMatthew G. Knepley MatFactorType ftype; 579d76b4d0SMatthew G. Knepley 589d76b4d0SMatthew G. Knepley ierr = MatGetFactorType(pc->pmat, &ftype);CHKERRQ(ierr); 599d76b4d0SMatthew G. Knepley if (ftype == MAT_FACTOR_NONE) { 60fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 61fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 62*2c7c0729SBarry Smith /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 63075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 6403c60df9SBarry Smith if (dir->row) { 653bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 663bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 6703c60df9SBarry Smith } 68075768bcSBarry Smith ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 6900e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 7000e125f8SBarry Smith if (err) { /* Factor() fails */ 7100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 726baea169SHong Zhang PetscFunctionReturn(0); 736baea169SHong Zhang } 749d76b4d0SMatthew G. Knepley } 75075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 769b54502bSHong Zhang } else { 779b54502bSHong Zhang MatInfo info; 7800e125f8SBarry Smith 799b54502bSHong Zhang if (!pc->setupcalled) { 80*2c7c0729SBarry Smith PetscBool useordering; 81*2c7c0729SBarry Smith if (!((PC_Factor*)dir)->fact) { 82*2c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 83*2c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 84*2c7c0729SBarry Smith } 85*2c7c0729SBarry Smith ierr = MatFactorGetUseOrdering(((PC_Factor*)dir)->fact,&useordering);CHKERRQ(ierr); 86*2c7c0729SBarry Smith if (useordering) { 87075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 889b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 899b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 909b54502bSHong Zhang } 913bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 923bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 9303c60df9SBarry Smith } 94075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 95075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 963d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 979b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 98*2c7c0729SBarry Smith PetscBool useordering; 993d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 100*2c7c0729SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 101*2c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 102*2c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 103*2c7c0729SBarry Smith ierr = MatFactorGetUseOrdering(((PC_Factor*)dir)->fact,&useordering);CHKERRQ(ierr); 104*2c7c0729SBarry Smith if (useordering) { 105fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 106fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 107075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 1089b54502bSHong Zhang if (dir->nonzerosalongdiagonal) { 1099b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 1109b54502bSHong Zhang } 1113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 1123bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 11303c60df9SBarry Smith } 1149b54502bSHong Zhang } 115075768bcSBarry Smith ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 116075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1173d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 11804545d6dSBarry Smith } else { 119b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 120160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 121b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 122b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 12304545d6dSBarry Smith } 1249b54502bSHong Zhang } 12500e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 12600e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12700e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1288c1cd74cSHong Zhang PetscFunctionReturn(0); 1298c1cd74cSHong Zhang } 1308c1cd74cSHong Zhang 131075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 13200e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 13300e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 13400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1358c1cd74cSHong Zhang } 136680c5173SHong Zhang 1379b54502bSHong Zhang } 13800c67f3bSHong Zhang 1393ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 14000c67f3bSHong Zhang if (!stype) { 141ea799195SBarry Smith MatSolverType solverpackage; 1423ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 1433ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 14400c67f3bSHong Zhang } 1459b54502bSHong Zhang PetscFunctionReturn(0); 1469b54502bSHong Zhang } 1479b54502bSHong Zhang 148574deadeSBarry Smith static PetscErrorCode PCReset_LU(PC pc) 1499b54502bSHong Zhang { 1509b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1519b54502bSHong Zhang PetscErrorCode ierr; 1529b54502bSHong Zhang 1539b54502bSHong Zhang PetscFunctionBegin; 1543d1c1ea0SBarry Smith if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 155fcfd50ebSBarry Smith if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 156fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 157574deadeSBarry Smith PetscFunctionReturn(0); 158574deadeSBarry Smith } 159574deadeSBarry Smith 160574deadeSBarry Smith static PetscErrorCode PCDestroy_LU(PC pc) 161574deadeSBarry Smith { 162574deadeSBarry Smith PC_LU *dir = (PC_LU*)pc->data; 163574deadeSBarry Smith PetscErrorCode ierr; 164574deadeSBarry Smith 165574deadeSBarry Smith PetscFunctionBegin; 166574deadeSBarry Smith ierr = PCReset_LU(pc);CHKERRQ(ierr); 167503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 168503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 169c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1709b54502bSHong Zhang PetscFunctionReturn(0); 1719b54502bSHong Zhang } 1729b54502bSHong Zhang 1739b54502bSHong Zhang static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 1749b54502bSHong Zhang { 1759b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 1769b54502bSHong Zhang PetscErrorCode ierr; 1779b54502bSHong Zhang 1789b54502bSHong Zhang PetscFunctionBegin; 1793d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1802fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1812fa5cd67SKarl Rupp } else { 1822fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1832fa5cd67SKarl Rupp } 1849b54502bSHong Zhang PetscFunctionReturn(0); 1859b54502bSHong Zhang } 1869b54502bSHong Zhang 1877b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y) 1887b6e2003SPierre Jolivet { 1897b6e2003SPierre Jolivet PC_LU *dir = (PC_LU*)pc->data; 1907b6e2003SPierre Jolivet PetscErrorCode ierr; 1917b6e2003SPierre Jolivet 1927b6e2003SPierre Jolivet PetscFunctionBegin; 1937b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1947b6e2003SPierre Jolivet ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr); 1957b6e2003SPierre Jolivet } else { 1967b6e2003SPierre Jolivet ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr); 1977b6e2003SPierre Jolivet } 1987b6e2003SPierre Jolivet PetscFunctionReturn(0); 1997b6e2003SPierre Jolivet } 2007b6e2003SPierre Jolivet 2019b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 2029b54502bSHong Zhang { 2039b54502bSHong Zhang PC_LU *dir = (PC_LU*)pc->data; 2049b54502bSHong Zhang PetscErrorCode ierr; 2059b54502bSHong Zhang 2069b54502bSHong Zhang PetscFunctionBegin; 2073d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2082fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2092fa5cd67SKarl Rupp } else { 2102fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2112fa5cd67SKarl Rupp } 2129b54502bSHong Zhang PetscFunctionReturn(0); 2139b54502bSHong Zhang } 2149b54502bSHong Zhang 2159b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2169b54502bSHong Zhang 2179b54502bSHong Zhang /*MC 2189b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 2199b54502bSHong Zhang 2209b54502bSHong Zhang Options Database Keys: 2212401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2223ca39a21SBarry Smith . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 2232401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 22455ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2252401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2262401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2272401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 2289b54502bSHong Zhang stability of factorization. 229145b9266SHong Zhang . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 230145b9266SHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 231e22d95b2SBarry Smith - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 232e22d95b2SBarry Smith diagonal. 2339b54502bSHong Zhang 23495452b02SPatrick Sanan Notes: 23595452b02SPatrick Sanan Not all options work for all matrix formats 2369b54502bSHong Zhang Run with -help to see additional options for particular matrix formats or factorization 2379b54502bSHong Zhang algorithms 2389b54502bSHong Zhang 2399b54502bSHong Zhang Level: beginner 2409b54502bSHong Zhang 24195452b02SPatrick Sanan Notes: 24295452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2439b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 2449b54502bSHong Zhang KSPSetType(ksp,KSPPREONLY) for the Krylov method 2459b54502bSHong Zhang 2469b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 247a4fd02acSBarry Smith PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 2488ff23777SHong Zhang PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 249145b9266SHong Zhang PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 2508ff23777SHong Zhang PCFactorReorderForNonzeroDiagonal() 2519b54502bSHong Zhang M*/ 2529b54502bSHong Zhang 2538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 2549b54502bSHong Zhang { 2559b54502bSHong Zhang PetscErrorCode ierr; 2569b54502bSHong Zhang PetscMPIInt size; 2579b54502bSHong Zhang PC_LU *dir; 2589b54502bSHong Zhang 2599b54502bSHong Zhang PetscFunctionBegin; 260b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 2613d1c1ea0SBarry Smith pc->data = (void*)dir; 2623d1c1ea0SBarry Smith ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 263d5f3da31SBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 2649b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE; 2659b54502bSHong Zhang 266075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 267075768bcSBarry Smith ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 268f4db908eSBarry Smith ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 2690a545947SLisandro Dalcin dir->col = NULL; 2700a545947SLisandro Dalcin dir->row = NULL; 2715c9eb25fSBarry Smith 272ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 2739b54502bSHong Zhang if (size == 1) { 27419fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 2759b54502bSHong Zhang } else { 27619fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 2779b54502bSHong Zhang } 2789b54502bSHong Zhang 279574deadeSBarry Smith pc->ops->reset = PCReset_LU; 2809b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU; 2819b54502bSHong Zhang pc->ops->apply = PCApply_LU; 2827b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_LU; 2839b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU; 2849b54502bSHong Zhang pc->ops->setup = PCSetUp_LU; 2859b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU; 28692e08861SBarry Smith pc->ops->view = PCView_Factor; 2870a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 288bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 2899b54502bSHong Zhang PetscFunctionReturn(0); 2909b54502bSHong Zhang } 291