xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1dba47a55SKris Buschelman 
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 */
7ee45ca4aSHong Zhang 
8c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/
99b54502bSHong Zhang 
109371c9d4SSatish Balay PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc, PetscReal z) {
119b54502bSHong Zhang   PC_LU *lu = (PC_LU *)pc->data;
129b54502bSHong Zhang 
139b54502bSHong Zhang   PetscFunctionBegin;
149b54502bSHong Zhang   lu->nonzerosalongdiagonal = PETSC_TRUE;
152fa5cd67SKarl Rupp   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
162fa5cd67SKarl Rupp   else lu->nonzerosalongdiagonaltol = z;
179b54502bSHong Zhang   PetscFunctionReturn(0);
189b54502bSHong Zhang }
199b54502bSHong Zhang 
209371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_LU(PC pc, PetscOptionItems *PetscOptionsObject) {
219b54502bSHong Zhang   PC_LU    *lu  = (PC_LU *)pc->data;
22ace3abfcSBarry Smith   PetscBool flg = PETSC_FALSE;
239b54502bSHong Zhang   PetscReal tol;
249b54502bSHong Zhang 
259b54502bSHong Zhang   PetscFunctionBegin;
26d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "LU options");
27dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
285c9eb25fSBarry Smith 
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
309b54502bSHong Zhang   if (flg) {
319b54502bSHong Zhang     tol = PETSC_DECIDE;
329566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", lu->nonzerosalongdiagonaltol, &tol, NULL));
339566063dSJacob Faibussowitsch     PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
349b54502bSHong Zhang   }
35d0609cedSBarry Smith   PetscOptionsHeadEnd();
369b54502bSHong Zhang   PetscFunctionReturn(0);
379b54502bSHong Zhang }
389b54502bSHong Zhang 
399371c9d4SSatish Balay static PetscErrorCode PCSetUp_LU(PC pc) {
409b54502bSHong Zhang   PC_LU         *dir = (PC_LU *)pc->data;
41ea799195SBarry Smith   MatSolverType  stype;
4200e125f8SBarry Smith   MatFactorError err;
43f023e1d5SPierre Jolivet   const char    *prefix;
443d1c1ea0SBarry Smith 
459b54502bSHong Zhang   PetscFunctionBegin;
46c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
473d1c1ea0SBarry Smith   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
489b54502bSHong Zhang 
4926cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
5026cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
5126cc229bSBarry Smith 
529566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
533d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
549d76b4d0SMatthew G. Knepley     MatFactorType ftype;
559d76b4d0SMatthew G. Knepley 
569566063dSJacob Faibussowitsch     PetscCall(MatGetFactorType(pc->pmat, &ftype));
579d76b4d0SMatthew G. Knepley     if (ftype == MAT_FACTOR_NONE) {
589566063dSJacob Faibussowitsch       if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
599566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&dir->col));
602c7c0729SBarry Smith       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
619566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
629566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
63*4dfa11a4SJacob Faibussowitsch       if (dir->row) { }
649566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
659566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat, &err));
6600e125f8SBarry Smith       if (err) { /* Factor() fails */
6700e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
686baea169SHong Zhang         PetscFunctionReturn(0);
696baea169SHong Zhang       }
709d76b4d0SMatthew G. Knepley     }
71075768bcSBarry Smith     ((PC_Factor *)dir)->fact = pc->pmat;
729b54502bSHong Zhang   } else {
739b54502bSHong Zhang     MatInfo info;
7400e125f8SBarry Smith 
759b54502bSHong Zhang     if (!pc->setupcalled) {
76f73b0415SBarry Smith       PetscBool canuseordering;
77*4dfa11a4SJacob Faibussowitsch       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact)); }
789566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
79f73b0415SBarry Smith       if (canuseordering) {
809566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
819566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
821baa6e33SBarry Smith         if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
8303c60df9SBarry Smith       }
849566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((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) {
88f73b0415SBarry Smith       PetscBool canuseordering;
893d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
909566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
919566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact));
929566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
93f73b0415SBarry Smith         if (canuseordering) {
949566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
959566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
969566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
979566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
981baa6e33SBarry Smith           if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
9903c60df9SBarry Smith         }
1009b54502bSHong Zhang       }
1019566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
1029566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1033d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
10404545d6dSBarry Smith     } else {
1059566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
106160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1079566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
108b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
10904545d6dSBarry Smith       }
1109b54502bSHong Zhang     }
1119566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
11200e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
11300e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1148c1cd74cSHong Zhang       PetscFunctionReturn(0);
1158c1cd74cSHong Zhang     }
1168c1cd74cSHong Zhang 
1179566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1189566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
11900e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12000e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1218c1cd74cSHong Zhang     }
1229b54502bSHong Zhang   }
12300c67f3bSHong Zhang 
1249566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
12500c67f3bSHong Zhang   if (!stype) {
126ea799195SBarry Smith     MatSolverType solverpackage;
1279566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1289566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
12900c67f3bSHong Zhang   }
1309b54502bSHong Zhang   PetscFunctionReturn(0);
1319b54502bSHong Zhang }
1329b54502bSHong Zhang 
1339371c9d4SSatish Balay static PetscErrorCode PCReset_LU(PC pc) {
1349b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1359b54502bSHong Zhang 
1369b54502bSHong Zhang   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1389566063dSJacob Faibussowitsch   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
140574deadeSBarry Smith   PetscFunctionReturn(0);
141574deadeSBarry Smith }
142574deadeSBarry Smith 
1439371c9d4SSatish Balay static PetscErrorCode PCDestroy_LU(PC pc) {
144574deadeSBarry Smith   PC_LU *dir = (PC_LU *)pc->data;
145574deadeSBarry Smith 
146574deadeSBarry Smith   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(PCReset_LU(pc));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1499566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1502e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1529b54502bSHong Zhang   PetscFunctionReturn(0);
1539b54502bSHong Zhang }
1549b54502bSHong Zhang 
1559371c9d4SSatish Balay static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y) {
1569b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1579b54502bSHong Zhang 
1589b54502bSHong Zhang   PetscFunctionBegin;
1593d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1609566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat, x, y));
1612fa5cd67SKarl Rupp   } else {
1629566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1632fa5cd67SKarl Rupp   }
1649b54502bSHong Zhang   PetscFunctionReturn(0);
1659b54502bSHong Zhang }
1669b54502bSHong Zhang 
1679371c9d4SSatish Balay static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y) {
1687b6e2003SPierre Jolivet   PC_LU *dir = (PC_LU *)pc->data;
1697b6e2003SPierre Jolivet 
1707b6e2003SPierre Jolivet   PetscFunctionBegin;
1717b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1729566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat, X, Y));
1737b6e2003SPierre Jolivet   } else {
1749566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1757b6e2003SPierre Jolivet   }
1767b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1777b6e2003SPierre Jolivet }
1787b6e2003SPierre Jolivet 
1799371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y) {
1809b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1819b54502bSHong Zhang 
1829b54502bSHong Zhang   PetscFunctionBegin;
1833d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1849566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat, x, y));
1852fa5cd67SKarl Rupp   } else {
1869566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
1872fa5cd67SKarl Rupp   }
1889b54502bSHong Zhang   PetscFunctionReturn(0);
1899b54502bSHong Zhang }
1909b54502bSHong Zhang 
1919b54502bSHong Zhang /*MC
1929b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
1939b54502bSHong Zhang 
1949b54502bSHong Zhang    Options Database Keys:
195f1580f4eSBarry Smith +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
196f1580f4eSBarry Smith .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
197f1580f4eSBarry Smith .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
19855ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
1992401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2002401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2012401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2029b54502bSHong Zhang                                          stability of factorization.
203f1580f4eSBarry Smith .  -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types
204f1580f4eSBarry Smith .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
20526cc229bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
206f1580f4eSBarry Smith .  -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities
20726cc229bSBarry Smith -  -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
2089b54502bSHong Zhang 
2099b54502bSHong Zhang    Level: beginner
2109b54502bSHong Zhang 
21195452b02SPatrick Sanan    Notes:
212f1580f4eSBarry Smith    Not all options work for all matrix formats
213f1580f4eSBarry Smith 
214f1580f4eSBarry Smith    Run with -help to see additional options for particular matrix formats or factorization algorithms
215f1580f4eSBarry Smith 
21695452b02SPatrick Sanan    Usually this will compute an "exact" solution in one iteration and does
2179b54502bSHong Zhang    not need a Krylov method (i.e. you can use -ksp_type preonly, or
218f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
2199b54502bSHong Zhang 
220f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`,
221db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
222db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
223c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
224db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2259b54502bSHong Zhang M*/
2269b54502bSHong Zhang 
2279371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) {
2289b54502bSHong Zhang   PC_LU *dir;
2299b54502bSHong Zhang 
2309b54502bSHong Zhang   PetscFunctionBegin;
231*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dir));
2323d1c1ea0SBarry Smith   pc->data = (void *)dir;
2339566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU));
2349b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
2359b54502bSHong Zhang 
236075768bcSBarry Smith   ((PC_Factor *)dir)->info.fill      = 5.0;
237075768bcSBarry Smith   ((PC_Factor *)dir)->info.dtcol     = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
238f4db908eSBarry Smith   ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
2390a545947SLisandro Dalcin   dir->col                           = NULL;
2400a545947SLisandro Dalcin   dir->row                           = NULL;
2415c9eb25fSBarry Smith 
242574deadeSBarry Smith   pc->ops->reset           = PCReset_LU;
2439b54502bSHong Zhang   pc->ops->destroy         = PCDestroy_LU;
2449b54502bSHong Zhang   pc->ops->apply           = PCApply_LU;
2457b6e2003SPierre Jolivet   pc->ops->matapply        = PCMatApply_LU;
2469b54502bSHong Zhang   pc->ops->applytranspose  = PCApplyTranspose_LU;
2479b54502bSHong Zhang   pc->ops->setup           = PCSetUp_LU;
2489b54502bSHong Zhang   pc->ops->setfromoptions  = PCSetFromOptions_LU;
24992e08861SBarry Smith   pc->ops->view            = PCView_Factor;
2500a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU));
2529b54502bSHong Zhang   PetscFunctionReturn(0);
2539b54502bSHong Zhang }
254