Actual source code: factimpl.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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  PCFactorGetLevels_Factor(PC pc,PetscInt *levels)
112: {
113:   PC_Factor      *ilu = (PC_Factor*)pc->data;

116:   *levels = ilu->info.levels;
117:   return(0);
118: }

122: PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
123: {
124:   PC_Factor      *ilu = (PC_Factor*)pc->data;

128:   if (!pc->setupcalled) ilu->info.levels = levels;
129:   else if (ilu->info.levels != levels) {
130:     (*pc->ops->reset)(pc); /* 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:   return(0);
135: }

139: PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc)
140: {
141:   PC_Factor *dir = (PC_Factor*)pc->data;

144:   dir->info.diagonal_fill = 1;
145:   return(0);
146: }

148: /* ------------------------------------------------------------------------------------------*/

152: PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
153: {
154:   PC_Factor *dir = (PC_Factor*)pc->data;

157:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
158:   return(0);
159: }

163: PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
164: {
165:   PC_Factor *ilu = (PC_Factor*)pc->data;

168:   if (!ilu->fact) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
169:   *mat = ilu->fact;
170:   return(0);
171: }

175: PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
176: {
178:   PC_Factor      *lu = (PC_Factor*)pc->data;

181:   if (lu->fact) {
182:     const MatSolverPackage ltype;
183:     PetscBool              flg;
184:     MatFactorGetSolverPackage(lu->fact,&ltype);
185:     PetscStrcmp(stype,ltype,&flg);
186:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
187:   } else {
188:     PetscFree(lu->solvertype);
189:     PetscStrallocpy(stype,&lu->solvertype);
190:   }
191:   return(0);
192: }

196: PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
197: {
198:   PC_Factor *lu = (PC_Factor*)pc->data;

201:   *stype = lu->solvertype;
202:   return(0);
203: }

207: PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
208: {
209:   PC_Factor *dir = (PC_Factor*)pc->data;

212:   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);
213:   dir->info.dtcol = dtcol;
214:   return(0);
215: }

219: PetscErrorCode  PCSetFromOptions_Factor(PC pc)
220: {
221:   PC_Factor         *factor = (PC_Factor*)pc->data;
222:   PetscErrorCode    ierr;
223:   PetscBool         flg = PETSC_FALSE,set;
224:   char              tname[256], solvertype[64];
225:   PetscFunctionList ordlist;
226:   PetscEnum         etmp;

229:   if (!MatOrderingRegisterAllCalled) {MatOrderingRegisterAll();}
230:   PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,NULL);
231:   if (flg) {
232:     PCFactorSetUseInPlace(pc);
233:   }
234:   PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);

236:   PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",
237:                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);
238:   if (flg) {
239:     PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);
240:   }
241:   PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);

243:   PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);
244:   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);

246:   flg  = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
247:   PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);
248:   if (set) {
249:     PCFactorSetPivotInBlocks(pc,flg);
250:   }

252:   flg  = PETSC_FALSE;
253:   PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,NULL);
254:   if (flg) {
255:     PCFactorSetReuseFill(pc,PETSC_TRUE);
256:   }
257:   flg  = PETSC_FALSE;
258:   PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,NULL);
259:   if (flg) {
260:     PCFactorSetReuseOrdering(pc,PETSC_TRUE);
261:   }

263:   MatGetOrderingList(&ordlist);
264:   PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);
265:   if (flg) {
266:     PCFactorSetMatOrderingType(pc,tname);
267:   }

269:   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
270:   PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);
271:   if (flg) {
272:     PCFactorSetMatSolverPackage(pc,solvertype);
273:   }
274:   return(0);
275: }

279: PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
280: {
281:   PC_Factor      *factor = (PC_Factor*)pc->data;
283:   PetscBool      isstring,iascii;

286:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
287:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
288:   if (iascii) {
289:     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
290:       if (factor->info.dt > 0) {
291:         PetscViewerASCIIPrintf(viewer,"  drop tolerance %g\n",(double)factor->info.dt);
292:         PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);
293:         PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %g\n",(double)factor->info.dtcol);
294:       } else if (factor->info.levels == 1) {
295:         PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);
296:       } else {
297:         PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);
298:       }
299:     }

301:     PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %g\n",(double)factor->info.zeropivot);
302:     if (MatFactorShiftTypesDetail[(int)factor->info.shifttype]) { /* Only print when using a nontrivial shift */
303:       PetscViewerASCIIPrintf(viewer,"  using %s [%s]\n",MatFactorShiftTypesDetail[(int)factor->info.shifttype],MatFactorShiftTypes[(int)factor->info.shifttype]);
304:     }

306:     PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);

308:     if (factor->fact) {
309:       MatInfo info;
310:       MatGetInfo(factor->fact,MAT_LOCAL,&info);
311:       PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);
312:       PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");
313:       PetscViewerASCIIPushTab(viewer);
314:       PetscViewerASCIIPushTab(viewer);
315:       PetscViewerASCIIPushTab(viewer);
316:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
317:       MatView(factor->fact,viewer);
318:       PetscViewerPopFormat(viewer);
319:       PetscViewerASCIIPopTab(viewer);
320:       PetscViewerASCIIPopTab(viewer);
321:       PetscViewerASCIIPopTab(viewer);
322:     }

324:   } else if (isstring) {
325:     MatFactorType t;
326:     MatGetFactorType(factor->fact,&t);
327:     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
328:       PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);
329:     }
330:   }
331:   return(0);
332: }