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