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 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 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) { 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; 9103e5aca4SStefano Zampini 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact)); 9303e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 949566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering)); 95f73b0415SBarry Smith if (canuseordering) { 969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->row)); 979566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 989566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col)); 999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col)); /* only use dir->row ordering in CholeskyFactor */ 1009d61c4f5SHong Zhang 10190d69ab7SBarry Smith flg = PETSC_FALSE; 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &flg, NULL)); 1039b54502bSHong Zhang if (flg) { 1049b54502bSHong Zhang PetscReal tol = 1.e-10; 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_factor_nonzeros_along_diagonal", &tol, NULL)); 1069566063dSJacob Faibussowitsch PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, tol, dir->row, dir->row)); 1079b54502bSHong Zhang } 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; 11304545d6dSBarry Smith } else { 1149566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 115160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 1169566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact)); 117b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR; 11804545d6dSBarry Smith } 1199b54502bSHong Zhang } 1209566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12100e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 12200e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1246baea169SHong Zhang } 1256baea169SHong Zhang 1269566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info)); 1279566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err)); 12800e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 12900e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1306baea169SHong Zhang } 1319b54502bSHong Zhang } 13200c67f3bSHong Zhang 1339566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype)); 13400c67f3bSHong Zhang if (!stype) { 135ea799195SBarry Smith MatSolverType solverpackage; 1369566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage)); 1379566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage)); 13800c67f3bSHong Zhang } 1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1409b54502bSHong Zhang } 1419b54502bSHong Zhang 142d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Cholesky(PC pc) 143d71ae5a4SJacob Faibussowitsch { 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)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15169d2c0f9SBarry Smith } 15269d2c0f9SBarry Smith 153d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Cholesky(PC pc) 154d71ae5a4SJacob Faibussowitsch { 15569d2c0f9SBarry Smith PC_Cholesky *dir = (PC_Cholesky *)pc->data; 15669d2c0f9SBarry Smith 15769d2c0f9SBarry Smith PetscFunctionBegin; 1589566063dSJacob Faibussowitsch PetscCall(PCReset_Cholesky(pc)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype)); 1612e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1649b54502bSHong Zhang } 1659b54502bSHong Zhang 166d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Cholesky(PC pc, Vec x, Vec y) 167d71ae5a4SJacob Faibussowitsch { 1689b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1699b54502bSHong Zhang 1709b54502bSHong Zhang PetscFunctionBegin; 1713d1c1ea0SBarry Smith if (dir->hdr.inplace) { 1729566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y)); 1732fa5cd67SKarl Rupp } else { 1749566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y)); 1752fa5cd67SKarl Rupp } 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1779b54502bSHong Zhang } 1789b54502bSHong Zhang 179d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Cholesky(PC pc, Mat X, Mat Y) 180d71ae5a4SJacob Faibussowitsch { 1817b6e2003SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1827b6e2003SPierre Jolivet 1837b6e2003SPierre Jolivet PetscFunctionBegin; 1847b6e2003SPierre Jolivet if (dir->hdr.inplace) { 1859566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y)); 1867b6e2003SPierre Jolivet } else { 1879566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y)); 1887b6e2003SPierre Jolivet } 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1907b6e2003SPierre Jolivet } 1917b6e2003SPierre Jolivet 192d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_Cholesky(PC pc, Vec x, Vec y) 193d71ae5a4SJacob Faibussowitsch { 1943d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 1953d72fe1eSJed Brown 1963d72fe1eSJed Brown PetscFunctionBegin; 1973d72fe1eSJed Brown if (dir->hdr.inplace) { 1989566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(pc->pmat, x, y)); 1993d72fe1eSJed Brown } else { 2009566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor *)dir)->fact, x, y)); 2013d72fe1eSJed Brown } 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2033d72fe1eSJed Brown } 2043d72fe1eSJed Brown 205d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_Cholesky(PC pc, Vec x, Vec y) 206d71ae5a4SJacob Faibussowitsch { 2073d72fe1eSJed Brown PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2083d72fe1eSJed Brown 2093d72fe1eSJed Brown PetscFunctionBegin; 2103d72fe1eSJed Brown if (dir->hdr.inplace) { 2119566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(pc->pmat, x, y)); 2123d72fe1eSJed Brown } else { 2139566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor *)dir)->fact, x, y)); 2143d72fe1eSJed Brown } 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2163d72fe1eSJed Brown } 2173d72fe1eSJed Brown 218d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Cholesky(PC pc, Vec x, Vec y) 219d71ae5a4SJacob Faibussowitsch { 2209b54502bSHong Zhang PC_Cholesky *dir = (PC_Cholesky *)pc->data; 2219b54502bSHong Zhang 2229b54502bSHong Zhang PetscFunctionBegin; 2233d1c1ea0SBarry Smith if (dir->hdr.inplace) { 2249566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y)); 2252fa5cd67SKarl Rupp } else { 2269566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y)); 2272fa5cd67SKarl Rupp } 2283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2299b54502bSHong Zhang } 2309b54502bSHong Zhang 231*4dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_Cholesky(PC pc, Mat X, Mat Y) 232*4dbf25a8SPierre Jolivet { 233*4dbf25a8SPierre Jolivet PC_Cholesky *dir = (PC_Cholesky *)pc->data; 234*4dbf25a8SPierre Jolivet 235*4dbf25a8SPierre Jolivet PetscFunctionBegin; 236*4dbf25a8SPierre Jolivet if (dir->hdr.inplace) { 237*4dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(pc->pmat, X, Y)); 238*4dbf25a8SPierre Jolivet } else { 239*4dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(((PC_Factor *)dir)->fact, X, Y)); 240*4dbf25a8SPierre Jolivet } 241*4dbf25a8SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS); 242*4dbf25a8SPierre Jolivet } 243*4dbf25a8SPierre Jolivet 2449b54502bSHong Zhang /*@ 2452401956bSBarry Smith PCFactorSetReuseOrdering - When similar matrices are factored, this 2469b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 2479b54502bSHong Zhang following factors. 2489b54502bSHong Zhang 249c3339decSBarry Smith Logically Collective 2509b54502bSHong Zhang 2519b54502bSHong Zhang Input Parameters: 2529b54502bSHong Zhang + pc - the preconditioner context 253f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 2549b54502bSHong Zhang 2559b54502bSHong Zhang Options Database Key: 256f1580f4eSBarry Smith . -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 2579b54502bSHong Zhang 2589b54502bSHong Zhang Level: intermediate 2599b54502bSHong Zhang 260562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetReuseFill()` 2619b54502bSHong Zhang @*/ 262d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseOrdering(PC pc, PetscBool flag) 263d71ae5a4SJacob Faibussowitsch { 2649b54502bSHong Zhang PetscFunctionBegin; 2650700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 266acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 267cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseOrdering_C", (PC, PetscBool), (pc, flag)); 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2699b54502bSHong Zhang } 2709b54502bSHong Zhang 2719b54502bSHong Zhang /*MC 27296fc60bcSBarry Smith PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner 2739b54502bSHong Zhang 2749b54502bSHong Zhang Options Database Keys: 275f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()` 276f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu 277f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 2780b4b7b1cSBarry Smith . -pc_factor_fill <fill> - Sets the explected fill amount 2792401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization 2800b4b7b1cSBarry 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 2810b4b7b1cSBarry Smith and thus be more effective 2829b54502bSHong Zhang 2839b54502bSHong Zhang Level: beginner 2849b54502bSHong Zhang 28595452b02SPatrick Sanan Notes: 2860b4b7b1cSBarry Smith The Cholesky factorization direct solver, `PCCHOLESKY` is only for symmetric positive-definite (SPD) matrices. For such 2870b4b7b1cSBarry Smith SPD matrices it is more efficient than using the LU factorization direct solver, `PCLU`. 2880b4b7b1cSBarry Smith 289f1580f4eSBarry Smith Not all options work for all matrix formats 290f1580f4eSBarry Smith 29195452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does 2929b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or 293f1580f4eSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method 2949b54502bSHong Zhang 295562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 296db781477SPatrick Sanan `PCILU`, `PCLU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`, 297ef959800SPierre Jolivet `PCFactorSetFill()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 298f1580f4eSBarry Smith `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()` 2999b54502bSHong Zhang M*/ 3009b54502bSHong Zhang 301d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Cholesky(PC pc) 302d71ae5a4SJacob Faibussowitsch { 3039b54502bSHong Zhang PC_Cholesky *dir; 3049b54502bSHong Zhang 3059b54502bSHong Zhang PetscFunctionBegin; 3064dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir)); 3073d1c1ea0SBarry Smith pc->data = (void *)dir; 3089566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_CHOLESKY)); 3092fa5cd67SKarl Rupp 310075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0; 3112fa5cd67SKarl Rupp 3129b54502bSHong Zhang pc->ops->destroy = PCDestroy_Cholesky; 31369d2c0f9SBarry Smith pc->ops->reset = PCReset_Cholesky; 3149b54502bSHong Zhang pc->ops->apply = PCApply_Cholesky; 3157b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Cholesky; 3163d72fe1eSJed Brown pc->ops->applysymmetricleft = PCApplySymmetricLeft_Cholesky; 3173d72fe1eSJed Brown pc->ops->applysymmetricright = PCApplySymmetricRight_Cholesky; 3189b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_Cholesky; 319*4dbf25a8SPierre Jolivet pc->ops->matapplytranspose = PCMatApplyTranspose_Cholesky; 3209b54502bSHong Zhang pc->ops->setup = PCSetUp_Cholesky; 3219b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_Cholesky; 32292e08861SBarry Smith pc->ops->view = PCView_Factor; 3230a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3259b54502bSHong Zhang } 326