19b54502bSHong Zhang 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 */ 7c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 89b54502bSHong Zhang 99b54502bSHong Zhang typedef struct { 10075768bcSBarry Smith PC_Factor hdr; 119b54502bSHong Zhang IS row,col; /* index sets used for reordering */ 129b54502bSHong Zhang } PC_Cholesky; 139b54502bSHong Zhang 144416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Cholesky(PetscOptionItems *PetscOptionsObject,PC pc) 159b54502bSHong Zhang { 169b54502bSHong Zhang PetscErrorCode ierr; 179b54502bSHong Zhang 189b54502bSHong Zhang PetscFunctionBegin; 19e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Cholesky options");CHKERRQ(ierr); 20e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 219b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 229b54502bSHong Zhang PetscFunctionReturn(0); 239b54502bSHong Zhang } 249b54502bSHong Zhang 259b54502bSHong Zhang static PetscErrorCode PCSetUp_Cholesky(PC pc) 269b54502bSHong Zhang { 279b54502bSHong Zhang PetscErrorCode ierr; 28ace3abfcSBarry Smith PetscBool flg; 299b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 30ea799195SBarry Smith MatSolverType stype; 3100e125f8SBarry Smith MatFactorError err; 329b54502bSHong Zhang 339b54502bSHong Zhang PetscFunctionBegin; 34c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 353d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill; 369b54502bSHong Zhang 3784d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 383d1c1ea0SBarry Smith if (dir->hdr.inplace) { 399b54502bSHong Zhang if (dir->row && dir->col && (dir->row != dir->col)) { 40fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 419b54502bSHong Zhang } 42fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 43075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 449b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 45fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 469b54502bSHong Zhang } 473bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 48075768bcSBarry Smith ierr = MatCholeskyFactor(pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 4900e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 5000e125f8SBarry Smith if (err) { /* Factor() fails */ 5100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 526baea169SHong Zhang PetscFunctionReturn(0); 536baea169SHong Zhang } 542fa5cd67SKarl Rupp 55075768bcSBarry Smith ((PC_Factor*)dir)->fact = pc->pmat; 569b54502bSHong Zhang } else { 579b54502bSHong Zhang MatInfo info; 5800e125f8SBarry Smith 599b54502bSHong Zhang if (!pc->setupcalled) { 60075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 619bfd6278SHong Zhang /* check if dir->row == dir->col */ 629bfd6278SHong Zhang ierr = ISEqual(dir->row,dir->col,&flg);CHKERRQ(ierr); 63e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row and column permutations must equal"); 64fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only pass one ordering into CholeskyFactor */ 659bfd6278SHong Zhang 6690d69ab7SBarry Smith flg = PETSC_FALSE; 67c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 689b54502bSHong Zhang if (flg) { 699b54502bSHong Zhang PetscReal tol = 1.e-10; 70c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 719b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 729b54502bSHong Zhang } 733bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 74d09a07f4SBarry Smith if (!((PC_Factor*)dir)->fact) { 75075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 76a1f19f5aSHong Zhang } 77075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 78075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 793d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 803bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 819b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 823d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 83fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 84075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 859d61c4f5SHong Zhang ierr = ISDestroy(&dir->col);CHKERRQ(ierr); /* only use dir->row ordering in CholeskyFactor */ 869d61c4f5SHong Zhang 8790d69ab7SBarry Smith flg = PETSC_FALSE; 88c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg,NULL);CHKERRQ(ierr); 899b54502bSHong Zhang if (flg) { 909b54502bSHong Zhang PetscReal tol = 1.e-10; 91c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,NULL);CHKERRQ(ierr); 929b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);CHKERRQ(ierr); 939b54502bSHong Zhang } 943bb1ff40SBarry Smith if (dir->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);} 959b54502bSHong Zhang } 966bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 97075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 98075768bcSBarry Smith ierr = MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 99075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1003d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 1013bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 10204545d6dSBarry Smith } else { 103b8b68cfdSBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 104160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 105b8b68cfdSBarry Smith ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 106b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 10704545d6dSBarry Smith } 1089b54502bSHong Zhang } 10900e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 11000e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 11100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1126baea169SHong Zhang PetscFunctionReturn(0); 1136baea169SHong Zhang } 1146baea169SHong Zhang 115075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 11600e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 11700e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 11800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1196baea169SHong Zhang } 1209b54502bSHong Zhang } 12100c67f3bSHong Zhang 1223ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 12300c67f3bSHong Zhang if (!stype) { 124ea799195SBarry Smith MatSolverType solverpackage; 1253ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 1263ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 12700c67f3bSHong Zhang } 1289b54502bSHong Zhang PetscFunctionReturn(0); 1299b54502bSHong Zhang } 1309b54502bSHong Zhang 13169d2c0f9SBarry Smith static PetscErrorCode PCReset_Cholesky(PC pc) 1329b54502bSHong Zhang { 1339b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1349b54502bSHong Zhang PetscErrorCode ierr; 1359b54502bSHong Zhang 1369b54502bSHong Zhang PetscFunctionBegin; 1373d1c1ea0SBarry Smith if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 138fcfd50ebSBarry Smith ierr = ISDestroy(&dir->row);CHKERRQ(ierr); 139fcfd50ebSBarry Smith ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 14069d2c0f9SBarry Smith PetscFunctionReturn(0); 14169d2c0f9SBarry Smith } 14269d2c0f9SBarry Smith 14369d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Cholesky(PC pc) 14469d2c0f9SBarry Smith { 14569d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky*)pc->data; 14669d2c0f9SBarry Smith PetscErrorCode ierr; 14769d2c0f9SBarry Smith 14869d2c0f9SBarry Smith PetscFunctionBegin; 14969d2c0f9SBarry Smith ierr = PCReset_Cholesky(pc);CHKERRQ(ierr); 150503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 151503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 152c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1539b54502bSHong Zhang PetscFunctionReturn(0); 1549b54502bSHong Zhang } 1559b54502bSHong Zhang 1569b54502bSHong Zhang static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y) 1579b54502bSHong Zhang { 1589b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1599b54502bSHong Zhang PetscErrorCode ierr; 1609b54502bSHong Zhang 1619b54502bSHong Zhang PetscFunctionBegin; 1623d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1632fa5cd67SKarl Rupp ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 1642fa5cd67SKarl Rupp } else { 1652fa5cd67SKarl Rupp ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1662fa5cd67SKarl Rupp } 1679b54502bSHong Zhang PetscFunctionReturn(0); 1689b54502bSHong Zhang } 1699b54502bSHong Zhang 1703d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y) 1713d72fe1eSJed Brown { 1723d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1733d72fe1eSJed Brown PetscErrorCode ierr; 1743d72fe1eSJed Brown 1753d72fe1eSJed Brown PetscFunctionBegin; 1763d72fe1eSJed Brown if (dir->hdr.inplace) { 1773d72fe1eSJed Brown ierr = MatForwardSolve(pc->pmat,x,y);CHKERRQ(ierr); 1783d72fe1eSJed Brown } else { 1793d72fe1eSJed Brown ierr = MatForwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1803d72fe1eSJed Brown } 1813d72fe1eSJed Brown PetscFunctionReturn(0); 1823d72fe1eSJed Brown } 1833d72fe1eSJed Brown 1843d72fe1eSJed Brown static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y) 1853d72fe1eSJed Brown { 1863d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky*)pc->data; 1873d72fe1eSJed Brown PetscErrorCode ierr; 1883d72fe1eSJed Brown 1893d72fe1eSJed Brown PetscFunctionBegin; 1903d72fe1eSJed Brown if (dir->hdr.inplace) { 1913d72fe1eSJed Brown ierr = MatBackwardSolve(pc->pmat,x,y);CHKERRQ(ierr); 1923d72fe1eSJed Brown } else { 1933d72fe1eSJed Brown ierr = MatBackwardSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 1943d72fe1eSJed Brown } 1953d72fe1eSJed Brown PetscFunctionReturn(0); 1963d72fe1eSJed Brown } 1973d72fe1eSJed Brown 1989b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y) 1999b54502bSHong Zhang { 2009b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky*)pc->data; 2019b54502bSHong Zhang PetscErrorCode ierr; 2029b54502bSHong Zhang 2039b54502bSHong Zhang PetscFunctionBegin; 2043d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2052fa5cd67SKarl Rupp ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 2062fa5cd67SKarl Rupp } else { 2072fa5cd67SKarl Rupp ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 2082fa5cd67SKarl Rupp } 2099b54502bSHong Zhang PetscFunctionReturn(0); 2109b54502bSHong Zhang } 2119b54502bSHong Zhang 2129b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2139b54502bSHong Zhang 2149b54502bSHong Zhang /* -----------------------------------------------------------------------------------*/ 2159b54502bSHong Zhang 2169b54502bSHong Zhang /*@ 2172401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2189b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2199b54502bSHong Zhang following factors. 2209b54502bSHong Zhang 221ad4df100SBarry Smith Logically Collective on PC 2229b54502bSHong Zhang 2239b54502bSHong Zhang Input Parameters: 2249b54502bSHong Zhang + pc - the preconditioner context 2259b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 2269b54502bSHong Zhang 2279b54502bSHong Zhang Options Database Key: 2282401956bSBarry Smith . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2299b54502bSHong Zhang 2309b54502bSHong Zhang Level: intermediate 2319b54502bSHong Zhang 2322401956bSBarry Smith .seealso: PCFactorSetReuseFill() 2339b54502bSHong Zhang @*/ 2347087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscBool flag) 2359b54502bSHong Zhang { 2364ac538c5SBarry Smith PetscErrorCode ierr; 2379b54502bSHong Zhang 2389b54502bSHong Zhang PetscFunctionBegin; 2390700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 240acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 2414ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCFactorSetReuseOrdering_C",(PC,PetscBool),(pc,flag));CHKERRQ(ierr); 2429b54502bSHong Zhang PetscFunctionReturn(0); 2439b54502bSHong Zhang } 2449b54502bSHong Zhang 2459b54502bSHong Zhang /*MC 24696fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2479b54502bSHong Zhang 2489b54502bSHong Zhang Options Database Keys: 2492401956bSBarry Smith + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 2503ca39a21SBarry Smith . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 2512401956bSBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 25255ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2532401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 254145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2559b54502bSHong Zhang 25695452b02SPatrick Sanan Notes: 25795452b02SPatrick Sanan Not all options work for all matrix formats 2589b54502bSHong Zhang 2599b54502bSHong Zhang Level: beginner 2609b54502bSHong Zhang 26195452b02SPatrick Sanan Notes: 26295452b02SPatrick Sanan 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, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 268145b9266SHong Zhang PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftType(), PCFactorSetShiftAmount() 2698e37d05fSBarry Smith PCFactorSetUseInPlace(), PCFactorGetUseInPlace(), PCFactorSetMatOrderingType() 2709b54502bSHong Zhang 2719b54502bSHong Zhang M*/ 2729b54502bSHong Zhang 2738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 2749b54502bSHong Zhang { 2759b54502bSHong Zhang PetscErrorCode ierr; 2769b54502bSHong Zhang PC_Cholesky *dir; 2779b54502bSHong Zhang 2789b54502bSHong Zhang PetscFunctionBegin; 279b00a9115SJed Brown ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 2803d1c1ea0SBarry Smith pc->data = (void*)dir; 2813d1c1ea0SBarry Smith ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 2822fa5cd67SKarl Rupp 283879e8a4dSBarry Smith ((PC_Factor*)dir)->factortype = MAT_FACTOR_CHOLESKY; 284075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = 5.0; 2852fa5cd67SKarl Rupp 2869b54502bSHong Zhang dir->col = 0; 2879b54502bSHong Zhang dir->row = 0; 2882fa5cd67SKarl Rupp 289*9bd791bbSBarry Smith /* MATORDERINGNATURAL_OR_ND allows selecting type based on matrix type sbaij or aij */ 290*9bd791bbSBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL_OR_ND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 2919b54502bSHong Zhang 2929b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 29369d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 2949b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 2953d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 2963d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 2979b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 2989b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 2999b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 30092e08861SBarry Smith pc->ops->view = PCView_Factor; 3019b54502bSHong Zhang pc->ops->applyrichardson = 0; 3029b54502bSHong Zhang PetscFunctionReturn(0); 3039b54502bSHong Zhang } 304