xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
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));
6303c60df9SBarry Smith       if (dir->row) {
649566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row));
659566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->col));
6603c60df9SBarry Smith       }
679566063dSJacob Faibussowitsch       PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
689566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(pc->pmat, &err));
6900e125f8SBarry Smith       if (err) { /* Factor() fails */
7000e125f8SBarry Smith         pc->failedreason = (PCFailedReason)err;
716baea169SHong Zhang         PetscFunctionReturn(0);
726baea169SHong Zhang       }
739d76b4d0SMatthew G. Knepley     }
74075768bcSBarry Smith     ((PC_Factor *)dir)->fact = pc->pmat;
759b54502bSHong Zhang   } else {
769b54502bSHong Zhang     MatInfo info;
7700e125f8SBarry Smith 
789b54502bSHong Zhang     if (!pc->setupcalled) {
79f73b0415SBarry Smith       PetscBool canuseordering;
802c7c0729SBarry Smith       if (!((PC_Factor *)dir)->fact) {
819566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact));
829566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)dir)->fact));
832c7c0729SBarry Smith       }
849566063dSJacob Faibussowitsch       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
85f73b0415SBarry Smith       if (canuseordering) {
869566063dSJacob Faibussowitsch         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
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));
899566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row));
909566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->col));
9103c60df9SBarry Smith       }
929566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
939566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
943d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
959b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
96f73b0415SBarry Smith       PetscBool canuseordering;
973d1c1ea0SBarry Smith       if (!dir->hdr.reuseordering) {
989566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
999566063dSJacob Faibussowitsch         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact));
1009566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)((PC_Factor *)dir)->fact));
1019566063dSJacob Faibussowitsch         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
102f73b0415SBarry Smith         if (canuseordering) {
1039566063dSJacob Faibussowitsch           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1049566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dir->col));
1059566063dSJacob Faibussowitsch           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
1069566063dSJacob Faibussowitsch           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
1071baa6e33SBarry Smith           if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
1089566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->row));
1099566063dSJacob Faibussowitsch           PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)dir->col));
11003c60df9SBarry Smith         }
1119b54502bSHong Zhang       }
1129566063dSJacob Faibussowitsch       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
1139566063dSJacob Faibussowitsch       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
1143d1c1ea0SBarry Smith       dir->hdr.actualfill = info.fill_ratio_needed;
11504545d6dSBarry Smith     } else {
1169566063dSJacob Faibussowitsch       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
117160a8794SBarry Smith       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
1189566063dSJacob Faibussowitsch         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
119b8b68cfdSBarry Smith         pc->failedreason = PC_NOERROR;
12004545d6dSBarry Smith       }
1219b54502bSHong Zhang     }
1229566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
12300e125f8SBarry Smith     if (err) { /* FactorSymbolic() fails */
12400e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1258c1cd74cSHong Zhang       PetscFunctionReturn(0);
1268c1cd74cSHong Zhang     }
1278c1cd74cSHong Zhang 
1289566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
1299566063dSJacob Faibussowitsch     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
13000e125f8SBarry Smith     if (err) { /* FactorNumeric() fails */
13100e125f8SBarry Smith       pc->failedreason = (PCFailedReason)err;
1328c1cd74cSHong Zhang     }
1339b54502bSHong Zhang   }
13400c67f3bSHong Zhang 
1359566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
13600c67f3bSHong Zhang   if (!stype) {
137ea799195SBarry Smith     MatSolverType solverpackage;
1389566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
1399566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
14000c67f3bSHong Zhang   }
1419b54502bSHong Zhang   PetscFunctionReturn(0);
1429b54502bSHong Zhang }
1439b54502bSHong Zhang 
1449371c9d4SSatish Balay static PetscErrorCode PCReset_LU(PC pc) {
1459b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1469b54502bSHong Zhang 
1479b54502bSHong Zhang   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
1499566063dSJacob Faibussowitsch   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
1509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&dir->col));
151574deadeSBarry Smith   PetscFunctionReturn(0);
152574deadeSBarry Smith }
153574deadeSBarry Smith 
1549371c9d4SSatish Balay static PetscErrorCode PCDestroy_LU(PC pc) {
155574deadeSBarry Smith   PC_LU *dir = (PC_LU *)pc->data;
156574deadeSBarry Smith 
157574deadeSBarry Smith   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(PCReset_LU(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));
1639b54502bSHong Zhang   PetscFunctionReturn(0);
1649b54502bSHong Zhang }
1659b54502bSHong Zhang 
1669371c9d4SSatish Balay static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y) {
1679b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1689b54502bSHong Zhang 
1699b54502bSHong Zhang   PetscFunctionBegin;
1703d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1719566063dSJacob Faibussowitsch     PetscCall(MatSolve(pc->pmat, x, y));
1722fa5cd67SKarl Rupp   } else {
1739566063dSJacob Faibussowitsch     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
1742fa5cd67SKarl Rupp   }
1759b54502bSHong Zhang   PetscFunctionReturn(0);
1769b54502bSHong Zhang }
1779b54502bSHong Zhang 
1789371c9d4SSatish Balay static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y) {
1797b6e2003SPierre Jolivet   PC_LU *dir = (PC_LU *)pc->data;
1807b6e2003SPierre Jolivet 
1817b6e2003SPierre Jolivet   PetscFunctionBegin;
1827b6e2003SPierre Jolivet   if (dir->hdr.inplace) {
1839566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(pc->pmat, X, Y));
1847b6e2003SPierre Jolivet   } else {
1859566063dSJacob Faibussowitsch     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
1867b6e2003SPierre Jolivet   }
1877b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1887b6e2003SPierre Jolivet }
1897b6e2003SPierre Jolivet 
1909371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y) {
1919b54502bSHong Zhang   PC_LU *dir = (PC_LU *)pc->data;
1929b54502bSHong Zhang 
1939b54502bSHong Zhang   PetscFunctionBegin;
1943d1c1ea0SBarry Smith   if (dir->hdr.inplace) {
1959566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(pc->pmat, x, y));
1962fa5cd67SKarl Rupp   } else {
1979566063dSJacob Faibussowitsch     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
1982fa5cd67SKarl Rupp   }
1999b54502bSHong Zhang   PetscFunctionReturn(0);
2009b54502bSHong Zhang }
2019b54502bSHong Zhang 
2029b54502bSHong Zhang /*MC
2039b54502bSHong Zhang    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
2049b54502bSHong Zhang 
2059b54502bSHong Zhang    Options Database Keys:
206*f1580f4eSBarry Smith +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
207*f1580f4eSBarry Smith .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
208*f1580f4eSBarry Smith .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
20955ba2a51SBarry Smith .  -pc_factor_fill <fill> - Sets fill amount
2102401956bSBarry Smith .  -pc_factor_in_place - Activates in-place factorization
2112401956bSBarry Smith .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
2122401956bSBarry Smith .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
2139b54502bSHong Zhang                                          stability of factorization.
214*f1580f4eSBarry Smith .  -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types
215*f1580f4eSBarry Smith .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
21626cc229bSBarry Smith .  -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
217*f1580f4eSBarry Smith .  -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities
21826cc229bSBarry Smith -  -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
2199b54502bSHong Zhang 
2209b54502bSHong Zhang    Level: beginner
2219b54502bSHong Zhang 
22295452b02SPatrick Sanan    Notes:
223*f1580f4eSBarry Smith    Not all options work for all matrix formats
224*f1580f4eSBarry Smith 
225*f1580f4eSBarry Smith    Run with -help to see additional options for particular matrix formats or factorization algorithms
226*f1580f4eSBarry Smith 
22795452b02SPatrick Sanan    Usually this will compute an "exact" solution in one iteration and does
2289b54502bSHong Zhang    not need a Krylov method (i.e. you can use -ksp_type preonly, or
229*f1580f4eSBarry Smith    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
2309b54502bSHong Zhang 
231*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`,
232db781477SPatrick Sanan           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
233db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
234c2e3fba1SPatrick Sanan           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
235db781477SPatrick Sanan           `PCFactorReorderForNonzeroDiagonal()`
2369b54502bSHong Zhang M*/
2379b54502bSHong Zhang 
2389371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) {
2399b54502bSHong Zhang   PC_LU *dir;
2409b54502bSHong Zhang 
2419b54502bSHong Zhang   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &dir));
2433d1c1ea0SBarry Smith   pc->data = (void *)dir;
2449566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU));
2459b54502bSHong Zhang   dir->nonzerosalongdiagonal = PETSC_FALSE;
2469b54502bSHong Zhang 
247075768bcSBarry Smith   ((PC_Factor *)dir)->info.fill      = 5.0;
248075768bcSBarry Smith   ((PC_Factor *)dir)->info.dtcol     = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
249f4db908eSBarry Smith   ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
2500a545947SLisandro Dalcin   dir->col                           = NULL;
2510a545947SLisandro Dalcin   dir->row                           = NULL;
2525c9eb25fSBarry Smith 
253574deadeSBarry Smith   pc->ops->reset           = PCReset_LU;
2549b54502bSHong Zhang   pc->ops->destroy         = PCDestroy_LU;
2559b54502bSHong Zhang   pc->ops->apply           = PCApply_LU;
2567b6e2003SPierre Jolivet   pc->ops->matapply        = PCMatApply_LU;
2579b54502bSHong Zhang   pc->ops->applytranspose  = PCApplyTranspose_LU;
2589b54502bSHong Zhang   pc->ops->setup           = PCSetUp_LU;
2599b54502bSHong Zhang   pc->ops->setfromoptions  = PCSetFromOptions_LU;
26092e08861SBarry Smith   pc->ops->view            = PCView_Factor;
2610a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
2629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU));
2639b54502bSHong Zhang   PetscFunctionReturn(0);
2649b54502bSHong Zhang }
265