Actual source code: cholesky.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    Defines a direct factorization preconditioner for any Mat implementation
  5:    Note: this need not be consided a preconditioner since it supplies
  6:          a direct solver.
  7: */
 8:  #include private/pcimpl.h

 10: typedef struct {
 11:   Mat             fact;             /* factored matrix */
 12:   PetscReal       actualfill;       /* actual fill in factor */
 13:   PetscTruth      inplace;          /* flag indicating in-place factorization */
 14:   IS              row,col;          /* index sets used for reordering */
 15:   MatOrderingType ordering;         /* matrix ordering */
 16:   PetscTruth      reuseordering;    /* reuses previous reordering computed */
 17:   PetscTruth      reusefill;        /* reuse fill from previous Cholesky */
 18:   MatFactorInfo   info;
 19: } PC_Cholesky;

 24: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_Cholesky(PC pc,PetscReal z)
 25: {
 26:   PC_Cholesky *ch;

 29:   ch                 = (PC_Cholesky*)pc->data;
 30:   ch->info.zeropivot = z;
 31:   return(0);
 32: }

 38: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_Cholesky(PC pc,PetscReal shift)
 39: {
 40:   PC_Cholesky *dir;

 43:   dir = (PC_Cholesky*)pc->data;
 44:   if (shift == (PetscReal) PETSC_DECIDE) {
 45:     dir->info.shiftnz = 1.e-12;
 46:   } else {
 47:     dir->info.shiftnz = shift;
 48:   }
 49:   return(0);
 50: }

 56: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_Cholesky(PC pc,PetscTruth shift)
 57: {
 58:   PC_Cholesky *dir;
 59: 
 61:   dir = (PC_Cholesky*)pc->data;
 62:   if (shift) {
 63:     dir->info.shift_fraction = 0.0;
 64:     dir->info.shiftpd = 1.0;
 65:   } else {
 66:     dir->info.shiftpd = 0.0;
 67:   }
 68:   return(0);
 69: }

 75: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
 76: {
 77:   PC_Cholesky *lu;
 78: 
 80:   lu               = (PC_Cholesky*)pc->data;
 81:   lu->reuseordering = flag;
 82:   return(0);
 83: }

 89: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag)
 90: {
 91:   PC_Cholesky *lu;
 92: 
 94:   lu = (PC_Cholesky*)pc->data;
 95:   lu->reusefill = flag;
 96:   return(0);
 97: }

102: static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
103: {
104:   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
106:   PetscTruth     flg;
107:   char           tname[256];
108:   PetscFList     ordlist;
109: 
111:   MatOrderingRegisterAll(PETSC_NULL);
112:   PetscOptionsHead("Cholesky options");
113:   PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);
114:   if (flg) {
115:     PCFactorSetUseInPlace(pc);
116:   }
117:   PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",lu->info.fill,&lu->info.fill,0);
118: 
119:   PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);
120:   if (flg) {
121:     PCFactorSetReuseFill(pc,PETSC_TRUE);
122:   }
123:   PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);
124:   if (flg) {
125:     PCFactorSetReuseOrdering(pc,PETSC_TRUE);
126:   }
127: 
128:   MatGetOrderingList(&ordlist);
129:   PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
130:   if (flg) {
131:     PCFactorSetMatOrdering(pc,tname);
132:   }
133:   PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);
134:   if (flg) {
135:     PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);
136:   }
137:   PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",lu->info.shiftnz,&lu->info.shiftnz,0);
138:   PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);
139:   if (flg) {
140:     PCFactorSetShiftPd(pc,PETSC_TRUE);
141:   }
142:   PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);

144:   PetscOptionsTail();
145:   return(0);
146: }

150: static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
151: {
152:   PC_Cholesky    *lu = (PC_Cholesky*)pc->data;
154:   PetscTruth     iascii,isstring;
155: 
157:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
158:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
159:   if (iascii) {
160:     MatInfo info;
161: 
162:     if (lu->inplace) {PetscViewerASCIIPrintf(viewer,"  Cholesky: in-place factorization\n");}
163:     else             {PetscViewerASCIIPrintf(viewer,"  Cholesky: out-of-place factorization\n");}
164:     PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);
165:     if (lu->fact) {
166:       MatGetInfo(lu->fact,MAT_LOCAL,&info);
167:       PetscViewerASCIIPrintf(viewer,"    Cholesky nonzeros %G\n",info.nz_used);
168:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
169:       MatView(lu->fact,viewer);
170:       PetscViewerPopFormat(viewer);
171:     }
172:     if (lu->reusefill)    {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");}
173:     if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");}
174:   } else if (isstring) {
175:     PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
176:   } else {
177:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name);
178:   }
179:   return(0);
180: }

