Actual source code: lu.c

  1: /*$Id: lu.c,v 1.144 2001/04/10 19:36:09 bsmith Exp bsmith $*/
  2: /*
  3:    Defines a direct factorization preconditioner for any Mat implementation
  4:    Note: this need not be consided a preconditioner since it supplies
  5:          a direct solver.
  6: */
 7:  #include src/sles/pc/pcimpl.h

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


 21: EXTERN_C_BEGIN
 22: int PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag)
 23: {
 24:   PC_LU *lu;

 27:   lu                = (PC_LU*)pc->data;
 28:   lu->reuseordering = flag;
 29:   return(0);
 30: }
 31: EXTERN_C_END

 33: EXTERN_C_BEGIN
 34: int PCLUSetReuseFill_LU(PC pc,PetscTruth flag)
 35: {
 36:   PC_LU *lu;

 39:   lu = (PC_LU*)pc->data;
 40:   lu->reusefill = flag;
 41:   return(0);
 42: }
 43: EXTERN_C_END

 45: static int PCSetFromOptions_LU(PC pc)
 46: {
 47:   PC_LU      *lu = (PC_LU*)pc->data;
 48:   int        ierr;
 49:   PetscTruth flg;
 50:   char       tname[256];
 51:   PetscFList ordlist;
 52:   double     tol;

 55:   MatOrderingRegisterAll(PETSC_NULL);
 56:   PetscOptionsHead("LU options");
 57:     PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);
 58:     if (flg) {
 59:       PCLUSetUseInPlace(pc);
 60:     }
 61:     PetscOptionsDouble("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);

 63:     PetscOptionsHasName(pc->prefix,"-pc_lu_damping",&flg);
 64:     if (flg) {
 65:         PCLUSetDamping(pc,0.0);
 66:     }
 67:     PetscOptionsDouble("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);

 69:     PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);
 70:     if (flg) {
 71:       PCLUSetReuseFill(pc,PETSC_TRUE);
 72:     }
 73:     PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);
 74:     if (flg) {
 75:       PCLUSetReuseOrdering(pc,PETSC_TRUE);
 76:     }

 78:     MatGetOrderingList(&ordlist);
 79:     PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
 80:     if (flg) {
 81:       PCLUSetMatOrdering(pc,tname);
 82:     }
 83:     PetscOptionsDouble("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);

 85:     PetscOptionsDouble("-pc_lu_column_pivoting","Column pivoting tolerance (not used)","PCLUSetColumnPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);
 86:   PetscOptionsTail();
 87:   return(0);
 88: }

 90: static int PCView_LU(PC pc,PetscViewer viewer)
 91: {
 92:   PC_LU      *lu = (PC_LU*)pc->data;
 93:   int        ierr;
 94:   PetscTruth isascii,isstring;

 97:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
 98:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
 99:   if (isascii) {
100:     MatInfo info;

102:     if (lu->inplace) {PetscViewerASCIIPrintf(viewer,"  LU: in-place factorizationn");}
103:     else             {PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorizationn");}
104:     PetscViewerASCIIPrintf(viewer,"    matrix ordering: %sn",lu->ordering);
105:     if (lu->fact) {
106:       MatGetInfo(lu->fact,MAT_LOCAL,&info);
107:       PetscViewerASCIIPrintf(viewer,"    LU nonzeros %gn",info.nz_used);
108:     }
109:     if (lu->reusefill)    {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorizationn");}
110:     if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorizationn");}
111:   } else if (isstring) {
112:     PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
113:   } else {
114:     SETERRQ1(1,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
115:   }
116:   return(0);
117: }

119: static int PCGetFactoredMatrix_LU(PC pc,Mat *mat)
120: {
121:   PC_LU *dir = (PC_LU*)pc->data;

124:   if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
125:   *mat = dir->fact;
126:   return(0);
127: }

129: static int PCSetUp_LU(PC pc)
130: {
131:   int        ierr;
132:   PetscTruth flg;
133:   PC_LU      *dir = (PC_LU*)pc->data;

136:   if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;

138:   if (dir->inplace) {
139:     if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
140:     if (dir->col) {ISDestroy(dir->col);}
141:     MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
142:     if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
143:     MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);
144:     dir->fact = pc->pmat;
145:   } else {
146:     MatInfo info;
147:     if (!pc->setupcalled) {
148:       MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
149:       PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
150:       if (flg) {
151:         PetscReal tol = 1.e-10;
152:         PetscOptionsGetDouble(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
153:         MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
154:       }
155:       if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
156:       MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
157:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
158:       dir->actualfill = info.fill_ratio_needed;
159:       PetscLogObjectParent(pc,dir->fact);
160:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
161:       if (!dir->reuseordering) {
162:         if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
163:         if (dir->col) {ISDestroy(dir->col);}
164:         MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
165:         PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
166:         if (flg) {
167:           PetscReal tol = 1.e-10;
168:           PetscOptionsGetDouble(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
169:           MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
170:         }
171:         if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
172:       }
173:       MatDestroy(dir->fact);
174:       MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
175:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
176:       dir->actualfill = info.fill_ratio_needed;
177:       PetscLogObjectParent(pc,dir->fact);
178:     }
179:     MatLUFactorNumeric(pc->pmat,&dir->fact);
180:   }
181:   return(0);
182: }

184: static int PCDestroy_LU(PC pc)
185: {
186:   PC_LU *dir = (PC_LU*)pc->data;
187:   int   ierr;

190:   if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
191:   if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
192:   if (dir->col) {ISDestroy(dir->col);}
193:   PetscStrfree(dir->ordering);
194:   PetscFree(dir);
195:   return(0);
196: }

198: static int PCApply_LU(PC pc,Vec x,Vec y)
199: {
200:   PC_LU *dir = (PC_LU*)pc->data;
201:   int   ierr;

204:   if (dir->inplace) {MatSolve(pc->pmat,x,y);}
205:   else              {MatSolve(dir->fact,x,y);}
206:   return(0);
207: }

209: static int PCApplyTranspose_LU(PC pc,Vec x,Vec y)
210: {
211:   PC_LU *dir = (PC_LU*)pc->data;
212:   int   ierr;

215:   if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
216:   else              {MatSolveTranspose(dir->fact,x,y);}
217:   return(0);
218: }

220: /* -----------------------------------------------------------------------------------*/

222: EXTERN_C_BEGIN
223: int PCLUSetFill_LU(PC pc,PetscReal fill)
224: {
225:   PC_LU *dir;

228:   dir = (PC_LU*)pc->data;
229:   dir->info.fill = fill;
230:   return(0);
231: }
232: EXTERN_C_END

234: EXTERN_C_BEGIN
235: int PCLUSetDamping_LU(PC pc,PetscReal damping)
236: {
237:   PC_LU *dir;

240:   dir = (PC_LU*)pc->data;
241:   dir->info.damping = damping;
242:   dir->info.damp    = 1.0;
243:   return(0);
244: }
245: EXTERN_C_END

247: EXTERN_C_BEGIN
248: int PCLUSetUseInPlace_LU(PC pc)
249: {
250:   PC_LU *dir;

253:   dir = (PC_LU*)pc->data;
254:   dir->inplace = 1;
255:   return(0);
256: }
257: EXTERN_C_END

259: EXTERN_C_BEGIN
260: int PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
261: {
262:   PC_LU *dir = (PC_LU*)pc->data;
263:   int   ierr;

266:   PetscStrfree(dir->ordering);
267:   PetscStrallocpy(ordering,&dir->ordering);
268:   return(0);
269: }
270: EXTERN_C_END

272: EXTERN_C_BEGIN
273: int PCLUSetColumnPivoting_LU(PC pc,PetscReal dtcol)
274: {
275:   PC_LU *dir = (PC_LU*)pc->data;

278:   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(1,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
279:   dir->info.dtcol = dtcol;
280:   return(0);
281: }
282: EXTERN_C_END

284: /* -----------------------------------------------------------------------------------*/

286: /*@
287:    PCLUSetReuseOrdering - When similar matrices are factored, this
288:    causes the ordering computed in the first factor to be used for all
289:    following factors; applies to both fill and drop tolerance LUs.

291:    Collective on PC

293:    Input Parameters:
294: +  pc - the preconditioner context
295: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

297:    Options Database Key:
298: .  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()

300:    Level: intermediate

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

304: .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
305: @*/
306: int PCLUSetReuseOrdering(PC pc,PetscTruth flag)
307: {
308:   int ierr,(*f)(PC,PetscTruth);

312:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)())&f);
313:   if (f) {
314:     (*f)(pc,flag);
315:   }
316:   return(0);
317: }

319: /*@
320:    PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
321:    this causes later ones to use the fill computed in the initial factorization.

323:    Collective on PC

325:    Input Parameters:
326: +  pc - the preconditioner context
327: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

329:    Options Database Key:
330: .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()

332:    Level: intermediate

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

336: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
337: @*/
338: int PCLUSetReuseFill(PC pc,PetscTruth flag)
339: {
340:   int ierr,(*f)(PC,PetscTruth);

344:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)())&f);
345:   if (f) {
346:     (*f)(pc,flag);
347:   }
348:   return(0);
349: }

