xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision e6e75211d226c622f451867f53ce5d558649ff4f)
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 #undef __FUNCT__
11 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU"
12 PetscErrorCode  PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
13 {
14   PC_LU *lu = (PC_LU*)pc->data;
15 
16   PetscFunctionBegin;
17   lu->nonzerosalongdiagonal = PETSC_TRUE;
18   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
19   else lu->nonzerosalongdiagonaltol = z;
20   PetscFunctionReturn(0);
21 }
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "PCFactorSetReuseOrdering_LU"
25 PetscErrorCode  PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag)
26 {
27   PC_LU *lu = (PC_LU*)pc->data;
28 
29   PetscFunctionBegin;
30   lu->reuseordering = flag;
31   PetscFunctionReturn(0);
32 }
33 
34 #undef __FUNCT__
35 #define __FUNCT__ "PCFactorSetReuseFill_LU"
36 PetscErrorCode  PCFactorSetReuseFill_LU(PC pc,PetscBool flag)
37 {
38   PC_LU *lu = (PC_LU*)pc->data;
39 
40   PetscFunctionBegin;
41   lu->reusefill = flag;
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "PCSetFromOptions_LU"
47 static PetscErrorCode PCSetFromOptions_LU(PetscOptions *PetscOptionsObject,PC pc)
48 {
49   PC_LU          *lu = (PC_LU*)pc->data;
50   PetscErrorCode ierr;
51   PetscBool      flg = PETSC_FALSE;
52   PetscReal      tol;
53 
54   PetscFunctionBegin;
55   ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr);
56   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
57 
58   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
59   if (flg) {
60     tol  = PETSC_DECIDE;
61     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
62     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
63   }
64   ierr = PetscOptionsTail();CHKERRQ(ierr);
65   PetscFunctionReturn(0);
66 }
67 
68 #undef __FUNCT__
69 #define __FUNCT__ "PCView_LU"
70 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
71 {
72   PC_LU          *lu = (PC_LU*)pc->data;
73   PetscErrorCode ierr;
74   PetscBool      iascii;
75 
76   PetscFunctionBegin;
77   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
78   if (iascii) {
79     if (lu->inplace) {
80       ierr = PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");CHKERRQ(ierr);
81     } else {
82       ierr = PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");CHKERRQ(ierr);
83     }
84 
85     if (lu->reusefill)    {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
86     if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
87   }
88   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "PCSetUp_LU"
94 static PetscErrorCode PCSetUp_LU(PC pc)
95 {
96   PetscErrorCode ierr;
97   PC_LU          *dir = (PC_LU*)pc->data;
98 
99   PetscFunctionBegin;
100   if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill;
101 
102   if (dir->inplace) {
103     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
104     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
105     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
106     if (dir->row) {
107       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
108       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
109     }
110     ierr = MatSetErrorIfFPE(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
111     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
112     ((PC_Factor*)dir)->fact = pc->pmat;
113   } else {
114     MatInfo info;
115     if (!pc->setupcalled) {
116       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
117       if (dir->nonzerosalongdiagonal) {
118         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
119       }
120       if (dir->row) {
121         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
122         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
123       }
124       if (!((PC_Factor*)dir)->fact) {
125         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
126       }
127       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
128       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
129       dir->actualfill = info.fill_ratio_needed;
130       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
131     } else if (pc->flag != SAME_NONZERO_PATTERN) {
132       if (!dir->reuseordering) {
133         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
134         ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
135         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
136         if (dir->nonzerosalongdiagonal) {
137           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
138         }
139         if (dir->row) {
140           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
141           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
142         }
143       }
144       ierr            = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
145       ierr            = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
146       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
147       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
148       dir->actualfill = info.fill_ratio_needed;
149       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
150     }
151     ierr = MatSetErrorIfFPE(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
152     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 #undef __FUNCT__
158 #define __FUNCT__ "PCReset_LU"
159 static PetscErrorCode PCReset_LU(PC pc)
160 {
161   PC_LU          *dir = (PC_LU*)pc->data;
162   PetscErrorCode ierr;
163 
164   PetscFunctionBegin;
165   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
166   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
167   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "PCDestroy_LU"
173 static PetscErrorCode PCDestroy_LU(PC pc)
174 {
175   PC_LU          *dir = (PC_LU*)pc->data;
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   ierr = PCReset_LU(pc);CHKERRQ(ierr);
180   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
181   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
182   ierr = PetscFree(pc->data);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PCApply_LU"
188 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
189 {
190   PC_LU          *dir = (PC_LU*)pc->data;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   if (dir->inplace) {
195     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
196   } else {
197     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "PCApplyTranspose_LU"
204 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
205 {
206   PC_LU          *dir = (PC_LU*)pc->data;
207   PetscErrorCode ierr;
208 
209   PetscFunctionBegin;
210   if (dir->inplace) {
211     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
212   } else {
213     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
214   }
215   PetscFunctionReturn(0);
216 }
217 
218 /* -----------------------------------------------------------------------------------*/
219 
220 #undef __FUNCT__
221 #define __FUNCT__ "PCFactorSetUseInPlace_LU"
222 PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
223 {
224   PC_LU *dir = (PC_LU*)pc->data;
225 
226   PetscFunctionBegin;
227   dir->inplace = flg;
228   PetscFunctionReturn(0);
229 }
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "PCFactorGetUseInPlace_LU"
233 PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
234 {
235   PC_LU *dir = (PC_LU*)pc->data;
236 
237   PetscFunctionBegin;
238   *flg = dir->inplace;
239   PetscFunctionReturn(0);
240 }
241 
242 /* ------------------------------------------------------------------------ */
243 
244 /*MC
245    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
246 
247    Options Database Keys:
248 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
249 .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
250 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
251 .  -pc_factor_fill <fill> - Sets fill amount
252 .  -pc_factor_in_place - Activates in-place factorization
253 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
254 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
255                                          stability of factorization.
256 .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
257 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
258 -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
259         diagonal.
260 
261    Notes: Not all options work for all matrix formats
262           Run with -help to see additional options for particular matrix formats or factorization
263           algorithms
264 
265    Level: beginner
266 
267    Concepts: LU factorization, direct solver
268 
269    Notes: Usually this will compute an "exact" solution in one iteration and does
270           not need a Krylov method (i.e. you can use -ksp_type preonly, or
271           KSPSetType(ksp,KSPPREONLY) for the Krylov method
272 
273 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
274            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
275            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
276            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
277            PCFactorReorderForNonzeroDiagonal()
278 M*/
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "PCCreate_LU"
282 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
283 {
284   PetscErrorCode ierr;
285   PetscMPIInt    size;
286   PC_LU          *dir;
287 
288   PetscFunctionBegin;
289   ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr);
290 
291   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
292 
293   ((PC_Factor*)dir)->fact       = NULL;
294   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
295   dir->inplace                  = PETSC_FALSE;
296   dir->nonzerosalongdiagonal    = PETSC_FALSE;
297 
298   ((PC_Factor*)dir)->info.fill          = 5.0;
299   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
300   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
301   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
302   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
303   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
304   dir->col                              = 0;
305   dir->row                              = 0;
306 
307   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */
308   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
309   if (size == 1) {
310     ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
311   } else {
312     ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
313   }
314   dir->reusefill     = PETSC_FALSE;
315   dir->reuseordering = PETSC_FALSE;
316   pc->data           = (void*)dir;
317 
318   pc->ops->reset             = PCReset_LU;
319   pc->ops->destroy           = PCDestroy_LU;
320   pc->ops->apply             = PCApply_LU;
321   pc->ops->applytranspose    = PCApplyTranspose_LU;
322   pc->ops->setup             = PCSetUp_LU;
323   pc->ops->setfromoptions    = PCSetFromOptions_LU;
324   pc->ops->view              = PCView_LU;
325   pc->ops->applyrichardson   = 0;
326   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
327 
328   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
329   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
330   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
331   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
332   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
333   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
334   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
335   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
336   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);CHKERRQ(ierr);
337   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
338   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
339   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr);
340   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
341   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
342   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
343   PetscFunctionReturn(0);
344 }
345