184: static PetscErrorCode PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
185: {
186:   PC_Cholesky *dir = (PC_Cholesky*)pc->data;
187: 
189:   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
190:   *mat = dir->fact;
191:   return(0);
192: }

196: static PetscErrorCode PCSetUp_Cholesky(PC pc)
197: {
199:   PetscTruth     flg;
200:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

203:   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
204: 
205:   if (dir->inplace) {
206:     if (dir->row && dir->col && (dir->row != dir->col)) {
207:       ISDestroy(dir->row);
208:       dir->row = 0;
209:     }
210:     if (dir->col) {
211:       ISDestroy(dir->col);
212:       dir->col = 0;
213:     }
214:     MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
215:     if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
216:       ISDestroy(dir->col);
217:       dir->col=0;
218:     }
219:     if (dir->row) {PetscLogObjectParent(pc,dir->row);}
220:     MatCholeskyFactor(pc->pmat,dir->row,&dir->info);
221:     dir->fact = pc->pmat;
222:   } else {
223:     MatInfo info;
224:     if (!pc->setupcalled) {
225:       MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
226:       if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
227:         ISDestroy(dir->col);
228:         dir->col=0;
229:       }
230:       PetscOptionsHasName(pc->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);
231:       if (flg) {
232:         PetscReal tol = 1.e-10;
233:         PetscOptionsGetReal(pc->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
234:         MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
235:       }
236:       if (dir->row) {PetscLogObjectParent(pc,dir->row);}
237:       MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
238:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
239:       dir->actualfill = info.fill_ratio_needed;
240:       PetscLogObjectParent(pc,dir->fact);
241:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
242:       if (!dir->reuseordering) {
243:         if (dir->row && dir->col && (dir->row != dir->col)) {
244:           ISDestroy(dir->row);
245:           dir->row = 0;
246:         }
247:         if (dir->col) {
248:           ISDestroy(dir->col);
249:           dir->col =0;
250:         }
251:         MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
252:         if (dir->col && (dir->row != dir->col)) {  /* only use row ordering for SBAIJ */
253:           ISDestroy(dir->col);
254:           dir->col=0;
255:         }
256:         PetscOptionsHasName(pc->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);
257:         if (flg) {
258:           PetscReal tol = 1.e-10;
259:           PetscOptionsGetReal(pc->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
260:           MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
261:         }
262:         if (dir->row) {PetscLogObjectParent(pc,dir->row);}
263:       }
264:       MatDestroy(dir->fact);
265:       MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
266:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
267:       dir->actualfill = info.fill_ratio_needed;
268:       PetscLogObjectParent(pc,dir->fact);
269:     }
270:     MatCholeskyFactorNumeric(pc->pmat,&dir->info,&dir->fact);
271:   }
272:   return(0);
273: }

277: static PetscErrorCode PCDestroy_Cholesky(PC pc)
278: {
279:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

283:   if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
284:   if (dir->row) {ISDestroy(dir->row);}
285:   if (dir->col) {ISDestroy(dir->col);}
286:   PetscStrfree(dir->ordering);
287:   PetscFree(dir);
288:   return(0);
289: }

293: static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
294: {
295:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;
297: 
299:   if (dir->inplace) {MatSolve(pc->pmat,x,y);}
300:   else              {MatSolve(dir->fact,x,y);}
301:   return(0);
302: }

306: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
307: {
308:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

312:   if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
313:   else              {MatSolveTranspose(dir->fact,x,y);}
314:   return(0);
315: }

317: /* -----------------------------------------------------------------------------------*/

322: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_Cholesky(PC pc,PetscReal fill)
323: {
324:   PC_Cholesky *dir;
325: 
327:   dir = (PC_Cholesky*)pc->data;
328:   dir->info.fill = fill;
329:   return(0);
330: }

336: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_Cholesky(PC pc)
337: {
338:   PC_Cholesky *dir;

341:   dir = (PC_Cholesky*)pc->data;
342:   dir->inplace = PETSC_TRUE;
343:   return(0);
344: }

