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