Actual source code: factimpl.c

petsc-3.7.5 2017-01-01
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->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"You can only call this routine after the matrix object has been provided to the solver, for example with KSPSetOperators() or SNESSetJacobian()");
 16:   if (!pc->setupcalled && !((PC_Factor*)icc)->fact) {
 17:     MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,((PC_Factor*)icc)->factortype,&((PC_Factor*)icc)->fact);
 18:   }
 19:   return(0);
 20: }

 24: PetscErrorCode  PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
 25: {
 26:   PC_Factor *ilu = (PC_Factor*)pc->data;

 29:   ilu->info.zeropivot = z;
 30:   return(0);
 31: }

 35: PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
 36: {
 37:   PC_Factor *dir = (PC_Factor*)pc->data;

 40:   if (shifttype == (MatFactorShiftType)PETSC_DECIDE) dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
 41:   else {
 42:     dir->info.shifttype = (PetscReal) shifttype;
 43:     if ((shifttype == MAT_SHIFT_NONZERO || shifttype ==  MAT_SHIFT_INBLOCKS) && dir->info.shiftamount == 0.0) {
 44:       dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; /* set default amount if user has not called PCFactorSetShiftAmount() yet */
 45:     }
 46:   }
 47:   return(0);
 48: }

 52: PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
 53: {
 54:   PC_Factor *dir = (PC_Factor*)pc->data;

 57:   if (shiftamount == (PetscReal) PETSC_DECIDE) dir->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
 58:   else dir->info.shiftamount = shiftamount;
 59:   return(0);
 60: }

 64: PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 65: {
 66:   PC_Factor *ilu = (PC_Factor*)pc->data;

 69:   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 70:     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
 71:   }
 72:   ilu->info.usedt   = PETSC_TRUE;
 73:   ilu->info.dt      = dt;
 74:   ilu->info.dtcol   = dtcol;
 75:   ilu->info.dtcount = dtcount;
 76:   ilu->info.fill    = PETSC_DEFAULT;
 77:   return(0);
 78: }

 82: PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
 83: {
 84:   PC_Factor *dir = (PC_Factor*)pc->data;

 87:   dir->info.fill = fill;
 88:   return(0);
 89: }

 93: PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,MatOrderingType ordering)
 94: {
 95:   PC_Factor      *dir = (PC_Factor*)pc->data;
 97:   PetscBool      flg;

100:   if (!pc->setupcalled) {
101:     PetscFree(dir->ordering);
102:     PetscStrallocpy(ordering,(char**)&dir->ordering);
103:   } else {
104:     PetscStrcmp(dir->ordering,ordering,&flg);
105:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
106:   }
107:   return(0);
108: }

112: PetscErrorCode  PCFactorGetLevels_Factor(PC pc,PetscInt *levels)
113: {
114:   PC_Factor      *ilu = (PC_Factor*)pc->data;

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

123: PetscErrorCode  PCFactorGetZeroPivot_Factor(PC pc,PetscReal *pivot)
124: {
125:   PC_Factor      *ilu = (PC_Factor*)pc->data;

128:   *pivot = ilu->info.zeropivot;
129:   return(0);
130: }

134: PetscErrorCode  PCFactorGetShiftAmount_Factor(PC pc,PetscReal *shift)
135: {
136:   PC_Factor      *ilu = (PC_Factor*)pc->data;

139:   *shift = ilu->info.shiftamount;
140:   return(0);
141: }

145: PetscErrorCode  PCFactorGetShiftType_Factor(PC pc,MatFactorShiftType *type)
146: {
147:   PC_Factor      *ilu = (PC_Factor*)pc->data;

150:   *type = (MatFactorShiftType) (int) ilu->info.shifttype;
151:   return(0);
152: }

156: PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
157: {
158:   PC_Factor      *ilu = (PC_Factor*)pc->data;

162:   if (!pc->setupcalled) ilu->info.levels = levels;
163:   else if (ilu->info.levels != levels) {
164:     (*pc->ops->reset)(pc); /* remove previous factored matrices */
165:     pc->setupcalled  = 0; /* force a complete rebuild of preconditioner factored matrices */
166:     ilu->info.levels = levels;
167:   } else if (ilu->info.usedt) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use with ILUdt");
168:   return(0);
169: }

173: PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc,PetscBool flg)
174: {
175:   PC_Factor *dir = (PC_Factor*)pc->data;

178:   dir->info.diagonal_fill = (PetscReal) flg;
179:   return(0);
180: }

