Actual source code: factimpl.c

  1: #define PETSCKSP_DLL

 3:  #include ../src/ksp/pc/impls/factor/factor.h

  5: /* ------------------------------------------------------------------------------------------*/

 10: PetscErrorCode  PCFactorSetZeroPivot_Factor(PC pc,PetscReal z)
 11: {
 12:   PC_Factor *ilu = (PC_Factor*)pc->data;

 15:   ilu->info.zeropivot = z;
 16:   return(0);
 17: }

 23: PetscErrorCode  PCFactorSetShiftType_Factor(PC pc,MatFactorShiftType shifttype)
 24: {
 25:   PC_Factor *dir = (PC_Factor*)pc->data;

 28:   if (shifttype == (MatFactorShiftType)PETSC_DECIDE){
 29:     dir->info.shifttype = (PetscReal) MAT_SHIFT_NONE;
 30:   } else {
 31:     dir->info.shifttype = (PetscReal) shifttype;
 32:   }
 33:   return(0);
 34: }

 38: PetscErrorCode  PCFactorSetShiftAmount_Factor(PC pc,PetscReal shiftamount)
 39: {
 40:   PC_Factor *dir = (PC_Factor*)pc->data;

 43:   if (shiftamount == (PetscReal) PETSC_DECIDE){
 44:     dir->info.shiftamount = 1.e-12;
 45:   } else {
 46:     dir->info.shiftamount = shiftamount;
 47:   }
 48:   return(0);
 49: }

 55: PetscErrorCode  PCFactorSetDropTolerance_Factor(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
 56: {
 57:   PC_Factor         *ilu = (PC_Factor*)pc->data;

 60:   if (pc->setupcalled && (!ilu->info.usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) {
 61:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change tolerance after use");
 62:   }
 63:   ilu->info.usedt   = PETSC_TRUE;
 64:   ilu->info.dt      = dt;
 65:   ilu->info.dtcol   = dtcol;
 66:   ilu->info.dtcount = dtcount;
 67:   ilu->info.fill    = PETSC_DEFAULT;
 68:   return(0);
 69: }

 75: PetscErrorCode  PCFactorSetFill_Factor(PC pc,PetscReal fill)
 76: {
 77:   PC_Factor *dir = (PC_Factor*)pc->data;

 80:   dir->info.fill = fill;
 81:   return(0);
 82: }

 88: PetscErrorCode  PCFactorSetMatOrderingType_Factor(PC pc,const MatOrderingType ordering)
 89: {
 90:   PC_Factor      *dir = (PC_Factor*)pc->data;
 92:   PetscTruth     flg;
 93: 
 95:   if (!pc->setupcalled) {
 96:      PetscStrfree(dir->ordering);
 97:      PetscStrallocpy(ordering,&dir->ordering);
 98:   } else {
 99:     PetscStrcmp(dir->ordering,ordering,&flg);
100:     if (!flg) {
101:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change ordering after use");
102:     }
103:   }
104:   return(0);
105: }

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

116:   if (!pc->setupcalled) {
117:     ilu->info.levels = levels;
118:   } else if (ilu->info.usedt || ilu->info.levels != levels) {
119:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change levels after use");
120:   }
121:   return(0);
122: }

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


139: /* ------------------------------------------------------------------------------------------*/

144: PetscErrorCode  PCFactorSetPivotInBlocks_Factor(PC pc,PetscTruth pivot)
145: {
146:   PC_Factor *dir = (PC_Factor*)pc->data;

149:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
150:   return(0);
151: }

157: PetscErrorCode  PCFactorGetMatrix_Factor(PC pc,Mat *mat)
158: {
159:   PC_Factor *ilu = (PC_Factor*)pc->data;

162:   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
163:   *mat = ilu->fact;
164:   return(0);
165: }

171: PetscErrorCode  PCFactorSetMatSolverPackage_Factor(PC pc,const MatSolverPackage stype)
172: {
174:   PC_Factor      *lu = (PC_Factor*)pc->data;

177:   if (lu->fact) {
178:     const MatSolverPackage ltype;
179:     PetscTruth             flg;
180:     MatFactorGetSolverPackage(lu->fact,&ltype);
181:     PetscStrcmp(stype,ltype,&flg);
182:     if (!flg) {
183:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot change solver matrix package after PC has been setup or used");
184:     }
185:   }
186:   PetscStrfree(lu->solvertype);
187:   PetscStrallocpy(stype,&lu->solvertype);
188:   return(0);
189: }

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

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

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

213:   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %G must be between 0 and 1",dtcol);
214:   dir->info.dtcol = dtcol;
215:   return(0);
216: }

