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