Actual source code: factimpl.c

  2: #include <../src/ksp/pc/impls/factor/factor.h>     /*I "petscpc.h"  I*/

  4: /* ------------------------------------------------------------------------------------------*/


 10: PetscErrorCode PCFactorSetUpMatSolverPackage_Factor(PC pc)
 11: {
 12:   PC_Factor      *icc = (PC_Factor*)pc->data;

 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: }

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

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

 39: PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
 40: {
 41:   PC_Factor *dir = (PC_Factor*)pc->data;

 44:   if (shifttype == (MatFactorShiftType)PETSC_DECIDE){
 45:     dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
 46:   } else {
 47:     dir->info.shifttype = (PetscReal) shifttype;
 48:     if (shifttype == MAT_SHIFT_NONZERO && dir->info.shiftamount == 0.0){
 49:       dir->info.shiftamount = 1.e-12; /* set default amount if user has not called PCFactorSetShiftAmount() yet */
 50:     }
 51:   }
 52:   return(0);
 53: }

 57: PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
 58: {
 59:   PC_Factor *dir = (PC_Factor*)pc->data;

 62:   if (shiftamount == (PetscReal) PETSC_DECIDE){
 63:     dir->info.shiftamount = 1.e-12;
 64:   } else {
 65:     dir->info.shiftamount = shiftamount;
 66:   }
 67:   return(0);
 68: }

 74: PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 75: {
 76:   PC_Factor         *ilu = (PC_Factor*)pc->data;

 79:   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 80:     SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
 81:   }
 82:   ilu->info.usedt   = PETSC_TRUE;
 83:   ilu->info.dt      = dt;
 84:   ilu->info.dtcol   = dtcol;
 85:   ilu->info.dtcount = dtcount;
 86:   ilu->info.fill    = PETSC_DEFAULT;
 87:   return(0);
 88: }

 94: PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
 95: {
 96:   PC_Factor *dir = (PC_Factor*)pc->data;

 99:   dir->info.fill = fill;
100:   return(0);
101: }

107: PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
108: {
109:   PC_Factor      *dir = (PC_Factor*)pc->data;
111:   PetscBool      flg;
112: 
114:   if (!pc->setupcalled) {
115:      PetscFree(dir->ordering);
116:      PetscStrallocpy(ordering,&dir->ordering);
117:   } else {
118:     PetscStrcmp(dir->ordering,ordering,&flg);
119:     if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
120:   }
121:   return(0);
122: }

128: PetscErrorCode  PCFactorSetLevels_Factor(PC pc,PetscInt levels)
129: {
130:   PC_Factor      *ilu = (PC_Factor*)pc->data;

133:   if (!pc->setupcalled) {
134:     ilu->info.levels = levels;
135:   } else if (ilu->info.usedt || ilu->info.levels != levels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
136:   return(0);
137: }

143: PetscErrorCode  PCFactorSetAllowDiagonalFill_Factor(PC pc)
144: {
145:   PC_Factor *dir = (PC_Factor*)pc->data;

148:   dir->info.diagonal_fill = 1;
149:   return(0);
150: }


154: /* ------------------------------------------------------------------------------------------*/

159: PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscBool  pivot)
160: {
161:   PC_Factor *dir = (PC_Factor*)pc->data;

164:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
165:   return(0);
166: }

171: PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
172: {
173:   PC_Factor *ilu = (PC_Factor*)pc->data;

176:   if (!ilu->fact) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
177:   *mat = ilu->fact;
178:   return(0);
179: }

184: PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
185: {
187:   PC_Factor      *lu = (PC_Factor*)pc->data;

190:   if (lu->fact) {
191:     const MatSolverPackage ltype;
192:     PetscBool              flg;
193:     MatFactorGetSolverPackage(lu->fact,&ltype);
194:     PetscStrcmp(stype,ltype,&flg);
195:     if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
196:   } else {
197:     PetscFree(lu->solvertype);
198:     PetscStrallocpy(stype,&lu->solvertype);
199:   }
200:   return(0);
201: }

207: PetscErrorCode  PCFactorGetMatSolverPackage_Factor(PC pc,const MatSolverPackage *stype)
208: {
209:   PC_Factor      *lu = (PC_Factor*)pc->data;

212:   *stype = lu->solvertype;
213:   return(0);
214: }

220: PetscErrorCode  PCFactorSetColumnPivot_Factor(PC pc,PetscReal dtcol)
221: {
222:   PC_Factor       *dir = (PC_Factor*)pc->data;

225:   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
226:   dir->info.dtcol = dtcol;
227:   return(0);
228: }

