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) { 40*fe143ba0SPablo Brubeck MatFactorType ftype; 41*fe143ba0SPablo Brubeck 42*fe143ba0SPablo Brubeck PetscCall(MatGetFactorType(pc->pmat, &ftype)); 43*fe143ba0SPablo 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)); 46*fe143ba0SPablo 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 } 58*fe143ba0SPablo 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; 654dfa11a4SJacob Faibussowitsch if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); } 669566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 67f73b0415SBarry Smith if (canuseordering) { 689566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 699566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 709bfd6278SHong Zhang /* check if dir->row == dir->col */ 714ac6704cSBarry Smith if (dir->row) { 729566063dSJacob Faibussowitsch PetscCall(ISEqual(dir->row, dir->col, &flg)); 7328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal"); 744ac6704cSBarry Smith } 759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */ 769bfd6278SHong Zhang 7790d69ab7SBarry Smith flg = PETSC_FALSE; 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 799b54502bSHong Zhang if (flg) { 809b54502bSHong Zhang PetscReal tol = 1.e-10; 819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 829566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 839b54502bSHong Zhang } 84a1f19f5aSHong Zhang } 859566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 869566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 873d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 889b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 893d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 90f73b0415SBarry Smith PetscBool canuseordering; 919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 929566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); 939566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 94f73b0415SBarry Smith if (canuseordering) { 959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 969566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 979566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */ 999d61c4f5SHong Zhang 10090d69ab7SBarry Smith flg = PETSC_FALSE; 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 1029b54502bSHong Zhang if (flg) { 1039b54502bSHong Zhang PetscReal tol = 1.e-10; 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 1069b54502bSHong Zhang } 1079b54502bSHong Zhang } 1082c7c0729SBarry Smith } 1099566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 1109566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 1113d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 11204545d6dSBarry Smith } else { 1139566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 114160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1159566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 116b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 11704545d6dSBarry Smith } 1189b54502bSHong Zhang } 1199566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12000e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1236baea169SHong Zhang } 1246baea169SHong Zhang 1259566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 1269566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12700e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 12800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1296baea169SHong Zhang } 1309b54502bSHong Zhang } 13100c67f3bSHong Zhang 1329566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 13300c67f3bSHong Zhang if (!stype) { 134ea799195SBarry Smith MatSolverType solverpackage; 1359566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 1369566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 13700c67f3bSHong Zhang } 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1399b54502bSHong Zhang } 1409b54502bSHong Zhang 141d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Cholesky(PC pc) 142d71ae5a4SJacob Faibussowitsch { 1439b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1449b54502bSHong Zhang 1459b54502bSHong Zhang PetscFunctionBegin; 1469566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 1479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 1489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15069d2c0f9SBarry Smith } 15169d2c0f9SBarry Smith 152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Cholesky(PC pc) 153d71ae5a4SJacob Faibussowitsch { 15469d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky *)pc->data; 15569d2c0f9SBarry Smith 15669d2c0f9SBarry Smith PetscFunctionBegin; 1579566063dSJacob Faibussowitsch PetscCall(PCReset_Cholesky(pc)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 1602e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1639b54502bSHong Zhang } 1649b54502bSHong Zhang 165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) 166d71ae5a4SJacob Faibussowitsch { 1679b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1689b54502bSHong Zhang 1699b54502bSHong Zhang PetscFunctionBegin; 1703d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1719566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1722fa5cd67SKarl Rupp } else { 1739566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1742fa5cd67SKarl Rupp } 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1769b54502bSHong Zhang } 1779b54502bSHong Zhang 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) 179d71ae5a4SJacob Faibussowitsch { 1807b6e2003SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1817b6e2003SPierre Jolivet 1827b6e2003SPierre Jolivet PetscFunctionBegin; 1837b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1849566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1857b6e2003SPierre Jolivet } else { 1869566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1877b6e2003SPierre Jolivet } 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1897b6e2003SPierre Jolivet } 1907b6e2003SPierre Jolivet 191d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) 192d71ae5a4SJacob Faibussowitsch { 1933d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1943d72fe1eSJed Brown 1953d72fe1eSJed Brown PetscFunctionBegin; 1963d72fe1eSJed Brown if (dir->hdr.inplace) { 1979566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(pc->pmat, x, y)); 1983d72fe1eSJed Brown } else { 1999566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y)); 2003d72fe1eSJed Brown } 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2023d72fe1eSJed Brown } 2033d72fe1eSJed Brown 204d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) 205d71ae5a4SJacob Faibussowitsch { 2063d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2073d72fe1eSJed Brown 2083d72fe1eSJed Brown PetscFunctionBegin; 2093d72fe1eSJed Brown if (dir->hdr.inplace) { 2109566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(pc->pmat, x, y)); 2113d72fe1eSJed Brown } else { 2129566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y)); 2133d72fe1eSJed Brown } 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2153d72fe1eSJed Brown } 2163d72fe1eSJed Brown 217d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) 218d71ae5a4SJacob Faibussowitsch { 2199b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2209b54502bSHong Zhang 2219b54502bSHong Zhang PetscFunctionBegin; 2223d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2239566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 2242fa5cd67SKarl Rupp } else { 2259566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 2262fa5cd67SKarl Rupp } 2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2289b54502bSHong Zhang } 2299b54502bSHong Zhang 2309b54502bSHong Zhang /*@ 2312401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2329b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2339b54502bSHong Zhang following factors. 2349b54502bSHong Zhang 235c3339decSBarry Smith Logically Collective 2369b54502bSHong Zhang 2379b54502bSHong Zhang Input Parameters: 2389b54502bSHong Zhang + pc - the preconditioner context 239f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 2409b54502bSHong Zhang 2419b54502bSHong Zhang Options Database Key: 242f1580f4eSBarry Smith . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 2439b54502bSHong Zhang 2449b54502bSHong Zhang Level: intermediate 2459b54502bSHong Zhang 246f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()` 2479b54502bSHong Zhang @*/ 248d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) 249d71ae5a4SJacob Faibussowitsch { 2509b54502bSHong Zhang PetscFunctionBegin; 2510700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 252acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 253cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag)); 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2559b54502bSHong Zhang } 2569b54502bSHong Zhang 2579b54502bSHong Zhang /*MC 25896fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2599b54502bSHong Zhang 2609b54502bSHong Zhang Options Database Keys: 261f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 262f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 263f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 26455ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2652401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 266145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2679b54502bSHong Zhang 2689b54502bSHong Zhang Level: beginner 2699b54502bSHong Zhang 27095452b02SPatrick Sanan Notes: 271f1580f4eSBarry Smith Not all options work for all matrix formats 272f1580f4eSBarry Smith 27395452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2749b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 275f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2769b54502bSHong Zhang 277db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 278db781477SPatrick Sanan `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 279db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 280f1580f4eSBarry Smith `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()` 2819b54502bSHong Zhang M*/ 2829b54502bSHong Zhang 283d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 284d71ae5a4SJacob Faibussowitsch { 2859b54502bSHong Zhang PC_Cholesky *dir; 2869b54502bSHong Zhang 2879b54502bSHong Zhang PetscFunctionBegin; 2884dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir)); 2893d1c1ea0SBarry Smith pc->data = (void *)dir; 2909566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY)); 2912fa5cd67SKarl Rupp 292075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 2932fa5cd67SKarl Rupp 2949b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 29569d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 2969b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 2977b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky; 2983d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 2993d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 3009b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 3019b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3029b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 30392e08861SBarry Smith pc->ops->view = PCView_Factor; 3040a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3069b54502bSHong Zhang } 307