xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 6baea16979a9c8ee8132e9efba43c7ba974ff0b3)
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(PetscOptionItems *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   ierr = MatSetErrorIfFPE(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
103   if (dir->inplace) {
104     if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
105     ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
106     ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
107     if (dir->row) {
108       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
109       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
110     }
111     ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
112     if (((PC_Factor*)dir)->info.errortype) { /* Factor() fails */
113       MatFactorInfo factinfo=((PC_Factor*)dir)->info;
114       pc->failedreason = (PCFailedReason)factinfo.errortype;
115       PetscFunctionReturn(0);
116     }
117 
118     ((PC_Factor*)dir)->fact = pc->pmat;
119   } else {
120     MatInfo info;
121     if (!pc->setupcalled) {
122       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
123       if (dir->nonzerosalongdiagonal) {
124         ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
125       }
126       if (dir->row) {
127         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
128         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
129       }
130       if (!((PC_Factor*)dir)->fact) {
131         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
132       }
133       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
134       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
135       dir->actualfill = info.fill_ratio_needed;
136       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
137     } else if (pc->flag != SAME_NONZERO_PATTERN) {
138       if (!dir->reuseordering) {
139         if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
140         ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
141         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
142         if (dir->nonzerosalongdiagonal) {
143           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
144         }
145         if (dir->row) {
146           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
147           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
148         }
149       }
150       ierr            = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
151       ierr            = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
152       ierr            = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
153       ierr            = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
154       dir->actualfill = info.fill_ratio_needed;
155       ierr            = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
156     }
157     if (((PC_Factor*)dir)->info.errortype) { /* FactorSymbolic() fails */
158       MatFactorInfo factinfo=((PC_Factor*)dir)->info;
159       pc->failedreason = (PCFailedReason)factinfo.errortype;
160       PetscFunctionReturn(0);
161     }
162 
163     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
164     if (((PC_Factor*)dir)->info.errortype) { /* FactorNumeric() fails */
165       MatFactorInfo factinfo=((PC_Factor*)dir)->info;
166       pc->failedreason = (PCFailedReason)factinfo.errortype;
167     }
168   }
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "PCReset_LU"
174 static PetscErrorCode PCReset_LU(PC pc)
175 {
176   PC_LU          *dir = (PC_LU*)pc->data;
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
181   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
182   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PCDestroy_LU"
188 static PetscErrorCode PCDestroy_LU(PC pc)
189 {
190   PC_LU          *dir = (PC_LU*)pc->data;
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   ierr = PCReset_LU(pc);CHKERRQ(ierr);
195   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
196   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
197   ierr = PetscFree(pc->data);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PCApply_LU"
203 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
204 {
205   PC_LU          *dir = (PC_LU*)pc->data;
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209   if (dir->inplace) {
210     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
211   } else {
212     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "PCApplyTranspose_LU"
219 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
220 {
221   PC_LU          *dir = (PC_LU*)pc->data;
222   PetscErrorCode ierr;
223 
224   PetscFunctionBegin;
225   if (dir->inplace) {
226     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
227   } else {
228     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
229   }
230   PetscFunctionReturn(0);
231 }
232 
233 /* -----------------------------------------------------------------------------------*/
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "PCFactorSetUseInPlace_LU"
237 PetscErrorCode  PCFactorSetUseInPlace_LU(PC pc,PetscBool flg)
238 {
239   PC_LU *dir = (PC_LU*)pc->data;
240 
241   PetscFunctionBegin;
242   dir->inplace = flg;
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "PCFactorGetUseInPlace_LU"
248 PetscErrorCode  PCFactorGetUseInPlace_LU(PC pc,PetscBool *flg)
249 {
250   PC_LU *dir = (PC_LU*)pc->data;
251 
252   PetscFunctionBegin;
253   *flg = dir->inplace;
254   PetscFunctionReturn(0);
255 }
256 
257 /* ------------------------------------------------------------------------ */
258 
259 /*MC
260    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
261 
262    Options Database Keys:
263 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
264 .  -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu
265 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
266 .  -pc_factor_fill <fill> - Sets fill amount
267 .  -pc_factor_in_place - Activates in-place factorization
268 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
269 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
270                                          stability of factorization.
271 .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
272 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
273 -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
274         diagonal.
275 
276    Notes: Not all options work for all matrix formats
277           Run with -help to see additional options for particular matrix formats or factorization
278           algorithms
279 
280    Level: beginner
281 
282    Concepts: LU factorization, direct solver
283 
284    Notes: Usually this will compute an "exact" solution in one iteration and does
285           not need a Krylov method (i.e. you can use -ksp_type preonly, or
286           KSPSetType(ksp,KSPPREONLY) for the Krylov method
287 
288 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
289            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
290            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
291            PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
292            PCFactorReorderForNonzeroDiagonal()
293 M*/
294 
295 #undef __FUNCT__
296 #define __FUNCT__ "PCCreate_LU"
297 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
298 {
299   PetscErrorCode ierr;
300   PetscMPIInt    size;
301   PC_LU          *dir;
302 
303   PetscFunctionBegin;
304   ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr);
305 
306   ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr);
307 
308   ((PC_Factor*)dir)->fact       = NULL;
309   ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU;
310   dir->inplace                  = PETSC_FALSE;
311   dir->nonzerosalongdiagonal    = PETSC_FALSE;
312 
313   ((PC_Factor*)dir)->info.fill          = 5.0;
314   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
315   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
316   ((PC_Factor*)dir)->info.shiftamount   = 0.0;
317   ((PC_Factor*)dir)->info.zeropivot     = 100.0*PETSC_MACHINE_EPSILON;
318   ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
319   dir->col                              = 0;
320   dir->row                              = 0;
321 
322   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */
323   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
324   if (size == 1) {
325     ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
326   } else {
327     ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
328   }
329   dir->reusefill     = PETSC_FALSE;
330   dir->reuseordering = PETSC_FALSE;
331   pc->data           = (void*)dir;
332 
333   pc->ops->reset             = PCReset_LU;
334   pc->ops->destroy           = PCDestroy_LU;
335   pc->ops->apply             = PCApply_LU;
336   pc->ops->applytranspose    = PCApplyTranspose_LU;
337   pc->ops->setup             = PCSetUp_LU;
338   pc->ops->setfromoptions    = PCSetFromOptions_LU;
339   pc->ops->view              = PCView_LU;
340   pc->ops->applyrichardson   = 0;
341   pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
342 
343   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
344   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
345   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
346   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
347   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
348   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
349   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
350   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_LU);CHKERRQ(ierr);
351   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_LU);CHKERRQ(ierr);
352   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
353   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_LU);CHKERRQ(ierr);
354   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_LU);CHKERRQ(ierr);
355   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetColumnPivot_C",PCFactorSetColumnPivot_Factor);CHKERRQ(ierr);
356   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr);
357   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
358   PetscFunctionReturn(0);
359 }
360