xref: /petsc/src/ksp/pc/impls/factor/cholesky/cholesky.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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();
209b54502bSHong Zhang   PetscFunctionReturn(0);
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) {
4048a46eb9SPierre Jolivet     if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row));
419566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&dir->col));
422c7c0729SBarry Smith     /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */
439566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
449566063dSJacob Faibussowitsch     PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
459b54502bSHong Zhang     if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
469566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
479b54502bSHong Zhang     }
489566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info));
499566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(pc->pmat, &err));
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) {
60f73b0415SBarry Smith       PetscBool canuseordering;
614dfa11a4SJacob Faibussowitsch       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); }
629566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
63f73b0415SBarry Smith       if (canuseordering) {
649566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
659566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
669bfd6278SHong Zhang         /* check if dir->row == dir->col */
674ac6704cSBarry Smith         if (dir->row) {
689566063dSJacob Faibussowitsch           PetscCall(ISEqual(dir->row, dir->col, &flg));
6928b400f6SJacob Faibussowitsch           PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal");
704ac6704cSBarry Smith         }
719566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
729bfd6278SHong Zhang 
7390d69ab7SBarry Smith         flg = PETSC_FALSE;
749566063dSJacob Faibussowitsch         PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
759b54502bSHong Zhang         if (flg) {
769b54502bSHong Zhang           PetscReal tol = 1.e-10;
779566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
789566063dSJacob Faibussowitsch           PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
799b54502bSHong Zhang         }
80a1f19f5aSHong Zhang       }
819566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
829566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
833d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
849b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
853d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
86f73b0415SBarry Smith         PetscBool canuseordering;
879566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
889566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact));
899566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
90f73b0415SBarry Smith         if (canuseordering) {
919566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->row));
929566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
939566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
949566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
959d61c4f5SHong Zhang 
9690d69ab7SBarry Smith           flg = PETSC_FALSE;
979566063dSJacob Faibussowitsch           PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL));
989b54502bSHong Zhang           if (flg) {
999b54502bSHong Zhang             PetscReal tol = 1.e-10;
1009566063dSJacob Faibussowitsch             PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
1019566063dSJacob Faibussowitsch             PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
1029b54502bSHong Zhang           }
1039b54502bSHong Zhang         }
1042c7c0729SBarry Smith       }
1059566063dSJacob Faibussowitsch       PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
1069566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1073d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
10804545d6dSBarry Smith     } else {
1099566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
110160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1119566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
112b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11304545d6dSBarry Smith       }
1149b54502bSHong Zhang     }
1159566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
11600e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
11700e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1186baea169SHong Zhang       PetscFunctionReturn(0);
1196baea169SHong Zhang     }
1206baea169SHong Zhang 
1219566063dSJacob Faibussowitsch     PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1229566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
12300e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1256baea169SHong Zhang     }
1269b54502bSHong Zhang   }
12700c67f3bSHong Zhang 
1289566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
12900c67f3bSHong Zhang   if (!stype) {
130ea799195SBarry Smith     MatSolverType solverpackage;
1319566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1329566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
13300c67f3bSHong Zhang   }
1349b54502bSHong Zhang   PetscFunctionReturn(0);
1359b54502bSHong Zhang }
1369b54502bSHong Zhang 
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Cholesky(PC pc)
138d71ae5a4SJacob Faibussowitsch {
1399b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1409b54502bSHong Zhang 
1419b54502bSHong Zhang   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->row));
1449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
14569d2c0f9SBarry Smith   PetscFunctionReturn(0);
14669d2c0f9SBarry Smith }
14769d2c0f9SBarry Smith 
148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Cholesky(PC pc)
149d71ae5a4SJacob Faibussowitsch {
15069d2c0f9SBarry Smith   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
15169d2c0f9SBarry Smith 
15269d2c0f9SBarry Smith   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(PCReset_Cholesky(pc));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1562e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1589b54502bSHong Zhang   PetscFunctionReturn(0);
1599b54502bSHong Zhang }
1609b54502bSHong Zhang 
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y)
162d71ae5a4SJacob Faibussowitsch {
1639b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1649b54502bSHong Zhang 
1659b54502bSHong Zhang   PetscFunctionBegin;
1663d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1679566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat, x, y));
1682fa5cd67SKarl Rupp   } else {
1699566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1702fa5cd67SKarl Rupp   }
1719b54502bSHong Zhang   PetscFunctionReturn(0);
1729b54502bSHong Zhang }
1739b54502bSHong Zhang 
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y)
175d71ae5a4SJacob Faibussowitsch {
1767b6e2003SPierre Jolivet   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1777b6e2003SPierre Jolivet 
1787b6e2003SPierre Jolivet   PetscFunctionBegin;
1797b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1809566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat, X, Y));
1817b6e2003SPierre Jolivet   } else {
1829566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1837b6e2003SPierre Jolivet   }
1847b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1857b6e2003SPierre Jolivet }
1867b6e2003SPierre Jolivet 
187d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y)
188d71ae5a4SJacob Faibussowitsch {
1893d72fe1eSJed Brown   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1903d72fe1eSJed Brown 
1913d72fe1eSJed Brown   PetscFunctionBegin;
1923d72fe1eSJed Brown   if (dir->hdr.inplace) {
1939566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(pc->pmat, x, y));
1943d72fe1eSJed Brown   } else {
1959566063dSJacob Faibussowitsch     PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y));
1963d72fe1eSJed Brown   }
1973d72fe1eSJed Brown   PetscFunctionReturn(0);
1983d72fe1eSJed Brown }
1993d72fe1eSJed Brown 
200d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y)
201d71ae5a4SJacob Faibussowitsch {
2023d72fe1eSJed Brown   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2033d72fe1eSJed Brown 
2043d72fe1eSJed Brown   PetscFunctionBegin;
2053d72fe1eSJed Brown   if (dir->hdr.inplace) {
2069566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(pc->pmat, x, y));
2073d72fe1eSJed Brown   } else {
2089566063dSJacob Faibussowitsch     PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y));
2093d72fe1eSJed Brown   }
2103d72fe1eSJed Brown   PetscFunctionReturn(0);
2113d72fe1eSJed Brown }
2123d72fe1eSJed Brown 
213d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y)
214d71ae5a4SJacob Faibussowitsch {
2159b54502bSHong Zhang   PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2169b54502bSHong Zhang 
2179b54502bSHong Zhang   PetscFunctionBegin;
2183d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
2199566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat, x, y));
2202fa5cd67SKarl Rupp   } else {
2219566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
2222fa5cd67SKarl Rupp   }
2239b54502bSHong Zhang   PetscFunctionReturn(0);
2249b54502bSHong Zhang }
2259b54502bSHong Zhang 
2269b54502bSHong Zhang /*@
2272401956bSBarry Smith    PCFactorSetReuseOrdering - When similar matrices are factored, this
2289b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
2299b54502bSHong Zhang    following factors.
2309b54502bSHong Zhang 
231*c3339decSBarry Smith    Logically Collective
2329b54502bSHong Zhang 
2339b54502bSHong Zhang    Input Parameters:
2349b54502bSHong Zhang +  pc - the preconditioner context
235f1580f4eSBarry Smith -  flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
2369b54502bSHong Zhang 
2379b54502bSHong Zhang    Options Database Key:
238f1580f4eSBarry Smith .  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
2399b54502bSHong Zhang 
2409b54502bSHong Zhang    Level: intermediate
2419b54502bSHong Zhang 
242f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()`
2439b54502bSHong Zhang @*/
244d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag)
245d71ae5a4SJacob Faibussowitsch {
2469b54502bSHong Zhang   PetscFunctionBegin;
2470700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
248acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, flag, 2);
249cac4c232SBarry Smith   PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag));
2509b54502bSHong Zhang   PetscFunctionReturn(0);
2519b54502bSHong Zhang }
2529b54502bSHong Zhang 
2539b54502bSHong Zhang /*MC
25496fc60bcSBarry Smith    PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2559b54502bSHong Zhang 
2569b54502bSHong Zhang    Options Database Keys:
257f1580f4eSBarry Smith +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
258f1580f4eSBarry Smith .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
259f1580f4eSBarry Smith .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
26055ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2612401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
262145b9266SHong Zhang -  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2639b54502bSHong Zhang 
2649b54502bSHong Zhang    Level: beginner
2659b54502bSHong Zhang 
26695452b02SPatrick Sanan    Notes:
267f1580f4eSBarry Smith    Not all options work for all matrix formats
268f1580f4eSBarry Smith 
26995452b02SPatrick Sanan    Usually this will compute an "exact" solution in one iteration and does
2709b54502bSHong Zhang    not need a Krylov method (i.e. you can use -ksp_type preonly, or
271f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
2729b54502bSHong Zhang 
273db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
274db781477SPatrick Sanan           `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
275db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
276f1580f4eSBarry Smith           `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`
2779b54502bSHong Zhang M*/
2789b54502bSHong Zhang 
279d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
280d71ae5a4SJacob Faibussowitsch {
2819b54502bSHong Zhang   PC_Cholesky *dir;
2829b54502bSHong Zhang 
2839b54502bSHong Zhang   PetscFunctionBegin;
2844dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dir));
2853d1c1ea0SBarry Smith   pc->data = (void *)dir;
2869566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY));
2872fa5cd67SKarl Rupp 
288075768bcSBarry Smith   ((PC_Factor *)dir)->info.fill = 5.0;
2892fa5cd67SKarl Rupp 
2909b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_Cholesky;
29169d2c0f9SBarry Smith   pc->ops->reset               = PCReset_Cholesky;
2929b54502bSHong Zhang   pc->ops->apply               = PCApply_Cholesky;
2937b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Cholesky;
2943d72fe1eSJed Brown   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_Cholesky;
2953d72fe1eSJed Brown   pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
2969b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_Cholesky;
2979b54502bSHong Zhang   pc->ops->setup               = PCSetUp_Cholesky;
2989b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_Cholesky;
29992e08861SBarry Smith   pc->ops->view                = PCView_Factor;
3000a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
3019b54502bSHong Zhang   PetscFunctionReturn(0);
3029b54502bSHong Zhang }
303