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