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