xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
19b54502bSHong Zhang /*
29b54502bSHong Zhang    Defines a direct factorization preconditioner for any Mat implementation
39b54502bSHong Zhang    Note: this need not be consided 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 
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 
20d71ae5a4SJacob Faibussowitsch 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 
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));
654dfa11a4SJacob Faibussowitsch       if (dir->row) { }
669566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
679566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat, &err));
6800e125f8SBarry Smith       if (err) { /* Factor() fails */
6900e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
703ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
716baea169SHong Zhang       }
729d76b4d0SMatthew G. Knepley     }
73075768bcSBarry Smith     ((PC_Factor *)dir)->fact = pc->pmat;
749b54502bSHong Zhang   } else {
759b54502bSHong Zhang     MatInfo info;
7600e125f8SBarry Smith 
779b54502bSHong Zhang     if (!pc->setupcalled) {
78f73b0415SBarry Smith       PetscBool canuseordering;
7903e5aca4SStefano Zampini 
8003e5aca4SStefano Zampini       PetscCall(PCFactorSetUpMatSolverType(pc));
819566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
82f73b0415SBarry Smith       if (canuseordering) {
839566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
849566063dSJacob Faibussowitsch         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
851baa6e33SBarry Smith         if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
8603c60df9SBarry Smith       }
879566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
889566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
893d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
909b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
91f73b0415SBarry Smith       PetscBool canuseordering;
9203e5aca4SStefano Zampini 
933d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
949566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
9503e5aca4SStefano Zampini         PetscCall(PCFactorSetUpMatSolverType(pc));
969566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
97f73b0415SBarry Smith         if (canuseordering) {
989566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
999566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
1009566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1019566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
1021baa6e33SBarry Smith           if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
10303c60df9SBarry Smith         }
1049b54502bSHong Zhang       }
1059566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
1069566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1073d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
10804545d6dSBarry Smith     } else {
1099566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
110160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1119566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
112b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
11304545d6dSBarry Smith       }
1149b54502bSHong Zhang     }
1159566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
11600e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
11700e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1183ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1198c1cd74cSHong Zhang     }
1208c1cd74cSHong Zhang 
1219566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1229566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
12300e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
12400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1258c1cd74cSHong Zhang     }
1269b54502bSHong Zhang   }
12700c67f3bSHong Zhang 
1289566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
12900c67f3bSHong Zhang   if (!stype) {
130ea799195SBarry Smith     MatSolverType solverpackage;
1319566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1329566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
13300c67f3bSHong Zhang   }
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1359b54502bSHong Zhang }
1369b54502bSHong Zhang 
137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_LU(PC pc)
138d71ae5a4SJacob Faibussowitsch {
1399b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1409b54502bSHong Zhang 
1419b54502bSHong Zhang   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1439566063dSJacob Faibussowitsch   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146574deadeSBarry Smith }
147574deadeSBarry Smith 
148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_LU(PC pc)
149d71ae5a4SJacob Faibussowitsch {
150574deadeSBarry Smith   PC_LU *dir = (PC_LU *)pc->data;
151574deadeSBarry Smith 
152574deadeSBarry Smith   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(PCReset_LU(pc));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
1562e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
1579566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1599b54502bSHong Zhang }
1609b54502bSHong Zhang 
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y)
162d71ae5a4SJacob Faibussowitsch {
1639b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1649b54502bSHong Zhang 
1659b54502bSHong Zhang   PetscFunctionBegin;
1663d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1679566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat, x, y));
1682fa5cd67SKarl Rupp   } else {
1699566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1702fa5cd67SKarl Rupp   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1729b54502bSHong Zhang }
1739b54502bSHong Zhang 
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y)
175d71ae5a4SJacob Faibussowitsch {
1767b6e2003SPierre Jolivet   PC_LU *dir = (PC_LU *)pc->data;
1777b6e2003SPierre Jolivet 
1787b6e2003SPierre Jolivet   PetscFunctionBegin;
1797b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1809566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat, X, Y));
1817b6e2003SPierre Jolivet   } else {
1829566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1837b6e2003SPierre Jolivet   }
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1857b6e2003SPierre Jolivet }
1867b6e2003SPierre Jolivet 
187d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y)
188d71ae5a4SJacob Faibussowitsch {
1899b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1909b54502bSHong Zhang 
1919b54502bSHong Zhang   PetscFunctionBegin;
1923d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1939566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat, x, y));
1942fa5cd67SKarl Rupp   } else {
1959566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
1962fa5cd67SKarl Rupp   }
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1989b54502bSHong Zhang }
1999b54502bSHong Zhang 
2009b54502bSHong Zhang /*MC
2019b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2029b54502bSHong Zhang 
2039b54502bSHong Zhang    Options Database Keys:
204f1580f4eSBarry Smith +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
205f1580f4eSBarry Smith .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
206f1580f4eSBarry Smith .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
20755ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2082401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2092401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2102401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2119b54502bSHong Zhang                                          stability of factorization.
212f1580f4eSBarry Smith .  -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types
213f1580f4eSBarry Smith .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
21426cc229bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
215f1580f4eSBarry Smith .  -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities
21626cc229bSBarry Smith -  -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
2179b54502bSHong Zhang 
2189b54502bSHong Zhang    Level: beginner
2199b54502bSHong Zhang 
22095452b02SPatrick Sanan    Notes:
221f1580f4eSBarry Smith    Not all options work for all matrix formats
222f1580f4eSBarry Smith 
223f1580f4eSBarry Smith    Run with -help to see additional options for particular matrix formats or factorization algorithms
224f1580f4eSBarry Smith 
22595452b02SPatrick Sanan    Usually this will compute an "exact" solution in one iteration and does
2269b54502bSHong Zhang    not need a Krylov method (i.e. you can use -ksp_type preonly, or
227f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
2289b54502bSHong Zhang 
229*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`,
230db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
231db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
232c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
233db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2349b54502bSHong Zhang M*/
2359b54502bSHong Zhang 
236d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
237d71ae5a4SJacob Faibussowitsch {
2389b54502bSHong Zhang   PC_LU *dir;
2399b54502bSHong Zhang 
2409b54502bSHong Zhang   PetscFunctionBegin;
2414dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&dir));
2423d1c1ea0SBarry Smith   pc->data = (void *)dir;
2439566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU));
2449b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
2459b54502bSHong Zhang 
246075768bcSBarry Smith   ((PC_Factor *)dir)->info.fill      = 5.0;
247075768bcSBarry Smith   ((PC_Factor *)dir)->info.dtcol     = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
248f4db908eSBarry Smith   ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
2490a545947SLisandro Dalcin   dir->col                           = NULL;
2500a545947SLisandro Dalcin   dir->row                           = NULL;
2515c9eb25fSBarry Smith 
252574deadeSBarry Smith   pc->ops->reset           = PCReset_LU;
2539b54502bSHong Zhang   pc->ops->destroy         = PCDestroy_LU;
2549b54502bSHong Zhang   pc->ops->apply           = PCApply_LU;
2557b6e2003SPierre Jolivet   pc->ops->matapply        = PCMatApply_LU;
2569b54502bSHong Zhang   pc->ops->applytranspose  = PCApplyTranspose_LU;
2579b54502bSHong Zhang   pc->ops->setup           = PCSetUp_LU;
2589b54502bSHong Zhang   pc->ops->setfromoptions  = PCSetFromOptions_LU;
25992e08861SBarry Smith   pc->ops->view            = PCView_Factor;
2600a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2639b54502bSHong Zhang }
264