Actual source code: factimpl.c
petsc-3.4.2 2013-07-02
2: #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/
4: /* ------------------------------------------------------------------------------------------*/
9: PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc)
10: {
11: PC_Factor *icc = (PC_Factor*)pc->data;
15: if (!pc->setupcalled && !((PC_Factor*)icc)->fact) {
16: MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);
17: }
18: return(0);
19: }
23: PetscErrorCode PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
24: {
25: PC_Factor *ilu = (PC_Factor*)pc->data;
28: ilu->info.zeropivot = z;
29: return(0);
30: }
34: PetscErrorCode PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
35: {
36: PC_Factor *dir = (PC_Factor*)pc->data;
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: return(0);
47: }
51: PetscErrorCode PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
52: {
53: PC_Factor *dir = (PC_Factor*)pc->data;
56: if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
57: else dir->info.shiftamount = shiftamount;
58: return(0);
59: }
63: PetscErrorCode PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
64: {
65: PC_Factor *ilu = (PC_Factor*)pc->data;
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: return(0);
77: }
81: PetscErrorCode PCFactorSetFill_Factor(PC pc,PetscReal fill)
82: {
83: PC_Factor *dir = (PC_Factor*)pc->data;
86: dir->info.fill = fill;
87: return(0);
88: }
92: PetscErrorCode PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering)
93: {
94: PC_Factor *dir = (PC_Factor*)pc->data;
96: PetscBool flg;
99: if (!pc->setupcalled) {
100: PetscFree(dir->ordering);
101: PetscStrallocpy(ordering,(char**)&dir->ordering);
102: } else {
103: PetscStrcmp(dir->ordering,ordering,&flg);
104: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
105: }
106: return(0);
107: }
111: PetscErrorCode PCFactorSetLevels_Factor(PC pc,PetscInt levels)
112: {
113: PC_Factor *ilu = (PC_Factor*)pc->data;
117: if (!pc->setupcalled) ilu->info.levels = levels;
118: else if (ilu->info.levels != levels) {
119: (*pc->ops->reset)(pc); /* 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: return(0);
124: }
128: PetscErrorCode PCFactorSetAllowDiagonalFill_Factor(PC pc)
129: {
130: PC_Factor *dir = (PC_Factor*)pc->data;
133: dir->info.diagonal_fill = 1;
134: return(0);
135: }
137: /* ------------------------------------------------------------------------------------------*/
141: PetscErrorCode PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
142: {
143: PC_Factor *dir = (PC_Factor*)pc->data;
146: dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
147: return(0);
148: }
152: PetscErrorCode PCFactorGetMatrix_Factor(PC pc,Mat *mat)
153: {
154: PC_Factor *ilu = (PC_Factor*)pc->data;
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: return(0);
160: }
164: PetscErrorCode PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
165: {
167: PC_Factor *lu = (PC_Factor*)pc->data;
170: if (lu->fact) {
171: const MatSolverPackage ltype;
172: PetscBool flg;
173: MatFactorGetSolverPackage(lu->fact,<ype);
174: PetscStrcmp(stype,ltype,&flg);
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: PetscFree(lu->solvertype);
178: PetscStrallocpy(stype,&lu->solvertype);
179: }
180: return(0);
181: }
185: PetscErrorCode PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
186: {
187: PC_Factor *lu = (PC_Factor*)pc->data;
190: *stype = lu->solvertype;
191: return(0);
192: }
196: PetscErrorCode PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
197: {
198: PC_Factor *dir = (PC_Factor*)pc->data;
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: return(0);
204: }
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;
218: if (!MatOrderingRegisterAllCalled) {MatOrderingRegisterAll();}
219: PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,NULL);
220: if (flg) {
221: PCFactorSetUseInPlace(pc);
222: }
223: PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);
225: PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",
226: MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);
227: if (flg) {
228: PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);
229: }
230: PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);
232: PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);
233: 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);
235: flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
236: PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);
237: if (set) {
238: PCFactorSetPivotInBlocks(pc,flg);
239: }
241: flg = PETSC_FALSE;
242: PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,NULL);
243: if (flg) {
244: PCFactorSetReuseFill(pc,PETSC_TRUE);
245: }
246: flg = PETSC_FALSE;
247: PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,NULL);
248: if (flg) {
249: PCFactorSetReuseOrdering(pc,PETSC_TRUE);
250: }
252: MatGetOrderingList(&ordlist);
253: PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);
254: if (flg) {
255: PCFactorSetMatOrderingType(pc,tname);
256: }
258: /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
259: PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);
260: if (flg) {
261: PCFactorSetMatSolverPackage(pc,solvertype);
262: }
263: return(0);
264: }
268: PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
269: {
270: PC_Factor *factor = (PC_Factor*)pc->data;
272: PetscBool isstring,iascii;
275: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
276: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
277: if (iascii) {
278: if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
279: if (factor->info.dt > 0) {
280: PetscViewerASCIIPrintf(viewer," drop tolerance %G\n",factor->info.dt);
281: PetscViewerASCIIPrintf(viewer," max nonzeros per row %D\n",factor->info.dtcount);
282: PetscViewerASCIIPrintf(viewer," column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);
283: } else if (factor->info.levels == 1) {
284: PetscViewerASCIIPrintf(viewer," %D level of fill\n",(PetscInt)factor->info.levels);
285: } else {
286: PetscViewerASCIIPrintf(viewer," %D levels of fill\n",(PetscInt)factor->info.levels);
287: }
288: }
290: PetscViewerASCIIPrintf(viewer," tolerance for zero pivot %G\n",factor->info.zeropivot);
291: if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
292: PetscViewerASCIIPrintf(viewer," using Manteuffel shift\n");
293: }
294: if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) {
295: PetscViewerASCIIPrintf(viewer," using diagonal shift to prevent zero pivot\n");
296: }
297: if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) {
298: PetscViewerASCIIPrintf(viewer," using diagonal shift on blocks to prevent zero pivot\n");
299: }
301: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",factor->ordering);
303: if (factor->fact) {
304: MatInfo info;
305: MatGetInfo(factor->fact,MAT_LOCAL,&info);
306: PetscViewerASCIIPrintf(viewer," factor fill ratio given %G, needed %G\n",info.fill_ratio_given,info.fill_ratio_needed);
307: PetscViewerASCIIPrintf(viewer," Factored matrix follows:\n");
308: PetscViewerASCIIPushTab(viewer);
309: PetscViewerASCIIPushTab(viewer);
310: PetscViewerASCIIPushTab(viewer);
311: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
312: MatView(factor->fact,viewer);
313: PetscViewerPopFormat(viewer);
314: PetscViewerASCIIPopTab(viewer);
315: PetscViewerASCIIPopTab(viewer);
316: PetscViewerASCIIPopTab(viewer);
317: }
319: } else if (isstring) {
320: MatFactorType t;
321: MatGetFactorType(factor->fact,&t);
322: if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
323: PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);
324: }
325: }
326: return(0);
327: }