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 14d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems *PetscOptionsObject) 15d71ae5a4SJacob Faibussowitsch { 169b54502bSHong Zhang PetscFunctionBegin; 17d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options"); 18dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject)); 19d0609cedSBarry Smith PetscOptionsHeadEnd(); 203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 219b54502bSHong Zhang } 229b54502bSHong Zhang 23d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Cholesky(PC pc) 24d71ae5a4SJacob Faibussowitsch { 25ace3abfcSBarry Smith PetscBool flg; 269b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 27ea799195SBarry Smith MatSolverType stype; 2800e125f8SBarry Smith MatFactorError err; 29f023e1d5SPierre Jolivet const char *prefix; 309b54502bSHong Zhang 319b54502bSHong Zhang PetscFunctionBegin; 32c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 333d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill; 349b54502bSHong Zhang 3526cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 3626cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 3726cc229bSBarry Smith 389566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 393d1c1ea0SBarry Smith if (dir->hdr.inplace) { 40fe143ba0SPablo Brubeck MatFactorType ftype; 41fe143ba0SPablo Brubeck 42fe143ba0SPablo Brubeck PetscCall(MatGetFactorType(pc->pmat, &ftype)); 43fe143ba0SPablo Brubeck if (ftype == MAT_FACTOR_NONE) { 4448a46eb9SPierre Jolivet if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row)); 459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 46fe143ba0SPablo Brubeck /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 479566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 489566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 499b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 519b54502bSHong Zhang } 529566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 539566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err)); 5400e125f8SBarry Smith if (err) { /* Factor() fails */ 5500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 576baea169SHong Zhang } 58fe143ba0SPablo Brubeck } 59075768bcSBarry Smith ((PC_Factor *)dir)->fact = pc->pmat; 609b54502bSHong Zhang } else { 619b54502bSHong Zhang MatInfo info; 6200e125f8SBarry Smith 639b54502bSHong Zhang if (!pc->setupcalled) { 64f73b0415SBarry Smith PetscBool canuseordering; 65*03e5aca4SStefano Zampini 66*03e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 679566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 68f73b0415SBarry Smith if (canuseordering) { 699566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 709566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 719bfd6278SHong Zhang /* check if dir->row == dir->col */ 724ac6704cSBarry Smith if (dir->row) { 739566063dSJacob Faibussowitsch PetscCall(ISEqual(dir->row, dir->col, &flg)); 7428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal"); 754ac6704cSBarry Smith } 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */ 779bfd6278SHong Zhang 7890d69ab7SBarry Smith flg = PETSC_FALSE; 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 809b54502bSHong Zhang if (flg) { 819b54502bSHong Zhang PetscReal tol = 1.e-10; 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 839566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 849b54502bSHong Zhang } 85a1f19f5aSHong Zhang } 869566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 879566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 883d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 899b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 903d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 91f73b0415SBarry Smith PetscBool canuseordering; 92*03e5aca4SStefano Zampini 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 94*03e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 959566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 96f73b0415SBarry Smith if (canuseordering) { 979566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 989566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 999566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 1009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */ 1019d61c4f5SHong Zhang 10290d69ab7SBarry Smith flg = PETSC_FALSE; 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 1049b54502bSHong Zhang if (flg) { 1059b54502bSHong Zhang PetscReal tol = 1.e-10; 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 1089b54502bSHong Zhang } 1099b54502bSHong Zhang } 1102c7c0729SBarry Smith } 1119566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 1129566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 1133d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 11404545d6dSBarry Smith } else { 1159566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 116160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1179566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 118b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 11904545d6dSBarry Smith } 1209b54502bSHong Zhang } 1219566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12200e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1256baea169SHong Zhang } 1266baea169SHong Zhang 1279566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 1289566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12900e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 13000e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1316baea169SHong Zhang } 1329b54502bSHong Zhang } 13300c67f3bSHong Zhang 1349566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 13500c67f3bSHong Zhang if (!stype) { 136ea799195SBarry Smith MatSolverType solverpackage; 1379566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 1389566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 13900c67f3bSHong Zhang } 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1419b54502bSHong Zhang } 1429b54502bSHong Zhang 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Cholesky(PC pc) 144d71ae5a4SJacob Faibussowitsch { 1459b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1469b54502bSHong Zhang 1479b54502bSHong Zhang PetscFunctionBegin; 1489566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 1499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 1509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15269d2c0f9SBarry Smith } 15369d2c0f9SBarry Smith 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Cholesky(PC pc) 155d71ae5a4SJacob Faibussowitsch { 15669d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky *)pc->data; 15769d2c0f9SBarry Smith 15869d2c0f9SBarry Smith PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(PCReset_Cholesky(pc)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 1622e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1659b54502bSHong Zhang } 1669b54502bSHong Zhang 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) 168d71ae5a4SJacob Faibussowitsch { 1699b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1709b54502bSHong Zhang 1719b54502bSHong Zhang PetscFunctionBegin; 1723d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1739566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1742fa5cd67SKarl Rupp } else { 1759566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1762fa5cd67SKarl Rupp } 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1789b54502bSHong Zhang } 1799b54502bSHong Zhang 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) 181d71ae5a4SJacob Faibussowitsch { 1827b6e2003SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1837b6e2003SPierre Jolivet 1847b6e2003SPierre Jolivet PetscFunctionBegin; 1857b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1869566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1877b6e2003SPierre Jolivet } else { 1889566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1897b6e2003SPierre Jolivet } 1903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1917b6e2003SPierre Jolivet } 1927b6e2003SPierre Jolivet 193d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) 194d71ae5a4SJacob Faibussowitsch { 1953d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1963d72fe1eSJed Brown 1973d72fe1eSJed Brown PetscFunctionBegin; 1983d72fe1eSJed Brown if (dir->hdr.inplace) { 1999566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(pc->pmat, x, y)); 2003d72fe1eSJed Brown } else { 2019566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y)); 2023d72fe1eSJed Brown } 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2043d72fe1eSJed Brown } 2053d72fe1eSJed Brown 206d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) 207d71ae5a4SJacob Faibussowitsch { 2083d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2093d72fe1eSJed Brown 2103d72fe1eSJed Brown PetscFunctionBegin; 2113d72fe1eSJed Brown if (dir->hdr.inplace) { 2129566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(pc->pmat, x, y)); 2133d72fe1eSJed Brown } else { 2149566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y)); 2153d72fe1eSJed Brown } 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2173d72fe1eSJed Brown } 2183d72fe1eSJed Brown 219d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) 220d71ae5a4SJacob Faibussowitsch { 2219b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2229b54502bSHong Zhang 2239b54502bSHong Zhang PetscFunctionBegin; 2243d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2259566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 2262fa5cd67SKarl Rupp } else { 2279566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 2282fa5cd67SKarl Rupp } 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2309b54502bSHong Zhang } 2319b54502bSHong Zhang 2329b54502bSHong Zhang /*@ 2332401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2349b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2359b54502bSHong Zhang following factors. 2369b54502bSHong Zhang 237c3339decSBarry Smith Logically Collective 2389b54502bSHong Zhang 2399b54502bSHong Zhang Input Parameters: 2409b54502bSHong Zhang + pc - the preconditioner context 241f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 2429b54502bSHong Zhang 2439b54502bSHong Zhang Options Database Key: 244f1580f4eSBarry Smith . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 2459b54502bSHong Zhang 2469b54502bSHong Zhang Level: intermediate 2479b54502bSHong Zhang 248f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()` 2499b54502bSHong Zhang @*/ 250d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) 251d71ae5a4SJacob Faibussowitsch { 2529b54502bSHong Zhang PetscFunctionBegin; 2530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 254acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 255cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2579b54502bSHong Zhang } 2589b54502bSHong Zhang 2599b54502bSHong Zhang /*MC 26096fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2619b54502bSHong Zhang 2629b54502bSHong Zhang Options Database Keys: 263f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 264f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 265f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 26655ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2672401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 268145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2699b54502bSHong Zhang 2709b54502bSHong Zhang Level: beginner 2719b54502bSHong Zhang 27295452b02SPatrick Sanan Notes: 273f1580f4eSBarry Smith Not all options work for all matrix formats 274f1580f4eSBarry Smith 27595452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2769b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 277f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2789b54502bSHong Zhang 279db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 280db781477SPatrick Sanan `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 281db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 282f1580f4eSBarry Smith `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()` 2839b54502bSHong Zhang M*/ 2849b54502bSHong Zhang 285d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 286d71ae5a4SJacob Faibussowitsch { 2879b54502bSHong Zhang PC_Cholesky *dir; 2889b54502bSHong Zhang 2899b54502bSHong Zhang PetscFunctionBegin; 2904dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir)); 2913d1c1ea0SBarry Smith pc->data = (void *)dir; 2929566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY)); 2932fa5cd67SKarl Rupp 294075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 2952fa5cd67SKarl Rupp 2969b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 29769d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 2989b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 2997b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky; 3003d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 3013d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 3029b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3039b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3049b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 30592e08861SBarry Smith pc->ops->view = PCView_Factor; 3060a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3089b54502bSHong Zhang } 309