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