351: /*@
352:    PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
353:    fill = number nonzeros in factor/number nonzeros in original matrix.

355:    Collective on PC
356:    
357:    Input Parameters:
358: +  pc - the preconditioner context
359: -  fill - amount of expected fill

361:    Options Database Key:
362: .  -pc_lu_fill <fill> - Sets fill amount

364:    Level: intermediate

366:    Note:
367:    For sparse matrix factorizations it is difficult to predict how much 
368:    fill to expect. By running with the option -log_info PETSc will print the 
369:    actual amount of fill used; allowing you to set the value accurately for
370:    future runs. Default PETSc uses a value of 5.0

372: .keywords: PC, set, factorization, direct, fill

374: .seealso: PCILUSetFill()
375: @*/
376: int PCLUSetFill(PC pc,PetscReal fill)
377: {
378:   int ierr,(*f)(PC,PetscReal);

382:   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
383:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)())&f);
384:   if (f) {
385:     (*f)(pc,fill);
386:   }
387:   return(0);
388: }

390: /*@
391:    PCLUSetDamping - adds this quantity to the diagonal of the matrix during the 
392:      LU numerical factorization

394:    Collective on PC
395:    
396:    Input Parameters:
397: +  pc - the preconditioner context
398: -  damping - amount of damping

400:    Options Database Key:
401: .  -pc_lu_damping <damping> - Sets damping amount

403:    Level: intermediate

405: .keywords: PC, set, factorization, direct, fill

407: .seealso: PCILUSetFill(), PCILUSetDamp()
408: @*/
409: int PCLUSetDamping(PC pc,PetscReal damping)
410: {
411:   int ierr,(*f)(PC,PetscReal);

415:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)())&f);
416:   if (f) {
417:     (*f)(pc,damping);
418:   }
419:   return(0);
420: }

