xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision e6e75211d226c622f451867f53ce5d558649ff4f)
1 
2 #include <../src/ksp/pc/impls/factor/factor.h>     /*I "petscpc.h"  I*/
3 
4 /* ------------------------------------------------------------------------------------------*/
5 
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "PCFactorSetUpMatSolverPackage_Factor"
9 PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc)
10 {
11   PC_Factor      *icc = (PC_Factor*)pc->data;
12   PetscErrorCode ierr;
13 
14   PetscFunctionBegin;
15   if (!pc->setupcalled && !((PC_Factor*)icc)->fact) {
16     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
17   }
18   PetscFunctionReturn(0);
19 }
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "PCFactorSetZeroPivot_Factor"
23 PetscErrorCode  PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
24 {
25   PC_Factor *ilu = (PC_Factor*)pc->data;
26 
27   PetscFunctionBegin;
28   ilu->info.zeropivot = z;
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "PCFactorSetShiftType_Factor"
34 PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
35 {
36   PC_Factor *dir = (PC_Factor*)pc->data;
37 
38   PetscFunctionBegin;
39   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
40   else {
41     dir->info.shifttype = (PetscReal) shifttype;
42     if ((shifttype == MAT_SHIFT_NONZERO || shifttype ==  MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) {
43       dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */
44     }
45   }
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "PCFactorSetShiftAmount_Factor"
51 PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
52 {
53   PC_Factor *dir = (PC_Factor*)pc->data;
54 
55   PetscFunctionBegin;
56   if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
57   else dir->info.shiftamount = shiftamount;
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "PCFactorSetDropTolerance_Factor"
63 PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
64 {
65   PC_Factor *ilu = (PC_Factor*)pc->data;
66 
67   PetscFunctionBegin;
68   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
69     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
70   }
71   ilu->info.usedt   = PETSC_TRUE;
72   ilu->info.dt      = dt;
73   ilu->info.dtcol   = dtcol;
74   ilu->info.dtcount = dtcount;
75   ilu->info.fill    = PETSC_DEFAULT;
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "PCFactorSetFill_Factor"
81 PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
82 {
83   PC_Factor *dir = (PC_Factor*)pc->data;
84 
85   PetscFunctionBegin;
86   dir->info.fill = fill;
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "PCFactorSetMatOrderingType_Factor"
92 PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering)
93 {
94   PC_Factor      *dir = (PC_Factor*)pc->data;
95   PetscErrorCode ierr;
96   PetscBool      flg;
97 
98   PetscFunctionBegin;
99   if (!pc->setupcalled) {
100     ierr = PetscFree(dir->ordering);CHKERRQ(ierr);
101     ierr = PetscStrallocpy(ordering,(char**)&dir->ordering);CHKERRQ(ierr);
102   } else {
103     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
104     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
105   }
106   PetscFunctionReturn(0);
107 }
108 
109 #undef __FUNCT__
110 #define __FUNCT__ "PCFactorGetLevels_Factor"
111 PetscErrorCode  PCFactorGetLevels_Factor(PC pc,PetscInt *levels)
112 {
113   PC_Factor      *ilu = (PC_Factor*)pc->data;
114 
115   PetscFunctionBegin;
116   *levels = ilu->info.levels;
117   PetscFunctionReturn(0);
118 }
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "PCFactorSetLevels_Factor"
122 PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
123 {
124   PC_Factor      *ilu = (PC_Factor*)pc->data;
125   PetscErrorCode ierr;
126 
127   PetscFunctionBegin;
128   if (!pc->setupcalled) ilu->info.levels = levels;
129   else if (ilu->info.levels != levels) {
130     ierr             = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */
131     pc->setupcalled  = 0; /* force a complete rebuild of preconditioner factored matrices */
132     ilu->info.levels = levels;
133   } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt");
134   PetscFunctionReturn(0);
135 }
136 
137 #undef __FUNCT__
138 #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
139 PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg)
140 {
141   PC_Factor *dir = (PC_Factor*)pc->data;
142 
143   PetscFunctionBegin;
144   dir->info.diagonal_fill = (PetscReal) flg;
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "PCFactorGetAllowDiagonalFill_Factor"
150 PetscErrorCode  PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg)
151 {
152   PC_Factor *dir = (PC_Factor*)pc->data;
153 
154   PetscFunctionBegin;
155   *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE;
156   PetscFunctionReturn(0);
157 }
158 
159 /* ------------------------------------------------------------------------------------------*/
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
163 PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
164 {
165   PC_Factor *dir = (PC_Factor*)pc->data;
166 
167   PetscFunctionBegin;
168   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "PCFactorGetMatrix_Factor"
174 PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
175 {
176   PC_Factor *ilu = (PC_Factor*)pc->data;
177 
178   PetscFunctionBegin;
179   if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
180   *mat = ilu->fact;
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
186 PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
187 {
188   PetscErrorCode ierr;
189   PC_Factor      *lu = (PC_Factor*)pc->data;
190 
191   PetscFunctionBegin;
192   if (lu->fact) {
193     const MatSolverPackage ltype;
194     PetscBool              flg;
195     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
196     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
197     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
198   } else {
199     ierr = PetscFree(lu->solvertype);CHKERRQ(ierr);
200     ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
201   }
202   PetscFunctionReturn(0);
203 }
204 
205 #undef __FUNCT__
206 #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
207 PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
208 {
209   PC_Factor *lu = (PC_Factor*)pc->data;
210 
211   PetscFunctionBegin;
212   *stype = lu->solvertype;
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "PCFactorSetColumnPivot_Factor"
218 PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
219 {
220   PC_Factor *dir = (PC_Factor*)pc->data;
221 
222   PetscFunctionBegin;
223   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",(double)dtcol);
224   dir->info.dtcol = dtcol;
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "PCSetFromOptions_Factor"
230 PetscErrorCode  PCSetFromOptions_Factor(PetscOptions *PetscOptionsObject,PC pc)
231 {
232   PC_Factor         *factor = (PC_Factor*)pc->data;
233   PetscErrorCode    ierr;
234   PetscBool         flg,set;
235   char              tname[256], solvertype[64];
236   PetscFunctionList ordlist;
237   PetscEnum         etmp;
238   PetscBool         inplace;
239 
240   PetscFunctionBegin;
241   ierr = PCFactorGetUseInPlace(pc,&inplace);CHKERRQ(ierr);
242   ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set);CHKERRQ(ierr);
243   if (set) {
244     ierr = PCFactorSetUseInPlace(pc,flg);CHKERRQ(ierr);
245   }
246   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,NULL);CHKERRQ(ierr);
247 
248   ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr);
249   if (flg) {
250     ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr);
251   }
252   ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);CHKERRQ(ierr);
253 
254   ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);CHKERRQ(ierr);
255   ierr = PetscOptionsReal("-pc_factor_column_pivot","Column pivot tolerance (used only for some factorization)","PCFactorSetColumnPivot",((PC_Factor*)factor)->info.dtcol,&((PC_Factor*)factor)->info.dtcol,&flg);CHKERRQ(ierr);
256 
257   ierr = PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
258   if (set) {
259     ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
260   }
261 
262   ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
263   if (set) {
264     ierr = PCFactorSetReuseFill(pc,flg);CHKERRQ(ierr);
265   }
266   ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
267   if (set) {
268     ierr = PCFactorSetReuseOrdering(pc,flg);CHKERRQ(ierr);
269   }
270 
271   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
272   ierr = PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr);
273   if (flg) {
274     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
275   }
276 
277   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
278   ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
279   if (flg) {
280     ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
281   }
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "PCView_Factor"
287 PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
288 {
289   PC_Factor      *factor = (PC_Factor*)pc->data;
290   PetscErrorCode ierr;
291   PetscBool      isstring,iascii;
292 
293   PetscFunctionBegin;
294   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
295   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
296   if (iascii) {
297     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
298       if (factor->info.dt > 0) {
299         ierr = PetscViewerASCIIPrintf(viewer,"  drop tolerance %g\n",(double)factor->info.dt);CHKERRQ(ierr);
300         ierr = PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr);
301         ierr = PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %g\n",(double)factor->info.dtcol);CHKERRQ(ierr);
302       } else if (factor->info.levels == 1) {
303         ierr = PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
304       } else {
305         ierr = PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
306       }
307     }
308 
309     ierr = PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %g\n",(double)factor->info.zeropivot);CHKERRQ(ierr);
310     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
311       ierr = PetscViewerASCIIPrintf(viewer,"  using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);CHKERRQ(ierr);
312     }
313 
314     ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr);
315 
316     if (factor->fact) {
317       MatInfo info;
318       ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
319       ierr = PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);CHKERRQ(ierr);
320       ierr = PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");CHKERRQ(ierr);
321       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
322       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
323       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
324       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
325       ierr = MatView(factor->fact,viewer);CHKERRQ(ierr);
326       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
327       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
328       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
329       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
330     }
331 
332   } else if (isstring) {
333     MatFactorType t;
334     ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr);
335     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
336       ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
337     }
338   }
339   PetscFunctionReturn(0);
340 }
341