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