422: /*@
423:    PCLUSetUseInPlace - Tells the system to do an in-place factorization.
424:    For dense matrices, this enables the solution of much larger problems. 
425:    For sparse matrices the factorization cannot be done truly in-place 
426:    so this does not save memory during the factorization, but after the matrix
427:    is factored, the original unfactored matrix is freed, thus recovering that
428:    space.

430:    Collective on PC

432:    Input Parameters:
433: .  pc - the preconditioner context

435:    Options Database Key:
436: .  -pc_lu_in_place - Activates in-place factorization

438:    Notes:
439:    PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 
440:    a different matrix is provided for the multiply and the preconditioner in 
441:    a call to SLESSetOperators().
442:    This is because the Krylov space methods require an application of the 
443:    matrix multiplication, which is not possible here because the matrix has 
444:    been factored in-place, replacing the original matrix.

446:    Level: intermediate

448: .keywords: PC, set, factorization, direct, inplace, in-place, LU

450: .seealso: PCILUSetUseInPlace()
451: @*/
452: int PCLUSetUseInPlace(PC pc)
453: {
454:   int ierr,(*f)(PC);

458:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)())&f);
459:   if (f) {
460:     (*f)(pc);
461:   }
462:   return(0);
463: }

465: /*@
466:     PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 
467:     be used in the LU factorization.

469:     Collective on PC

471:     Input Parameters:
472: +   pc - the preconditioner context
473: -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM

475:     Options Database Key:
476: .   -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine

478:     Level: intermediate

480:     Notes: nested dissection is used by default

482: .seealso: PCILUSetMatOrdering()
483: @*/
484: int PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
485: {
486:   int ierr,(*f)(PC,MatOrderingType);

489:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)())&f);
490:   if (f) {
491:     (*f)(pc,ordering);
492:   }
493:   return(0);
494: }

