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