350: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
351: {
352:   PC_Cholesky    *dir = (PC_Cholesky*)pc->data;

356:   PetscStrfree(dir->ordering);
357:   PetscStrallocpy(ordering,&dir->ordering);
358:   return(0);
359: }

362: /* -----------------------------------------------------------------------------------*/

366: /*@
367:    PCFactorSetReuseOrdering - When similar matrices are factored, this
368:    causes the ordering computed in the first factor to be used for all
369:    following factors.

371:    Collective on PC

373:    Input Parameters:
374: +  pc - the preconditioner context
375: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

377:    Options Database Key:
378: .  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()

380:    Level: intermediate

382: .keywords: PC, levels, reordering, factorization, incomplete, LU

384: .seealso: PCFactorSetReuseFill()
385: @*/
386: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
387: {
388:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

392:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);
393:   if (f) {
394:     (*f)(pc,flag);
395:   }
396:   return(0);
397: }

401: /*@
402:    PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
403:    this causes later ones to use the fill ratio computed in the initial factorization.

405:    Collective on PC

407:    Input Parameters:
408: +  pc - the preconditioner context
409: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

411:    Options Database Key:
412: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()

414:    Level: intermediate

416: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky

418: .seealso: PCFactorSetReuseOrdering()
419: @*/
420: PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill(PC pc,PetscTruth flag)
421: {
422:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

426:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);
427:   if (f) {
428:     (*f)(pc,flag);
429:   }
430:   return(0);
431: }

433: /*MC
434:    PCCholesky - Uses a direct solver, based on Cholesky factorization, as a preconditioner

436:    Options Database Keys:
437: +  -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
438: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
439: .  -pc_factor_fill <fill> - Sets fill amount
440: .  -pc_factor_in_place - Activates in-place factorization
441: .  -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
442: .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
443: -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
444:    is optional with PETSC_TRUE being the default

446:    Notes: Not all options work for all matrix formats

448:    Level: beginner

450:    Concepts: Cholesky factorization, direct solver

452:    Notes: Usually this will compute an "exact" solution in one iteration and does 
453:           not need a Krylov method (i.e. you can use -ksp_type preonly, or 
454:           KSPSetType(ksp,KSPPREONLY) for the Krylov method

456: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
457:            PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCGetFactoredMatrix(),
458:            PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
459:            PCFactorSetUseInPlace(), PCFactorSetMatOrdering(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()

461: M*/

466: PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Cholesky(PC pc)
467: {
469:   PC_Cholesky    *dir;

472:   PetscNew(PC_Cholesky,&dir);
473:   PetscLogObjectMemory(pc,sizeof(PC_Cholesky));

475:   dir->fact                   = 0;
476:   dir->inplace                = PETSC_FALSE;
477:   MatFactorInfoInitialize(&dir->info);
478:   dir->info.fill              = 5.0;
479:   dir->info.shiftnz           = 0.0;
480:   dir->info.shiftpd           = 0.0; /* false */
481:   dir->info.shift_fraction    = 0.0;
482:   dir->info.pivotinblocks     = 1.0;
483:   dir->col                    = 0;
484:   dir->row                    = 0;
485:   PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
486:   dir->reusefill        = PETSC_FALSE;
487:   dir->reuseordering    = PETSC_FALSE;
488:   pc->data              = (void*)dir;

490:   pc->ops->destroy           = PCDestroy_Cholesky;
491:   pc->ops->apply             = PCApply_Cholesky;
492:   pc->ops->applytranspose    = PCApplyTranspose_Cholesky;
493:   pc->ops->setup             = PCSetUp_Cholesky;
494:   pc->ops->setfromoptions    = PCSetFromOptions_Cholesky;
495:   pc->ops->view              = PCView_Cholesky;
496:   pc->ops->applyrichardson   = 0;
497:   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;

499:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Cholesky",
500:                     PCFactorSetZeroPivot_Cholesky);
501:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Cholesky",
502:                     PCFactorSetShiftNonzero_Cholesky);
503:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Cholesky",
504:                     PCFactorSetShiftPd_Cholesky);

506:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Cholesky",
507:                     PCFactorSetFill_Cholesky);
508:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
509:                     PCFactorSetUseInPlace_Cholesky);
510:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrdering_C","PCFactorSetMatOrdering_Cholesky",
511:                     PCFactorSetMatOrdering_Cholesky);
512:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
513:                     PCFactorSetReuseOrdering_Cholesky);
514:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
515:                     PCFactorSetReuseFill_Cholesky);
516:   return(0);
517: }