184: PetscErrorCode  PCFactorGetAllowDiagonalFill_Factor(PC pc,PetscBool *flg)
185: {
186:   PC_Factor *dir = (PC_Factor*)pc->data;

189:   *flg = dir->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE;
190:   return(0);
191: }

193: /* ------------------------------------------------------------------------------------------*/

197: PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool pivot)
198: {
199:   PC_Factor *dir = (PC_Factor*)pc->data;

202:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
203:   return(0);
204: }

208: PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
209: {
210:   PC_Factor *ilu = (PC_Factor*)pc->data;

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

220: PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
221: {
223:   PC_Factor      *lu = (PC_Factor*)pc->data;

226:   if (lu->fact) {
227:     const MatSolverPackage ltype;
228:     PetscBool              flg;
229:     MatFactorGetSolverPackage(lu->fact,&ltype);
230:     PetscStrcmp(stype,ltype,&flg);
231:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
232:   }
233: 
234:   PetscFree(lu->solvertype);
235:   PetscStrallocpy(stype,&lu->solvertype);
236:   return(0);
237: }

241: PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
242: {
243:   PC_Factor *lu = (PC_Factor*)pc->data;

246:   *stype = lu->solvertype;
247:   return(0);
248: }

252: PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
253: {
254:   PC_Factor *dir = (PC_Factor*)pc->data;

257:   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);
258:   dir->info.dtcol = dtcol;
259:   return(0);
260: }

264: PetscErrorCode  PCSetFromOptions_Factor(PetscOptionItems *PetscOptionsObject,PC pc)
265: {
266:   PC_Factor         *factor = (PC_Factor*)pc->data;
267:   PetscErrorCode    ierr;
268:   PetscBool         flg,set;
269:   char              tname[256], solvertype[64];
270:   PetscFunctionList ordlist;
271:   PetscEnum         etmp;
272:   PetscBool         inplace;

275:   PCFactorGetUseInPlace(pc,&inplace);
276:   PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",inplace,&flg,&set);
277:   if (set) {
278:     PCFactorSetUseInPlace(pc,flg);
279:   }
280:   PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,NULL);

282:   PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);
283:   if (flg) {
284:     PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);
285:   }
286:   PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);

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

291:   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);
292:   if (set) {
293:     PCFactorSetPivotInBlocks(pc,flg);
294:   }

296:   PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",PETSC_FALSE,&flg,&set);
297:   if (set) {
298:     PCFactorSetReuseFill(pc,flg);
299:   }
300:   PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",PETSC_FALSE,&flg,&set);
301:   if (set) {
302:     PCFactorSetReuseOrdering(pc,flg);
303:   }

305:   MatGetOrderingList(&ordlist);
306:   PetscOptionsFList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);
307:   if (flg) {
308:     PCFactorSetMatOrderingType(pc,tname);
309:   }

311:   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
312:   PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);
313:   if (flg) {
314:     PCFactorSetMatSolverPackage(pc,solvertype);
315:   }
316:   return(0);
317: }

321: PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
322: {
323:   PC_Factor      *factor = (PC_Factor*)pc->data;
325:   PetscBool      isstring,iascii;

328:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
329:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
330:   if (iascii) {
331:     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC) {
332:       if (factor->info.dt > 0) {
333:         PetscViewerASCIIPrintf(viewer,"  drop tolerance %g\n",(double)factor->info.dt);
334:         PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);
335:         PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %g\n",(double)factor->info.dtcol);
336:       } else if (factor->info.levels == 1) {
337:         PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);
338:       } else {
339:         PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);
340:       }
341:     }

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

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

350:     if (factor->fact) {
351:       MatInfo info;
352:       MatGetInfo(factor->fact,MAT_LOCAL,&info);
353:       PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %g, needed %g\n",(double)info.fill_ratio_given,(double)info.fill_ratio_needed);
354:       PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");
355:       PetscViewerASCIIPushTab(viewer);
356:       PetscViewerASCIIPushTab(viewer);
357:       PetscViewerASCIIPushTab(viewer);
358:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
359:       MatView(factor->fact,viewer);
360:       PetscViewerPopFormat(viewer);
361:       PetscViewerASCIIPopTab(viewer);
362:       PetscViewerASCIIPopTab(viewer);
363:       PetscViewerASCIIPopTab(viewer);
364:     }

366:   } else if (isstring) {
367:     MatFactorType t;
368:     MatGetFactorType(factor->fact,&t);
369:     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC) {
370:       PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);
371:     }
372:   }
373:   return(0);
374: }