234: PetscErrorCode  PCSetFromOptions_Factor(PC pc)
235: {
236:   PC_Factor       *factor = (PC_Factor*)pc->data;
237:   PetscErrorCode  ierr;
238:   PetscBool       flg = PETSC_FALSE,set;
239:   char            tname[256], solvertype[64];
240:   PetscFList      ordlist;
241:   PetscEnum       etmp;

244:   if (!MatOrderingRegisterAllCalled) {MatOrderingRegisterAll(PETSC_NULL);}
245:   PetscOptionsBool("-pc_factor_in_place","Form factored matrix in the same memory as the matrix","PCFactorSetUseInPlace",flg,&flg,PETSC_NULL);
246:   if (flg) {
247:     PCFactorSetUseInPlace(pc);
248:   }
249:   PetscOptionsReal("-pc_factor_fill","Expected non-zeros in factored matrix","PCFactorSetFill",((PC_Factor*)factor)->info.fill,&((PC_Factor*)factor)->info.fill,0);

251:   PetscOptionsEnum("-pc_factor_shift_type","Type of shift to add to diagonal","PCFactorSetShiftType",
252:                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);
253:   if (flg) {
254:     PCFactorSetShiftType(pc,(MatFactorShiftType)etmp);
255:   }
256:   PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);
257: 
258:   PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);
259:   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);

261:   flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
262:   PetscOptionsBool("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);
263:   if (set) {
264:     PCFactorSetPivotInBlocks(pc,flg);
265:   }
266: 
267:   flg  = PETSC_FALSE;
268:   PetscOptionsBool("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);
269:   if (flg) {
270:     PCFactorSetReuseFill(pc,PETSC_TRUE);
271:   }
272:   flg  = PETSC_FALSE;
273:   PetscOptionsBool("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);
274:   if (flg) {
275:     PCFactorSetReuseOrdering(pc,PETSC_TRUE);
276:   }

278:   MatGetOrderingList(&ordlist);
279:   PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in factored matrix","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)factor)->ordering,tname,256,&flg);
280:   if (flg) {
281:     PCFactorSetMatOrderingType(pc,tname);
282:   }

284:   /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
285:   PetscOptionsString("-pc_factor_mat_solver_package","Specific direct solver to use","MatGetFactor",((PC_Factor*)factor)->solvertype,solvertype,64,&flg);
286:   if (flg) {
287:     PCFactorSetMatSolverPackage(pc,solvertype);
288:   }
289:   return(0);
290: }

296: PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
297: {
298:   PC_Factor       *factor = (PC_Factor*)pc->data;
299:   PetscErrorCode  ierr;
300:   PetscBool       isstring,iascii;

303:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
304:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
305:   if (iascii) {
306:     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC){
307:       if (factor->info.dt > 0) {
308:         PetscViewerASCIIPrintf(viewer,"  drop tolerance %G\n",factor->info.dt);
309:         PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);
310:         PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);
311:       } else if (factor->info.levels == 1) {
312:         PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);
313:       } else {
314:         PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);
315:       }
316:     }
317: 
318:     PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %G\n",factor->info.zeropivot);
319:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
320:       PetscViewerASCIIPrintf(viewer,"  using Manteuffel shift\n");
321:     }
322:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) {
323:       PetscViewerASCIIPrintf(viewer,"  using diagonal shift to prevent zero pivot\n");
324:     }
325:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) {
326:       PetscViewerASCIIPrintf(viewer,"  using diagonal shift on blocks to prevent zero pivot\n");
327:     }
328: 
329:     PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);
330: 
331:     if (factor->fact) {
332:       MatInfo info;
333:       MatGetInfo(factor->fact,MAT_LOCAL,&info);
334:       PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %G, needed %G\n",info.fill_ratio_given,info.fill_ratio_needed);
335:       PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");
336:       PetscViewerASCIIPushTab(viewer);
337:       PetscViewerASCIIPushTab(viewer);
338:       PetscViewerASCIIPushTab(viewer);
339:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
340:       MatView(factor->fact,viewer);
341:       PetscViewerPopFormat(viewer);
342:       PetscViewerASCIIPopTab(viewer);
343:       PetscViewerASCIIPopTab(viewer);
344:       PetscViewerASCIIPopTab(viewer);
345:     }
346: 
347:   } else if (isstring) {
348:     MatFactorType t;
349:     MatGetFactorType(factor->fact,&t);
350:     if (t == MAT_FACTOR_ILU || t == MAT_FACTOR_ICC){
351:       PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);
352:     }
353:   } else {
354:     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC_Factor",((PetscObject)viewer)->type_name);
355:   }
356:   return(0);
357: }