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 149371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems *PetscOptionsObject) { 159b54502bSHong Zhang PetscFunctionBegin; 16d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options"); 17dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject)); 18d0609cedSBarry Smith PetscOptionsHeadEnd(); 199b54502bSHong Zhang PetscFunctionReturn(0); 209b54502bSHong Zhang } 219b54502bSHong Zhang 229371c9d4SSatish Balay static PetscErrorCode PCSetUp_Cholesky(PC pc) { 23ace3abfcSBarry Smith PetscBool flg; 249b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 25ea799195SBarry Smith MatSolverType stype; 2600e125f8SBarry Smith MatFactorError err; 27f023e1d5SPierre Jolivet const char *prefix; 289b54502bSHong Zhang 299b54502bSHong Zhang PetscFunctionBegin; 30c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 313d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill; 329b54502bSHong Zhang 3326cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 3426cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 3526cc229bSBarry Smith 369566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 373d1c1ea0SBarry Smith if (dir->hdr.inplace) { 3848a46eb9SPierre Jolivet if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row)); 399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 402c7c0729SBarry Smith /* should only get reordering if the factor matrix uses it but cannot determine because MatGetFactor() not called */ 419566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 429566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 439b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */ 449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 459b54502bSHong Zhang } 469566063dSJacob Faibussowitsch if (dir->row) PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row)); 479566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 489566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err)); 4900e125f8SBarry Smith if (err) { /* Factor() fails */ 5000e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 516baea169SHong Zhang PetscFunctionReturn(0); 526baea169SHong Zhang } 532fa5cd67SKarl Rupp 54075768bcSBarry Smith ((PC_Factor *)dir)->fact = pc->pmat; 559b54502bSHong Zhang } else { 569b54502bSHong Zhang MatInfo info; 5700e125f8SBarry Smith 589b54502bSHong Zhang if (!pc->setupcalled) { 59f73b0415SBarry Smith PetscBool canuseordering; 602c7c0729SBarry Smith if (!((PC_Factor *)dir)->fact) { 619566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); 629566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)dir)->fact)); 632c7c0729SBarry Smith } 649566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 65f73b0415SBarry Smith if (canuseordering) { 669566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 679566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 689bfd6278SHong Zhang /* check if dir->row == dir->col */ 694ac6704cSBarry Smith if (dir->row) { 709566063dSJacob Faibussowitsch PetscCall(ISEqual(dir->row, dir->col, &flg)); 7128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal"); 724ac6704cSBarry Smith } 739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */ 749bfd6278SHong Zhang 7590d69ab7SBarry Smith flg = PETSC_FALSE; 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 779b54502bSHong Zhang if (flg) { 789b54502bSHong Zhang PetscReal tol = 1.e-10; 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 809566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 819b54502bSHong Zhang } 829566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row)); 83a1f19f5aSHong Zhang } 849566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 859566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 863d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 879b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 883d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) { 89f73b0415SBarry Smith PetscBool canuseordering; 909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 919566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_CHOLESKY, &((PC_Factor *)dir)->fact)); 929566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((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 } 1079566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row)); 1089b54502bSHong Zhang } 1092c7c0729SBarry Smith } 1109566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info)); 1119566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info)); 1123d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed; 1139566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)dir)->fact)); 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; 1246baea169SHong Zhang PetscFunctionReturn(0); 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 } 1409b54502bSHong Zhang PetscFunctionReturn(0); 1419b54502bSHong Zhang } 1429b54502bSHong Zhang 1439371c9d4SSatish Balay static PetscErrorCode PCReset_Cholesky(PC pc) { 1449b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1459b54502bSHong Zhang 1469b54502bSHong Zhang PetscFunctionBegin; 1479566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 1489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 1499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); 15069d2c0f9SBarry Smith PetscFunctionReturn(0); 15169d2c0f9SBarry Smith } 15269d2c0f9SBarry Smith 1539371c9d4SSatish Balay static PetscErrorCode PCDestroy_Cholesky(PC pc) { 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)); 1629b54502bSHong Zhang PetscFunctionReturn(0); 1639b54502bSHong Zhang } 1649b54502bSHong Zhang 1659371c9d4SSatish Balay static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) { 1669b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1679b54502bSHong Zhang 1689b54502bSHong Zhang PetscFunctionBegin; 1693d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1709566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1712fa5cd67SKarl Rupp } else { 1729566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1732fa5cd67SKarl Rupp } 1749b54502bSHong Zhang PetscFunctionReturn(0); 1759b54502bSHong Zhang } 1769b54502bSHong Zhang 1779371c9d4SSatish Balay static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) { 1787b6e2003SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1797b6e2003SPierre Jolivet 1807b6e2003SPierre Jolivet PetscFunctionBegin; 1817b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1829566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1837b6e2003SPierre Jolivet } else { 1849566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1857b6e2003SPierre Jolivet } 1867b6e2003SPierre Jolivet PetscFunctionReturn(0); 1877b6e2003SPierre Jolivet } 1887b6e2003SPierre Jolivet 1899371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) { 1903d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1913d72fe1eSJed Brown 1923d72fe1eSJed Brown PetscFunctionBegin; 1933d72fe1eSJed Brown if (dir->hdr.inplace) { 1949566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(pc->pmat, x, y)); 1953d72fe1eSJed Brown } else { 1969566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y)); 1973d72fe1eSJed Brown } 1983d72fe1eSJed Brown PetscFunctionReturn(0); 1993d72fe1eSJed Brown } 2003d72fe1eSJed Brown 2019371c9d4SSatish Balay static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) { 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 2139371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) { 2149b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2159b54502bSHong Zhang 2169b54502bSHong Zhang PetscFunctionBegin; 2173d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2189566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 2192fa5cd67SKarl Rupp } else { 2209566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 2212fa5cd67SKarl Rupp } 2229b54502bSHong Zhang PetscFunctionReturn(0); 2239b54502bSHong Zhang } 2249b54502bSHong Zhang 2259b54502bSHong Zhang /*@ 2262401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2279b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2289b54502bSHong Zhang following factors. 2299b54502bSHong Zhang 230*f1580f4eSBarry Smith Logically Collective on pc 2319b54502bSHong Zhang 2329b54502bSHong Zhang Input Parameters: 2339b54502bSHong Zhang + pc - the preconditioner context 234*f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 2359b54502bSHong Zhang 2369b54502bSHong Zhang Options Database Key: 237*f1580f4eSBarry Smith . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 2389b54502bSHong Zhang 2399b54502bSHong Zhang Level: intermediate 2409b54502bSHong Zhang 241*f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()` 2429b54502bSHong Zhang @*/ 2439371c9d4SSatish Balay PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) { 2449b54502bSHong Zhang PetscFunctionBegin; 2450700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 246acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 247cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag)); 2489b54502bSHong Zhang PetscFunctionReturn(0); 2499b54502bSHong Zhang } 2509b54502bSHong Zhang 2519b54502bSHong Zhang /*MC 25296fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2539b54502bSHong Zhang 2549b54502bSHong Zhang Options Database Keys: 255*f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 256*f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 257*f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 25855ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 2592401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 260145b9266SHong Zhang - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2619b54502bSHong Zhang 2629b54502bSHong Zhang Level: beginner 2639b54502bSHong Zhang 26495452b02SPatrick Sanan Notes: 265*f1580f4eSBarry Smith Not all options work for all matrix formats 266*f1580f4eSBarry Smith 26795452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2689b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 269*f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2709b54502bSHong Zhang 271db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 272db781477SPatrick Sanan `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 273db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 274*f1580f4eSBarry Smith `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()` 2759b54502bSHong Zhang M*/ 2769b54502bSHong Zhang 2779371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) { 2789b54502bSHong Zhang PC_Cholesky *dir; 2799b54502bSHong Zhang 2809b54502bSHong Zhang PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &dir)); 2823d1c1ea0SBarry Smith pc->data = (void *)dir; 2839566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY)); 2842fa5cd67SKarl Rupp 285075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 2862fa5cd67SKarl Rupp 2879b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 28869d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 2899b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 2907b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky; 2913d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 2923d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 2939b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 2949b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 2959b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 29692e08861SBarry Smith pc->ops->view = PCView_Factor; 2970a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2989b54502bSHong Zhang PetscFunctionReturn(0); 2999b54502bSHong Zhang } 300