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