xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 750b007cd8d816cecd9de99077bb0a703b4cf61a)
1 
2 /*
3    Defines a direct factorization preconditioner for any Mat implementation
4    Note: this need not be consided a preconditioner since it supplies
5          a direct solver.
6 */
7 
8 #include <../src/ksp/pc/impls/factor/lu/lu.h> /*I "petscpc.h" I*/
9 
10 PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc, PetscReal z) {
11   PC_LU *lu = (PC_LU *)pc->data;
12 
13   PetscFunctionBegin;
14   lu->nonzerosalongdiagonal = PETSC_TRUE;
15   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
16   else lu->nonzerosalongdiagonaltol = z;
17   PetscFunctionReturn(0);
18 }
19 
20 static PetscErrorCode PCSetFromOptions_LU(PC pc, PetscOptionItems *PetscOptionsObject) {
21   PC_LU    *lu  = (PC_LU *)pc->data;
22   PetscBool flg = PETSC_FALSE;
23   PetscReal tol;
24 
25   PetscFunctionBegin;
26   PetscOptionsHeadBegin(PetscOptionsObject, "LU options");
27   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
28 
29   PetscCall(PetscOptionsName("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", &flg));
30   if (flg) {
31     tol = PETSC_DECIDE;
32     PetscCall(PetscOptionsReal("-pc_factor_nonzeros_along_diagonal", "Reorder to remove zeros from diagonal", "PCFactorReorderForNonzeroDiagonal", lu->nonzerosalongdiagonaltol, &tol, NULL));
33     PetscCall(PCFactorReorderForNonzeroDiagonal(pc, tol));
34   }
35   PetscOptionsHeadEnd();
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode PCSetUp_LU(PC pc) {
40   PC_LU         *dir = (PC_LU *)pc->data;
41   MatSolverType  stype;
42   MatFactorError err;
43   const char    *prefix;
44 
45   PetscFunctionBegin;
46   pc->failedreason = PC_NOERROR;
47   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor *)dir)->info.fill = dir->hdr.actualfill;
48 
49   PetscCall(PCGetOptionsPrefix(pc, &prefix));
50   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
51 
52   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
53   if (dir->hdr.inplace) {
54     MatFactorType ftype;
55 
56     PetscCall(MatGetFactorType(pc->pmat, &ftype));
57     if (ftype == MAT_FACTOR_NONE) {
58       if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
59       PetscCall(ISDestroy(&dir->col));
60       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
61       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
62       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
63       if (dir->row) { }
64       PetscCall(MatLUFactor(pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
65       PetscCall(MatFactorGetError(pc->pmat, &err));
66       if (err) { /* Factor() fails */
67         pc->failedreason = (PCFailedReason)err;
68         PetscFunctionReturn(0);
69       }
70     }
71     ((PC_Factor *)dir)->fact = pc->pmat;
72   } else {
73     MatInfo info;
74 
75     if (!pc->setupcalled) {
76       PetscBool canuseordering;
77       if (!((PC_Factor *)dir)->fact) { PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact)); }
78       PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
79       if (canuseordering) {
80         PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
81         PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
82         if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
83       }
84       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
85       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
86       dir->hdr.actualfill = info.fill_ratio_needed;
87     } else if (pc->flag != SAME_NONZERO_PATTERN) {
88       PetscBool canuseordering;
89       if (!dir->hdr.reuseordering) {
90         PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
91         PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)dir)->solvertype, MAT_FACTOR_LU, &((PC_Factor *)dir)->fact));
92         PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)dir)->fact, &canuseordering));
93         if (canuseordering) {
94           if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
95           PetscCall(ISDestroy(&dir->col));
96           PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
97           PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)dir)->ordering, &dir->row, &dir->col));
98           if (dir->nonzerosalongdiagonal) PetscCall(MatReorderForNonzeroDiagonal(pc->pmat, dir->nonzerosalongdiagonaltol, dir->row, dir->col));
99         }
100       }
101       PetscCall(MatLUFactorSymbolic(((PC_Factor *)dir)->fact, pc->pmat, dir->row, dir->col, &((PC_Factor *)dir)->info));
102       PetscCall(MatGetInfo(((PC_Factor *)dir)->fact, MAT_LOCAL, &info));
103       dir->hdr.actualfill = info.fill_ratio_needed;
104     } else {
105       PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
106       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
107         PetscCall(MatFactorClearError(((PC_Factor *)dir)->fact));
108         pc->failedreason = PC_NOERROR;
109       }
110     }
111     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
112     if (err) { /* FactorSymbolic() fails */
113       pc->failedreason = (PCFailedReason)err;
114       PetscFunctionReturn(0);
115     }
116 
117     PetscCall(MatLUFactorNumeric(((PC_Factor *)dir)->fact, pc->pmat, &((PC_Factor *)dir)->info));
118     PetscCall(MatFactorGetError(((PC_Factor *)dir)->fact, &err));
119     if (err) { /* FactorNumeric() fails */
120       pc->failedreason = (PCFailedReason)err;
121     }
122   }
123 
124   PetscCall(PCFactorGetMatSolverType(pc, &stype));
125   if (!stype) {
126     MatSolverType solverpackage;
127     PetscCall(MatFactorGetSolverType(((PC_Factor *)dir)->fact, &solverpackage));
128     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
129   }
130   PetscFunctionReturn(0);
131 }
132 
133 static PetscErrorCode PCReset_LU(PC pc) {
134   PC_LU *dir = (PC_LU *)pc->data;
135 
136   PetscFunctionBegin;
137   if (!dir->hdr.inplace && ((PC_Factor *)dir)->fact) PetscCall(MatDestroy(&((PC_Factor *)dir)->fact));
138   if (dir->row && dir->col && dir->row != dir->col) PetscCall(ISDestroy(&dir->row));
139   PetscCall(ISDestroy(&dir->col));
140   PetscFunctionReturn(0);
141 }
142 
143 static PetscErrorCode PCDestroy_LU(PC pc) {
144   PC_LU *dir = (PC_LU *)pc->data;
145 
146   PetscFunctionBegin;
147   PetscCall(PCReset_LU(pc));
148   PetscCall(PetscFree(((PC_Factor *)dir)->ordering));
149   PetscCall(PetscFree(((PC_Factor *)dir)->solvertype));
150   PetscCall(PCFactorClearComposedFunctions(pc));
151   PetscCall(PetscFree(pc->data));
152   PetscFunctionReturn(0);
153 }
154 
155 static PetscErrorCode PCApply_LU(PC pc, Vec x, Vec y) {
156   PC_LU *dir = (PC_LU *)pc->data;
157 
158   PetscFunctionBegin;
159   if (dir->hdr.inplace) {
160     PetscCall(MatSolve(pc->pmat, x, y));
161   } else {
162     PetscCall(MatSolve(((PC_Factor *)dir)->fact, x, y));
163   }
164   PetscFunctionReturn(0);
165 }
166 
167 static PetscErrorCode PCMatApply_LU(PC pc, Mat X, Mat Y) {
168   PC_LU *dir = (PC_LU *)pc->data;
169 
170   PetscFunctionBegin;
171   if (dir->hdr.inplace) {
172     PetscCall(MatMatSolve(pc->pmat, X, Y));
173   } else {
174     PetscCall(MatMatSolve(((PC_Factor *)dir)->fact, X, Y));
175   }
176   PetscFunctionReturn(0);
177 }
178 
179 static PetscErrorCode PCApplyTranspose_LU(PC pc, Vec x, Vec y) {
180   PC_LU *dir = (PC_LU *)pc->data;
181 
182   PetscFunctionBegin;
183   if (dir->hdr.inplace) {
184     PetscCall(MatSolveTranspose(pc->pmat, x, y));
185   } else {
186     PetscCall(MatSolveTranspose(((PC_Factor *)dir)->fact, x, y));
187   }
188   PetscFunctionReturn(0);
189 }
190 
191 /*MC
192    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
193 
194    Options Database Keys:
195 +  -pc_factor_reuse_ordering - Activate `PCFactorSetReuseOrdering()`
196 .  -pc_factor_mat_solver_type - Actives `PCFactorSetMatSolverType()` to choose the direct solver, like superlu
197 .  -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()`
198 .  -pc_factor_fill <fill> - Sets fill amount
199 .  -pc_factor_in_place - Activates in-place factorization
200 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
201 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
202                                          stability of factorization.
203 .  -pc_factor_shift_type <shifttype> - Sets shift type or -1 for the default; use '-help' for a list of available types
204 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default
205 .  -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the diagonal.
206 .  -pc_factor_mat_solver_type <packagename> - use an external package for the solve, see `MatSolverType` for possibilities
207 -  -mat_solvertype_optionname - options for a specific solver package, for example -mat_mumps_cntl_1
208 
209    Level: beginner
210 
211    Notes:
212    Not all options work for all matrix formats
213 
214    Run with -help to see additional options for particular matrix formats or factorization algorithms
215 
216    Usually this will compute an "exact" solution in one iteration and does
217    not need a Krylov method (i.e. you can use -ksp_type preonly, or
218    `KSPSetType`(ksp,`KSPPREONLY`) for the Krylov method
219 
220 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MatSolverType`, `MatGetFactor()`, `PCQR`, `PCSVD`,
221           `PCILU`, `PCCHOLESKY`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetReuseFill()`, `PCFactorGetMatrix()`,
222           `PCFactorSetFill()`, `PCFactorSetUseInPlace()`, `PCFactorSetMatOrderingType()`, `PCFactorSetColumnPivot()`,
223           `PCFactorSetPivotInBlocks()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`
224           `PCFactorReorderForNonzeroDiagonal()`
225 M*/
226 
227 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) {
228   PC_LU *dir;
229 
230   PetscFunctionBegin;
231   PetscCall(PetscNew(&dir));
232   pc->data = (void *)dir;
233   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_LU));
234   dir->nonzerosalongdiagonal = PETSC_FALSE;
235 
236   ((PC_Factor *)dir)->info.fill      = 5.0;
237   ((PC_Factor *)dir)->info.dtcol     = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
238   ((PC_Factor *)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
239   dir->col                           = NULL;
240   dir->row                           = NULL;
241 
242   pc->ops->reset           = PCReset_LU;
243   pc->ops->destroy         = PCDestroy_LU;
244   pc->ops->apply           = PCApply_LU;
245   pc->ops->matapply        = PCMatApply_LU;
246   pc->ops->applytranspose  = PCApplyTranspose_LU;
247   pc->ops->setup           = PCSetUp_LU;
248   pc->ops->setfromoptions  = PCSetFromOptions_LU;
249   pc->ops->view            = PCView_Factor;
250   pc->ops->applyrichardson = NULL;
251   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", PCFactorReorderForNonzeroDiagonal_LU));
252   PetscFunctionReturn(0);
253 }
254