19b54502bSHong Zhang /*
29b54502bSHong Zhang Defines a direct factorization preconditioner for any Mat implementation
3c14e9f00SDavid Andrs Note: this need not be considered a preconditioner since it supplies
49b54502bSHong Zhang a direct solver.
59b54502bSHong Zhang */
6c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
79b54502bSHong Zhang
89b54502bSHong Zhang typedef struct {
9075768bcSBarry Smith PC_Factor hdr;
109b54502bSHong Zhang IS row, col; /* index sets used for reordering */
119b54502bSHong Zhang } PC_Cholesky;
129b54502bSHong Zhang
PCSetFromOptions_Cholesky(PC pc,PetscOptionItems PetscOptionsObject)13ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Cholesky(PC pc, PetscOptionItems PetscOptionsObject)
14d71ae5a4SJacob Faibussowitsch {
159b54502bSHong Zhang PetscFunctionBegin;
16d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Cholesky options");
17dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
18d0609cedSBarry Smith PetscOptionsHeadEnd();
193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
209b54502bSHong Zhang }
219b54502bSHong Zhang
PCSetUp_Cholesky(PC pc)22d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Cholesky(PC pc)
23d71ae5a4SJacob Faibussowitsch {
24ace3abfcSBarry Smith PetscBool flg;
259b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data;
26ea799195SBarry Smith MatSolverType stype;
2700e125f8SBarry Smith MatFactorError err;
28f023e1d5SPierre Jolivet const char *prefix;
299b54502bSHong Zhang
309b54502bSHong Zhang PetscFunctionBegin;
31c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR;
323d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
339b54502bSHong Zhang
3426cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix));
3526cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
3626cc229bSBarry Smith
379566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
383d1c1ea0SBarry Smith if (dir->hdr.inplace) {
39fe143ba0SPablo Brubeck MatFactorType ftype;
40fe143ba0SPablo Brubeck
41fe143ba0SPablo Brubeck PetscCall(MatGetFactorType(pc->pmat, &ftype));
42fe143ba0SPablo Brubeck if (ftype == MAT_FACTOR_NONE) {
4348a46eb9SPierre Jolivet if (dir->row && dir->col && (dir->row != dir->col)) PetscCall(ISDestroy(&dir->row));
449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
45fe143ba0SPablo Brubeck /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
469566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
479566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
489b54502bSHong Zhang if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
509b54502bSHong Zhang }
519566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(pc->pmat, dir->row, &((PC_Factor *)dir)->info));
529566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err));
5300e125f8SBarry Smith if (err) { /* Factor() fails */
5400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
566baea169SHong Zhang }
57fe143ba0SPablo Brubeck }
58075768bcSBarry Smith ((PC_Factor *)dir)->fact = pc->pmat;
599b54502bSHong Zhang } else {
609b54502bSHong Zhang MatInfo info;
6100e125f8SBarry Smith
629b54502bSHong Zhang if (!pc->setupcalled) {
63f73b0415SBarry Smith PetscBool canuseordering;
6403e5aca4SStefano Zampini
6503e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc));
669566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
67f73b0415SBarry Smith if (canuseordering) {
68c294e005SBarry Smith PetscBool external;
69c294e005SBarry Smith
709566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
71c294e005SBarry Smith PetscCall(PetscStrcmp(((PC_Factor *)dir)->ordering, MATORDERINGEXTERNAL, &external));
72c294e005SBarry Smith if (!external) {
739566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
749bfd6278SHong Zhang /* check if dir->row == dir->col */
754ac6704cSBarry Smith if (dir->row) {
769566063dSJacob Faibussowitsch PetscCall(ISEqual(dir->row, dir->col, &flg));
7728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row and column permutations must be equal");
784ac6704cSBarry Smith }
799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only pass one ordering into CholeskyFactor */
809bfd6278SHong Zhang
8190d69ab7SBarry Smith flg = PETSC_FALSE;
82c294e005SBarry Smith PetscCall(PetscOptionsHasName(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg));
839b54502bSHong Zhang if (flg) {
849b54502bSHong Zhang PetscReal tol = 1.e-10;
859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
869566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
879b54502bSHong Zhang }
88a1f19f5aSHong Zhang }
89c294e005SBarry Smith }
909566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
919566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
923d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed;
939b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) {
94*da34b7cdSBarry Smith PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
95*da34b7cdSBarry Smith PetscCall(PCFactorSetUpMatSolverType(pc));
963d1c1ea0SBarry Smith if (!dir->hdr.reuseordering) {
97f73b0415SBarry Smith PetscBool canuseordering;
9803e5aca4SStefano Zampini
999566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
100f73b0415SBarry Smith if (canuseordering) {
101c294e005SBarry Smith PetscBool external;
102c294e005SBarry Smith
1039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row));
1049566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
105c294e005SBarry Smith PetscCall(PetscStrcmp(((PC_Factor *)dir)->ordering, MATORDERINGEXTERNAL, &external));
106c294e005SBarry Smith if (!external) {
1079566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
1089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */
1099d61c4f5SHong Zhang
11090d69ab7SBarry Smith flg = PETSC_FALSE;
111c294e005SBarry Smith PetscCall(PetscOptionsHasName(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg));
1129b54502bSHong Zhang if (flg) {
1139b54502bSHong Zhang PetscReal tol = 1.e-10;
1149566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL));
1159566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row));
1169b54502bSHong Zhang }
1179b54502bSHong Zhang }
1182c7c0729SBarry Smith }
119c294e005SBarry Smith }
1209566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, &((PC_Factor *)dir)->info));
1219566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1223d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed;
12304545d6dSBarry Smith } else {
1249566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
125160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1269566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
127b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR;
12804545d6dSBarry Smith }
1299b54502bSHong Zhang }
1309566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
13100e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */
13200e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1346baea169SHong Zhang }
1356baea169SHong Zhang
1369566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1379566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
13800e125f8SBarry Smith if (err) { /* FactorNumeric() fails */
13900e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
1406baea169SHong Zhang }
1419b54502bSHong Zhang }
14200c67f3bSHong Zhang
1439566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype));
14400c67f3bSHong Zhang if (!stype) {
145ea799195SBarry Smith MatSolverType solverpackage;
1469566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1479566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
14800c67f3bSHong Zhang }
1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1509b54502bSHong Zhang }
1519b54502bSHong Zhang
PCReset_Cholesky(PC pc)152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Cholesky(PC pc)
153d71ae5a4SJacob Faibussowitsch {
1549b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1559b54502bSHong Zhang
1569b54502bSHong Zhang PetscFunctionBegin;
1579566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row));
1599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16169d2c0f9SBarry Smith }
16269d2c0f9SBarry Smith
PCDestroy_Cholesky(PC pc)163d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Cholesky(PC pc)
164d71ae5a4SJacob Faibussowitsch {
16569d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky *)pc->data;
16669d2c0f9SBarry Smith
16769d2c0f9SBarry Smith PetscFunctionBegin;
1689566063dSJacob Faibussowitsch PetscCall(PCReset_Cholesky(pc));
1699566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1709566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1712e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc));
1729566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1749b54502bSHong Zhang }
1759b54502bSHong Zhang
PCApply_Cholesky(PC pc,Vec x,Vec y)176d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y)
177d71ae5a4SJacob Faibussowitsch {
1789b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1799b54502bSHong Zhang
1809b54502bSHong Zhang PetscFunctionBegin;
1813d1c1ea0SBarry Smith if (dir->hdr.inplace) {
1829566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y));
1832fa5cd67SKarl Rupp } else {
1849566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1852fa5cd67SKarl Rupp }
1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1879b54502bSHong Zhang }
1889b54502bSHong Zhang
PCMatApply_Cholesky(PC pc,Mat X,Mat Y)189d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y)
190d71ae5a4SJacob Faibussowitsch {
1917b6e2003SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data;
1927b6e2003SPierre Jolivet
1937b6e2003SPierre Jolivet PetscFunctionBegin;
1947b6e2003SPierre Jolivet if (dir->hdr.inplace) {
1959566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y));
1967b6e2003SPierre Jolivet } else {
1979566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1987b6e2003SPierre Jolivet }
1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2007b6e2003SPierre Jolivet }
2017b6e2003SPierre Jolivet
PCApplySymmetricLeft_Cholesky(PC pc,Vec x,Vec y)202d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y)
203d71ae5a4SJacob Faibussowitsch {
2043d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2053d72fe1eSJed Brown
2063d72fe1eSJed Brown PetscFunctionBegin;
2073d72fe1eSJed Brown if (dir->hdr.inplace) {
2089566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(pc->pmat, x, y));
2093d72fe1eSJed Brown } else {
2109566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y));
2113d72fe1eSJed Brown }
2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2133d72fe1eSJed Brown }
2143d72fe1eSJed Brown
PCApplySymmetricRight_Cholesky(PC pc,Vec x,Vec y)215d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y)
216d71ae5a4SJacob Faibussowitsch {
2173d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2183d72fe1eSJed Brown
2193d72fe1eSJed Brown PetscFunctionBegin;
2203d72fe1eSJed Brown if (dir->hdr.inplace) {
2219566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(pc->pmat, x, y));
2223d72fe1eSJed Brown } else {
2239566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y));
2243d72fe1eSJed Brown }
2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2263d72fe1eSJed Brown }
2273d72fe1eSJed Brown
PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)228d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y)
229d71ae5a4SJacob Faibussowitsch {
2309b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2319b54502bSHong Zhang
2329b54502bSHong Zhang PetscFunctionBegin;
2333d1c1ea0SBarry Smith if (dir->hdr.inplace) {
2349566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y));
2352fa5cd67SKarl Rupp } else {
2369566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
2372fa5cd67SKarl Rupp }
2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2399b54502bSHong Zhang }
2409b54502bSHong Zhang
PCMatApplyTranspose_Cholesky(PC pc,Mat X,Mat Y)2414dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_Cholesky(PC pc, Mat X, Mat Y)
2424dbf25a8SPierre Jolivet {
2434dbf25a8SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data;
2444dbf25a8SPierre Jolivet
2454dbf25a8SPierre Jolivet PetscFunctionBegin;
2464dbf25a8SPierre Jolivet if (dir->hdr.inplace) {
2474dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(pc->pmat, X, Y));
2484dbf25a8SPierre Jolivet } else {
2494dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(((PC_Factor *)dir)->fact, X, Y));
2504dbf25a8SPierre Jolivet }
2514dbf25a8SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS);
2524dbf25a8SPierre Jolivet }
2534dbf25a8SPierre Jolivet
2549b54502bSHong Zhang /*@
2552401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this
2569b54502bSHong Zhang causes the ordering computed in the first factor to be used for all
2579b54502bSHong Zhang following factors.
2589b54502bSHong Zhang
259c3339decSBarry Smith Logically Collective
2609b54502bSHong Zhang
2619b54502bSHong Zhang Input Parameters:
2629b54502bSHong Zhang + pc - the preconditioner context
263f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE`
2649b54502bSHong Zhang
2659b54502bSHong Zhang Options Database Key:
266f1580f4eSBarry Smith . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
2679b54502bSHong Zhang
2689b54502bSHong Zhang Level: intermediate
2699b54502bSHong Zhang
270562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()`
2719b54502bSHong Zhang @*/
PCFactorSetReuseOrdering(PC pc,PetscBool flag)272d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag)
273d71ae5a4SJacob Faibussowitsch {
2749b54502bSHong Zhang PetscFunctionBegin;
2750700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
276acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2);
277cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag));
2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2799b54502bSHong Zhang }
2809b54502bSHong Zhang
2819b54502bSHong Zhang /*MC
28296fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
2839b54502bSHong Zhang
2849b54502bSHong Zhang Options Database Keys:
285f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
286f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
287f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
2880b4b7b1cSBarry Smith . -pc_factor_fill <fill> - Sets the explected fill amount
2892401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization
2900b4b7b1cSBarry Smith - -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine used to determine the order the rows are used in the factorization to reduce fill
2910b4b7b1cSBarry Smith and thus be more effective
2929b54502bSHong Zhang
2939b54502bSHong Zhang Level: beginner
2949b54502bSHong Zhang
29595452b02SPatrick Sanan Notes:
2960b4b7b1cSBarry Smith The Cholesky factorization direct solver, `PCCHOLESKY` is only for symmetric positive-definite (SPD) matrices. For such
2970b4b7b1cSBarry Smith SPD matrices it is more efficient than using the LU factorization direct solver, `PCLU`.
2980b4b7b1cSBarry Smith
299f1580f4eSBarry Smith Not all options work for all matrix formats
300f1580f4eSBarry Smith
30195452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does
3029b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or
303f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
3049b54502bSHong Zhang
305562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
306db781477SPatrick Sanan `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
307ef959800SPierre Jolivet `PCFactorSetFill()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
308f1580f4eSBarry Smith `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`
3099b54502bSHong Zhang M*/
3109b54502bSHong Zhang
PCCreate_Cholesky(PC pc)311d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc)
312d71ae5a4SJacob Faibussowitsch {
3139b54502bSHong Zhang PC_Cholesky *dir;
3149b54502bSHong Zhang
3159b54502bSHong Zhang PetscFunctionBegin;
3164dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir));
3173d1c1ea0SBarry Smith pc->data = (void *)dir;
3189566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY));
3192fa5cd67SKarl Rupp
320075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0;
3212fa5cd67SKarl Rupp
3229b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky;
32369d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky;
3249b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky;
3257b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky;
3263d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky;
3273d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky;
3289b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky;
3294dbf25a8SPierre Jolivet pc->ops->matapplytranspose = PCMatApplyTranspose_Cholesky;
3309b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky;
3319b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky;
33292e08861SBarry Smith pc->ops->view = PCView_Factor;
3330a545947SLisandro Dalcin pc->ops->applyrichardson = NULL;
3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3359b54502bSHong Zhang }
336