496: /*@
497:     PCLUSetColumnPivoting - Determines when column pivoting is done during LU. 
498:       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
499:       it is never done. For the Matlab factorization this is used.

501:     Collective on PC

503:     Input Parameters:
504: +   pc - the preconditioner context
505: -   dtcol - 0.0 implies no pivoting, 1.0 complete column pivoting (slower, requires more memory but more stable)

507:     Options Database Key:
508: .   -pc_lu_column_pivoting - dttol

510:     Level: intermediate

512: .seealso: PCILUSetMatOrdering()
513: @*/
514: int PCLUSetColumnPivoting(PC pc,PetscReal dtcol)
515: {
516:   int ierr,(*f)(PC,PetscReal);

519:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetColumnPivoting_C",(void (**)())&f);
520:   if (f) {
521:     (*f)(pc,dtcol);
522:   }
523:   return(0);
524: }

526: /* ------------------------------------------------------------------------ */

528: EXTERN_C_BEGIN
529: int PCCreate_LU(PC pc)
530: {
531:   int   ierr;
532:   PC_LU *dir;

535:   PetscNew(PC_LU,&dir);
536:   PetscLogObjectMemory(pc,sizeof(PC_LU));

538:   dir->fact             = 0;
539:   dir->inplace          = 0;
540:   dir->info.fill        = 5.0;
541:   dir->info.dtcol       = 0.0; /* default to no pivoting; this is only thing PETSc LU supports */
542:   dir->info.damping     = 0.0;
543:   dir->info.damp        = 0.0;
544:   dir->col              = 0;
545:   dir->row              = 0;
546:   PetscStrallocpy(MATORDERING_ND,&dir->ordering);
547:   dir->reusefill        = PETSC_FALSE;
548:   dir->reuseordering    = PETSC_FALSE;
549:   pc->data              = (void*)dir;

551:   pc->ops->destroy           = PCDestroy_LU;
552:   pc->ops->apply             = PCApply_LU;
553:   pc->ops->applytranspose    = PCApplyTranspose_LU;
554:   pc->ops->setup             = PCSetUp_LU;
555:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
556:   pc->ops->view              = PCView_LU;
557:   pc->ops->applyrichardson   = 0;
558:   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;

560:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
561:                     PCLUSetFill_LU);
562:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
563:                     PCLUSetDamping_LU);
564:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
565:                     PCLUSetUseInPlace_LU);
566:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
567:                     PCLUSetMatOrdering_LU);
568:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
569:                     PCLUSetReuseOrdering_LU);
570:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
571:                     PCLUSetReuseFill_LU);
572:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetColumnPivoting_C","PCLUSetColumnPivoting_LU",
573:                     PCLUSetColumnPivoting_LU);
574:   return(0);
575: }
576: EXTERN_C_END