xref: /petsc/src/ksp/pc/impls/factor/lu/lu.c (revision 7d5fd1e4d9337468ad3f05b65b7facdcd2dfd2a4)
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 PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z)
11 {
12   PC_LU *lu = (PC_LU*)pc->data;
13 
14   PetscFunctionBegin;
15   lu->nonzerosalongdiagonal = PETSC_TRUE;
16   if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10;
17   else lu->nonzerosalongdiagonaltol = z;
18   PetscFunctionReturn(0);
19 }
20 
21 static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc)
22 {
23   PC_LU          *lu = (PC_LU*)pc->data;
24   PetscErrorCode ierr;
25   PetscBool      flg = PETSC_FALSE;
26   PetscReal      tol;
27 
28   PetscFunctionBegin;
29   ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr);
30   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
31 
32   ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
33   if (flg) {
34     tol  = PETSC_DECIDE;
35     ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL);CHKERRQ(ierr);
36     ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
37   }
38   ierr = PetscOptionsTail();CHKERRQ(ierr);
39   PetscFunctionReturn(0);
40 }
41 
42 static PetscErrorCode PCSetUp_LU(PC pc)
43 {
44   PetscErrorCode         ierr;
45   PC_LU                  *dir = (PC_LU*)pc->data;
46   MatSolverType          stype;
47   MatFactorError         err;
48 
49   PetscFunctionBegin;
50   pc->failedreason = PC_NOERROR;
51   if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill;
52 
53   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
54   if (dir->hdr.inplace) {
55     MatFactorType ftype;
56 
57     ierr = MatGetFactorType(pc->pmat, &ftype);CHKERRQ(ierr);
58     if (ftype == MAT_FACTOR_NONE) {
59       if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
60       ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
61       /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */
62       ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr);
63       ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
64       if (dir->row) {
65         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
66         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
67       }
68       ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
69       ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr);
70       if (err) { /* Factor() fails */
71         pc->failedreason = (PCFailedReason)err;
72         PetscFunctionReturn(0);
73       }
74     }
75     ((PC_Factor*)dir)->fact = pc->pmat;
76   } else {
77     MatInfo info;
78 
79     if (!pc->setupcalled) {
80       PetscBool canuseordering;
81       if (!((PC_Factor*)dir)->fact) {
82         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
83         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
84       }
85       ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr);
86       if (canuseordering) {
87         ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr);
88         ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
89         if (dir->nonzerosalongdiagonal) {
90           ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
91         }
92         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
93         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
94       }
95       ierr                = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
96       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
97       dir->hdr.actualfill = info.fill_ratio_needed;
98     } else if (pc->flag != SAME_NONZERO_PATTERN) {
99       PetscBool canuseordering;
100       if (!dir->hdr.reuseordering) {
101         ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
102         ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr);
103         ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr);
104         ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr);
105         if (canuseordering) {
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 = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr);
109           ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr);
110           if (dir->nonzerosalongdiagonal) {
111             ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr);
112           }
113           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr);
114           ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr);
115         }
116       }
117       ierr                = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
118       ierr                = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
119       dir->hdr.actualfill = info.fill_ratio_needed;
120     } else {
121       ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
122       if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) {
123         ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr);
124         pc->failedreason = PC_NOERROR;
125       }
126     }
127     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
128     if (err) { /* FactorSymbolic() fails */
129       pc->failedreason = (PCFailedReason)err;
130       PetscFunctionReturn(0);
131     }
132 
133     ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr);
134     ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr);
135     if (err) { /* FactorNumeric() fails */
136       pc->failedreason = (PCFailedReason)err;
137     }
138 
139   }
140 
141   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
142   if (!stype) {
143     MatSolverType solverpackage;
144     ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr);
145     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 static PetscErrorCode PCReset_LU(PC pc)
151 {
152   PC_LU          *dir = (PC_LU*)pc->data;
153   PetscErrorCode ierr;
154 
155   PetscFunctionBegin;
156   if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);}
157   if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);}
158   ierr = ISDestroy(&dir->col);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 static PetscErrorCode PCDestroy_LU(PC pc)
163 {
164   PC_LU          *dir = (PC_LU*)pc->data;
165   PetscErrorCode ierr;
166 
167   PetscFunctionBegin;
168   ierr = PCReset_LU(pc);CHKERRQ(ierr);
169   ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
170   ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr);
171   ierr = PetscFree(pc->data);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
176 {
177   PC_LU          *dir = (PC_LU*)pc->data;
178   PetscErrorCode ierr;
179 
180   PetscFunctionBegin;
181   if (dir->hdr.inplace) {
182     ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);
183   } else {
184     ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
185   }
186   PetscFunctionReturn(0);
187 }
188 
189 static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y)
190 {
191   PC_LU          *dir = (PC_LU*)pc->data;
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   if (dir->hdr.inplace) {
196     ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr);
197   } else {
198     ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr);
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
204 {
205   PC_LU          *dir = (PC_LU*)pc->data;
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209   if (dir->hdr.inplace) {
210     ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);
211   } else {
212     ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 /* -----------------------------------------------------------------------------------*/
218 
219 /*MC
220    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
221 
222    Options Database Keys:
223 +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
224 .  -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu
225 .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
226 .  -pc_factor_fill <fill> - Sets fill amount
227 .  -pc_factor_in_place - Activates in-place factorization
228 .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
229 .  -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
230                                          stability of factorization.
231 .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
232 .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
233 -   -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the
234         diagonal.
235 
236    Notes:
237     Not all options work for all matrix formats
238           Run with -help to see additional options for particular matrix formats or factorization
239           algorithms
240 
241    Level: beginner
242 
243    Notes:
244     Usually this will compute an "exact" solution in one iteration and does
245           not need a Krylov method (i.e. you can use -ksp_type preonly, or
246           KSPSetType(ksp,KSPPREONLY) for the Krylov method
247 
248 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
249            PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
250            PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(),
251            PCFactorSetPivotInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount()
252            PCFactorReorderForNonzeroDiagonal()
253 M*/
254 
255 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc)
256 {
257   PetscErrorCode ierr;
258   PC_LU          *dir;
259 
260   PetscFunctionBegin;
261   ierr     = PetscNewLog(pc,&dir);CHKERRQ(ierr);
262   pc->data = (void*)dir;
263   ierr     = PCFactorInitialize(pc,MAT_FACTOR_LU);CHKERRQ(ierr);
264   dir->nonzerosalongdiagonal    = PETSC_FALSE;
265 
266   ((PC_Factor*)dir)->info.fill          = 5.0;
267   ((PC_Factor*)dir)->info.dtcol         = 1.e-6;  /* default to pivoting; this is only thing PETSc LU supports */
268   ((PC_Factor*)dir)->info.shifttype     = (PetscReal)MAT_SHIFT_NONE;
269   dir->col                              = NULL;
270   dir->row                              = NULL;
271 
272   pc->ops->reset             = PCReset_LU;
273   pc->ops->destroy           = PCDestroy_LU;
274   pc->ops->apply             = PCApply_LU;
275   pc->ops->matapply          = PCMatApply_LU;
276   pc->ops->applytranspose    = PCApplyTranspose_LU;
277   pc->ops->setup             = PCSetUp_LU;
278   pc->ops->setfromoptions    = PCSetFromOptions_LU;
279   pc->ops->view              = PCView_Factor;
280   pc->ops->applyrichardson   = NULL;
281   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284