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 */
6ee45ca4aSHong Zhang
7c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/
89b54502bSHong Zhang
PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)966976f2fSJacob Faibussowitsch static PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc, PetscReal z)
10d71ae5a4SJacob Faibussowitsch {
119b54502bSHong Zhang PC_LU *lu = (PC_LU *)pc->data;
129b54502bSHong Zhang
139b54502bSHong Zhang PetscFunctionBegin;
149b54502bSHong Zhang lu->nonzerosalongdiagonal = PETSC_TRUE;
1513bcc0bdSJacob Faibussowitsch if (z == (PetscReal)PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
162fa5cd67SKarl Rupp else lu->nonzerosalongdiagonaltol = z;
173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
189b54502bSHong Zhang }
199b54502bSHong Zhang
PCSetFromOptions_LU(PC pc,PetscOptionItems PetscOptionsObject)20ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_LU(PC pc, PetscOptionItems PetscOptionsObject)
21d71ae5a4SJacob Faibussowitsch {
229b54502bSHong Zhang PC_LU *lu = (PC_LU *)pc->data;
23ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE;
249b54502bSHong Zhang PetscReal tol;
259b54502bSHong Zhang
269b54502bSHong Zhang PetscFunctionBegin;
27d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "LU options");
28dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
295c9eb25fSBarry Smith
309566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
319b54502bSHong Zhang if (flg) {
329b54502bSHong Zhang tol = PETSC_DECIDE;
339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", lu->nonzerosalongdiagonaltol, &tol, NULL));
349566063dSJacob Faibussowitsch PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
359b54502bSHong Zhang }
36d0609cedSBarry Smith PetscOptionsHeadEnd();
373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
389b54502bSHong Zhang }
399b54502bSHong Zhang
PCSetUp_LU(PC pc)40d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_LU(PC pc)
41d71ae5a4SJacob Faibussowitsch {
429b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data;
43ea799195SBarry Smith MatSolverType stype;
4400e125f8SBarry Smith MatFactorError err;
45f023e1d5SPierre Jolivet const char *prefix;
463d1c1ea0SBarry Smith
479b54502bSHong Zhang PetscFunctionBegin;
48c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR;
493d1c1ea0SBarry Smith if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
509b54502bSHong Zhang
5126cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix));
5226cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
5326cc229bSBarry Smith
549566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
553d1c1ea0SBarry Smith if (dir->hdr.inplace) {
569d76b4d0SMatthew G. Knepley MatFactorType ftype;
579d76b4d0SMatthew G. Knepley
589566063dSJacob Faibussowitsch PetscCall(MatGetFactorType(pc->pmat, &ftype));
599d76b4d0SMatthew G. Knepley if (ftype == MAT_FACTOR_NONE) {
609566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
622c7c0729SBarry Smith /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
639566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
649566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
659566063dSJacob Faibussowitsch PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
669566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, &err));
6700e125f8SBarry Smith if (err) { /* Factor() fails */
6800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
706baea169SHong Zhang }
719d76b4d0SMatthew G. Knepley }
72075768bcSBarry Smith ((PC_Factor *)dir)->fact = pc->pmat;
739b54502bSHong Zhang } else {
749b54502bSHong Zhang MatInfo info;
7500e125f8SBarry Smith
769b54502bSHong Zhang if (!pc->setupcalled) {
77f73b0415SBarry Smith PetscBool canuseordering;
7803e5aca4SStefano Zampini
7903e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc));
809566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
81f73b0415SBarry Smith if (canuseordering) {
82c294e005SBarry Smith PetscBool external;
83c294e005SBarry Smith
849566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
85c294e005SBarry Smith PetscCall(PetscStrcmp(((PC_Factor *)dir)->ordering, MATORDERINGEXTERNAL, &external));
86c294e005SBarry Smith if (!external) {
879566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
881baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
8903c60df9SBarry Smith }
90c294e005SBarry Smith }
919566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
929566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
933d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed;
949b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) {
95f73b0415SBarry Smith PetscBool canuseordering;
9603e5aca4SStefano Zampini
979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
9803e5aca4SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc));
99*da34b7cdSBarry Smith if (!dir->hdr.reuseordering) {
1009566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
101f73b0415SBarry Smith if (canuseordering) {
102c294e005SBarry Smith PetscBool external;
103c294e005SBarry Smith
1049566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
1069566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
107c294e005SBarry Smith PetscCall(PetscStrcmp(((PC_Factor *)dir)->ordering, MATORDERINGEXTERNAL, &external));
108c294e005SBarry Smith if (!external) {
1099566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
1101baa6e33SBarry Smith if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
11103c60df9SBarry Smith }
1129b54502bSHong Zhang }
113c294e005SBarry Smith }
1149566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
1159566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1163d1c1ea0SBarry Smith dir->hdr.actualfill = info.fill_ratio_needed;
11704545d6dSBarry Smith } else {
1189566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
119160a8794SBarry Smith if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1209566063dSJacob Faibussowitsch PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
121b8b68cfdSBarry Smith pc->failedreason = PC_NOERROR;
12204545d6dSBarry Smith }
1239b54502bSHong Zhang }
1249566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
12500e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */
12600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1288c1cd74cSHong Zhang }
1298c1cd74cSHong Zhang
1309566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1319566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
13200e125f8SBarry Smith if (err) { /* FactorNumeric() fails */
13300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err;
1348c1cd74cSHong Zhang }
1359b54502bSHong Zhang }
13600c67f3bSHong Zhang
1379566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc, &stype));
13800c67f3bSHong Zhang if (!stype) {
139ea799195SBarry Smith MatSolverType solverpackage;
1409566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1419566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
14200c67f3bSHong Zhang }
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1449b54502bSHong Zhang }
1459b54502bSHong Zhang
PCReset_LU(PC pc)146d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LU(PC pc)
147d71ae5a4SJacob Faibussowitsch {
1489b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data;
1499b54502bSHong Zhang
1509b54502bSHong Zhang PetscFunctionBegin;
1519566063dSJacob Faibussowitsch if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1529566063dSJacob Faibussowitsch if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dir->col));
1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
155574deadeSBarry Smith }
156574deadeSBarry Smith
PCDestroy_LU(PC pc)157d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LU(PC pc)
158d71ae5a4SJacob Faibussowitsch {
159574deadeSBarry Smith PC_LU *dir = (PC_LU *)pc->data;
160574deadeSBarry Smith
161574deadeSBarry Smith PetscFunctionBegin;
1629566063dSJacob Faibussowitsch PetscCall(PCReset_LU(pc));
1639566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1649566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1652e956fe4SStefano Zampini PetscCall(PCFactorClearComposedFunctions(pc));
1669566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1689b54502bSHong Zhang }
1699b54502bSHong Zhang
PCApply_LU(PC pc,Vec x,Vec y)170d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y)
171d71ae5a4SJacob Faibussowitsch {
1729b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data;
1739b54502bSHong Zhang
1749b54502bSHong Zhang PetscFunctionBegin;
1753d1c1ea0SBarry Smith if (dir->hdr.inplace) {
1769566063dSJacob Faibussowitsch PetscCall(MatSolve(pc->pmat, x, y));
1772fa5cd67SKarl Rupp } else {
1789566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1792fa5cd67SKarl Rupp }
1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1819b54502bSHong Zhang }
1829b54502bSHong Zhang
PCMatApply_LU(PC pc,Mat X,Mat Y)183d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y)
184d71ae5a4SJacob Faibussowitsch {
1857b6e2003SPierre Jolivet PC_LU *dir = (PC_LU *)pc->data;
1867b6e2003SPierre Jolivet
1877b6e2003SPierre Jolivet PetscFunctionBegin;
1887b6e2003SPierre Jolivet if (dir->hdr.inplace) {
1899566063dSJacob Faibussowitsch PetscCall(MatMatSolve(pc->pmat, X, Y));
1907b6e2003SPierre Jolivet } else {
1919566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1927b6e2003SPierre Jolivet }
1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1947b6e2003SPierre Jolivet }
1957b6e2003SPierre Jolivet
PCApplyTranspose_LU(PC pc,Vec x,Vec y)196d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y)
197d71ae5a4SJacob Faibussowitsch {
1989b54502bSHong Zhang PC_LU *dir = (PC_LU *)pc->data;
1999b54502bSHong Zhang
2009b54502bSHong Zhang PetscFunctionBegin;
2013d1c1ea0SBarry Smith if (dir->hdr.inplace) {
2029566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(pc->pmat, x, y));
2032fa5cd67SKarl Rupp } else {
2049566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
2052fa5cd67SKarl Rupp }
2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2079b54502bSHong Zhang }
2089b54502bSHong Zhang
PCMatApplyTranspose_LU(PC pc,Mat X,Mat Y)2094dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_LU(PC pc, Mat X, Mat Y)
2104dbf25a8SPierre Jolivet {
2114dbf25a8SPierre Jolivet PC_LU *dir = (PC_LU *)pc->data;
2124dbf25a8SPierre Jolivet
2134dbf25a8SPierre Jolivet PetscFunctionBegin;
2144dbf25a8SPierre Jolivet if (dir->hdr.inplace) {
2154dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(pc->pmat, X, Y));
2164dbf25a8SPierre Jolivet } else {
2174dbf25a8SPierre Jolivet PetscCall(MatMatSolveTranspose(((PC_Factor *)dir)->fact, X, Y));
2184dbf25a8SPierre Jolivet }
2194dbf25a8SPierre Jolivet PetscFunctionReturn(PETSC_SUCCESS);
2204dbf25a8SPierre Jolivet }
2214dbf25a8SPierre Jolivet
2229b54502bSHong Zhang /*MC
2239b54502bSHong Zhang PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2249b54502bSHong Zhang
2259b54502bSHong Zhang Options Database Keys:
226f1580f4eSBarry Smith + -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
227f1580f4eSBarry Smith . -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
228f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
22955ba2a51SBarry Smith . -pc_factor_fill <fill> - Sets fill amount
2302401956bSBarry Smith . -pc_factor_in_place - Activates in-place factorization
2312401956bSBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2322401956bSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2339b54502bSHong Zhang stability of factorization.
234f1580f4eSBarry Smith . -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types
235f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
23626cc229bSBarry Smith . -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
237f1580f4eSBarry Smith . -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities
23826cc229bSBarry Smith - -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
2399b54502bSHong Zhang
2409b54502bSHong Zhang Level: beginner
2419b54502bSHong Zhang
24295452b02SPatrick Sanan Notes:
243f1580f4eSBarry Smith Not all options work for all matrix formats
244f1580f4eSBarry Smith
2450b4b7b1cSBarry Smith Run with `-help` to see additional options for particular matrix formats or factorization algorithms
2460b4b7b1cSBarry Smith
2470b4b7b1cSBarry Smith The Cholesky factorization direct solver, `PCCHOLESKY` will be more efficient than `PCLU` for symmetric positive-definite (SPD) matrices
248f1580f4eSBarry Smith
24995452b02SPatrick Sanan Usually this will compute an "exact" solution in one iteration and does
2509b54502bSHong Zhang not need a Krylov method (i.e. you can use -ksp_type preonly, or
2510b4b7b1cSBarry Smith `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method.
2529b54502bSHong Zhang
253562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`,
254db781477SPatrick Sanan `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
255db781477SPatrick Sanan `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
256c2e3fba1SPatrick Sanan `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
257db781477SPatrick Sanan `PCFactorReorderForNonzeroDiagonal()`
2589b54502bSHong Zhang M*/
2599b54502bSHong Zhang
PCCreate_LU(PC pc)260d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
261d71ae5a4SJacob Faibussowitsch {
2629b54502bSHong Zhang PC_LU *dir;
2639b54502bSHong Zhang
2649b54502bSHong Zhang PetscFunctionBegin;
2654dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dir));
2663d1c1ea0SBarry Smith pc->data = (void *)dir;
2679566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU));
2689b54502bSHong Zhang dir->nonzerosalongdiagonal = PETSC_FALSE;
2699b54502bSHong Zhang
270075768bcSBarry Smith ((PC_Factor *)dir)->info.fill = 5.0;
271075768bcSBarry Smith ((PC_Factor *)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
272f4db908eSBarry Smith ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
2730a545947SLisandro Dalcin dir->col = NULL;
2740a545947SLisandro Dalcin dir->row = NULL;
2755c9eb25fSBarry Smith
276574deadeSBarry Smith pc->ops->reset = PCReset_LU;
2779b54502bSHong Zhang pc->ops->destroy = PCDestroy_LU;
2789b54502bSHong Zhang pc->ops->apply = PCApply_LU;
2797b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_LU;
2809b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_LU;
2814dbf25a8SPierre Jolivet pc->ops->matapplytranspose = PCMatApplyTranspose_LU;
2829b54502bSHong Zhang pc->ops->setup = PCSetUp_LU;
2839b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_LU;
28492e08861SBarry Smith pc->ops->view = PCView_Factor;
2850a545947SLisandro Dalcin pc->ops->applyrichardson = NULL;
2869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU));
2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2889b54502bSHong Zhang }
289