222: PetscErrorCode  PCSetFromOptions_Factor(PC pc)
223: {
224:   PC_Factor       *factor = (PC_Factor*)pc->data;
225:   PetscErrorCode  ierr;
226:   PetscTruth      flg = PETSC_FALSE,set;
227:   char            tname[256], solvertype[64];
228:   PetscFList      ordlist;
229:   PetscEnum       etmp;

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

239:   PetscOptionsEnum("-pc_factor_shift_type","Shift added to diagonal","PCFactorSetShiftType",
240:                           MatFactorShiftTypes,(PetscEnum)(int)((PC_Factor*)factor)->info.shifttype,&etmp,&flg);
241:   if (flg) {
242:     ((PC_Factor*)factor)->info.shifttype = (PetscReal)(etmp);
243:   }
244:   PetscOptionsReal("-pc_factor_shift_amount","Shift added to diagonal","PCFactorSetShiftAmount",((PC_Factor*)factor)->info.shiftamount,&((PC_Factor*)factor)->info.shiftamount,0);
245: 
246:   PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)factor)->info.zeropivot,&((PC_Factor*)factor)->info.zeropivot,0);
247:   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);

249:   flg = ((PC_Factor*)factor)->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
250:   PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix dense blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);
251:   if (set) {
252:     PCFactorSetPivotInBlocks(pc,flg);
253:   }
254: 
255:   flg  = PETSC_FALSE;
256:   PetscOptionsTruth("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",flg,&flg,PETSC_NULL);
257:   if (flg) {
258:     PCFactorSetReuseFill(pc,PETSC_TRUE);
259:   }
260:   flg  = PETSC_FALSE;
261:   PetscOptionsTruth("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",flg,&flg,PETSC_NULL);
262:   if (flg) {
263:     PCFactorSetReuseOrdering(pc,PETSC_TRUE);
264:   }

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

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

282:  #include private/matimpl.h
285: PetscErrorCode PCView_Factor(PC pc,PetscViewer viewer)
286: {
287:   PC_Factor       *factor = (PC_Factor*)pc->data;
288:   PetscErrorCode  ierr;
289:   PetscTruth      isstring,iascii;

292:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
293:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
294:   if (iascii) {
295:     if (factor->factortype == MAT_FACTOR_ILU || factor->factortype == MAT_FACTOR_ICC){
296:       if (factor->info.dt > 0) {
297:         PetscViewerASCIIPrintf(viewer,"  drop tolerance %G\n",factor->info.dt);
298:         PetscViewerASCIIPrintf(viewer,"  max nonzeros per row %D\n",factor->info.dtcount);
299:         PetscViewerASCIIPrintf(viewer,"  column permutation tolerance %G\n",(PetscInt)factor->info.dtcol);
300:       } else if (factor->info.levels == 1) {
301:         PetscViewerASCIIPrintf(viewer,"  %D level of fill\n",(PetscInt)factor->info.levels);
302:       } else {
303:         PetscViewerASCIIPrintf(viewer,"  %D levels of fill\n",(PetscInt)factor->info.levels);
304:       }
305:     }
306: 
307:     PetscViewerASCIIPrintf(viewer,"  tolerance for zero pivot %G\n",factor->info.zeropivot);
308:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
309:       PetscViewerASCIIPrintf(viewer,"  using Manteuffel shift\n");
310:     }
311:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_NONZERO) {
312:       PetscViewerASCIIPrintf(viewer,"  using diagonal shift to prevent zero pivot\n");
313:     }
314:     if (factor->info.shifttype==(PetscReal)MAT_SHIFT_INBLOCKS) {
315:       PetscViewerASCIIPrintf(viewer,"  using diagonal shift on blocks to prevent zero pivot\n");
316:     }
317: 
318:     PetscViewerASCIIPrintf(viewer,"  matrix ordering: %s\n",factor->ordering);
319: 
320:     if (factor->fact) {
321:       PetscViewerASCIIPrintf(viewer,"  factor fill ratio given %G, needed %G\n",factor->fact->info.fill_ratio_given,factor->fact->info.fill_ratio_needed);
322:       PetscViewerASCIIPrintf(viewer,"    Factored matrix follows:\n");
323:       PetscViewerASCIIPushTab(viewer);
324:       PetscViewerASCIIPushTab(viewer);
325:       PetscViewerASCIIPushTab(viewer);
326:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
327:       MatView(factor->fact,viewer);
328:       PetscViewerPopFormat(viewer);
329:       PetscViewerASCIIPopTab(viewer);
330:       PetscViewerASCIIPopTab(viewer);
331:       PetscViewerASCIIPopTab(viewer);
332:     }
333: 
334:   } else if (isstring) {
335:     if (factor->fact->factor == MAT_FACTOR_ILU || factor->fact->factor == MAT_FACTOR_ICC){
336:       PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)factor->info.levels,factor->ordering);
337:     }
338:   } else {
339:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC_Factor",((PetscObject)viewer)->type_name);
340:   }
341:   return(0);
342: }