xref: /petsc/src/ksp/pc/impls/factor/factimpl.c (revision 0005d66cf2be2b0fca19463d69605dfad837b879)
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__ "PCFactorSetLevels_Factor"
111 PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
112 {
113   PC_Factor      *ilu = (PC_Factor*)pc->data;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   if (!pc->setupcalled) ilu->info.levels = levels;
118   else if (ilu->info.levels != levels) {
119     ierr             = (*pc->ops->reset)(pc);CHKERRQ(ierr); /* remove previous factored matrices */
120     pc->setupcalled  = 0; /* force a complete rebuild of preconditioner factored matrices */
121     ilu->info.levels = levels;
122   } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt");
123   PetscFunctionReturn(0);
124 }
125 
126 #undef __FUNCT__
127 #define __FUNCT__ "PCFactorSetAllowDiagonalFill_Factor"
128 PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc)
129 {
130   PC_Factor *dir = (PC_Factor*)pc->data;
131 
132   PetscFunctionBegin;
133   dir->info.diagonal_fill = 1;
134   PetscFunctionReturn(0);
135 }
136 
137 /* ------------------------------------------------------------------------------------------*/
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "PCFactorSetPivotInBlocks_Factor"
141 PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
142 {
143   PC_Factor *dir = (PC_Factor*)pc->data;
144 
145   PetscFunctionBegin;
146   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "PCFactorGetMatrix_Factor"
152 PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
153 {
154   PC_Factor *ilu = (PC_Factor*)pc->data;
155 
156   PetscFunctionBegin;
157   if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
158   *mat = ilu->fact;
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "PCFactorSetMatSolverPackage_Factor"
164 PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
165 {
166   PetscErrorCode ierr;
167   PC_Factor      *lu = (PC_Factor*)pc->data;
168 
169   PetscFunctionBegin;
170   if (lu->fact) {
171     const MatSolverPackage ltype;
172     PetscBool              flg;
173     ierr = MatFactorGetSolverPackage(lu->fact,&ltype);CHKERRQ(ierr);
174     ierr = PetscStrcmp(stype,ltype,&flg);CHKERRQ(ierr);
175     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
176   } else {
177     ierr = PetscFree(lu->solvertype);CHKERRQ(ierr);
178     ierr = PetscStrallocpy(stype,&lu->solvertype);CHKERRQ(ierr);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "PCFactorGetMatSolverPackage_Factor"
185 PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
186 {
187   PC_Factor *lu = (PC_Factor*)pc->data;
188 
189   PetscFunctionBegin;
190   *stype = lu->solvertype;
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "PCFactorSetColumnPivot_Factor"
196 PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
197 {
198   PC_Factor *dir = (PC_Factor*)pc->data;
199 
200   PetscFunctionBegin;
201   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",dtcol);
202   dir->info.dtcol = dtcol;
203   PetscFunctionReturn(0);
204 }
205 
206 #undef __FUNCT__
207 #define __FUNCT__ "PCSetFromOptions_Factor"
208 PetscErrorCode  PCSetFromOptions_Factor(PC pc)
209 {
210   PC_Factor         *factor = (PC_Factor*)pc->data;
211   PetscErrorCode    ierr;
212   PetscBool         flg = PETSC_FALSE,set;
213   char              tname[256], solvertype[64];
214   PetscFunctionList ordlist;
215   PetscEnum         etmp;
216 
217   PetscFunctionBegin;
218   if (!MatOrderingRegisterAllCalled) {ierr = MatOrderingRegisterAll();CHKERRQ(ierr);}
219   ierr = PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,NULL);CHKERRQ(ierr);
220   if (flg) {
221     ierr = PCFactorSetUseInPlace(pc);CHKERRQ(ierr);
222   }
223   ierr = PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);CHKERRQ(ierr);
224 
225   ierr = PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",
226                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);CHKERRQ(ierr);
227   if (flg) {
228     ierr = PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);CHKERRQ(ierr);
229   }
230   ierr = PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);CHKERRQ(ierr);
231 
232   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);
233   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);
234 
235   flg  = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
236   ierr = PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
237   if (set) {
238     ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
239   }
240 
241   flg  = PETSC_FALSE;
242   ierr = PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,NULL);CHKERRQ(ierr);
243   if (flg) {
244     ierr = PCFactorSetReuseFill(pc,PETSC_TRUE);CHKERRQ(ierr);
245   }
246   flg  = PETSC_FALSE;
247   ierr = PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,NULL);CHKERRQ(ierr);
248   if (flg) {
249     ierr = PCFactorSetReuseOrdering(pc,PETSC_TRUE);CHKERRQ(ierr);
250   }
251 
252   ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
253   ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);CHKERRQ(ierr);
254   if (flg) {
255     ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
256   }
257 
258   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
259   ierr = PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);CHKERRQ(ierr);
260   if (flg) {
261     ierr = PCFactorSetMatSolverPackage(pc,solvertype);CHKERRQ(ierr);
262   }
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "PCView_Factor"
268 PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
269 {
270   PC_Factor      *factor = (PC_Factor*)pc->data;
271   PetscErrorCode ierr;
272   PetscBool      isstring,iascii;
273 
274   PetscFunctionBegin;
275   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
276   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
277   if (iascii) {
278     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
279       if (factor->info.dt > 0) {
280         ierr = PetscViewerASCIIPrintf(viewer,"  drop tolerance %G\n",factor->info.dt);CHKERRQ(ierr);
281         ierr = PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);CHKERRQ(ierr);
282         ierr = PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);CHKERRQ(ierr);
283       } else if (factor->info.levels == 1) {
284         ierr = PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
285       } else {
286         ierr = PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);CHKERRQ(ierr);
287       }
288     }
289 
290     ierr = PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %G\n",factor->info.zeropivot);CHKERRQ(ierr);
291     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
292       ierr = PetscViewerASCIIPrintf(viewer,"  using Manteuffel shift\n");CHKERRQ(ierr);
293     }
294     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) {
295       ierr = PetscViewerASCIIPrintf(viewer,"  using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);
296     }
297     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) {
298       ierr = PetscViewerASCIIPrintf(viewer,"  using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);
299     }
300 
301     ierr = PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);CHKERRQ(ierr);
302 
303     if (factor->fact) {
304       MatInfo info;
305       ierr = MatGetInfo(factor->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
306       ierr = PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %G, needed %G\n",info.fill_ratio_given,info.fill_ratio_needed);CHKERRQ(ierr);
307       ierr = PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");CHKERRQ(ierr);
308       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
309       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
310       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
311       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
312       ierr = MatView(factor->fact,viewer);CHKERRQ(ierr);
313       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
314       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
315       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
316       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
317     }
318 
319   } else if (isstring) {
320     MatFactorType t;
321     ierr = MatGetFactorType(factor->fact,&t);CHKERRQ(ierr);
322     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
323       ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
324     }
325   }
326   PetscFunctionReturn(0);
327 }
328