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