Actual source code: lu.c

  1: /*$Id: lu.c,v 1.149 2001/08/07 21:30:26 bsmith Exp $*/
  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:   PetscReal  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:     PetscOptionsReal("-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:     PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);
 68:     PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);

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

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

 86:     PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for SuperLU)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);
 87:   PetscOptionsTail();
 88:   return(0);
 89: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

222: /* -----------------------------------------------------------------------------------*/

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

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

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

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

249: EXTERN_C_BEGIN
250: int PCLUSetUseInPlace_LU(PC pc)
251: {
252:   PC_LU *dir;

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

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

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

274: EXTERN_C_BEGIN
275: int PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
276: {
277:   PC_LU *dir = (PC_LU*)pc->data;

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

286: /* -----------------------------------------------------------------------------------*/

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

293:    Collective on PC

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

299:    Options Database Key:
300: .  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()

302:    Level: intermediate

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

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

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

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

325:    Collective on PC

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

331:    Options Database Key:
332: .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()

334:    Level: intermediate

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

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

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

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

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

363:    Options Database Key:
364: .  -pc_lu_fill <fill> - Sets fill amount

366:    Level: intermediate

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

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

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

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

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

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

402:    Options Database Key:
403: .  -pc_lu_damping <damping> - Sets damping amount

405:    Level: intermediate

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

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

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

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

432:    Collective on PC

434:    Input Parameters:
435: .  pc - the preconditioner context

437:    Options Database Key:
438: .  -pc_lu_in_place - Activates in-place factorization

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

448:    Level: intermediate

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

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

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

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

471:     Collective on PC

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

477:     Options Database Key:
478: .   -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine

480:     Level: intermediate

482:     Notes: nested dissection is used by default

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

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

498: /*@
499:     PCLUSetPivoting - Determines when pivoting is done during LU. 
500:       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
501:       it is never done. For the Matlab and SuperLU factorization this is used.

503:     Collective on PC

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

509:     Options Database Key:
510: .   -pc_lu_pivoting - dttol

512:     Level: intermediate

514: .seealso: PCILUSetMatOrdering()
515: @*/
516: int PCLUSetPivoting(PC pc,PetscReal dtcol)
517: {
518:   int ierr,(*f)(PC,PetscReal);

521:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);
522:   if (f) {
523:     (*f)(pc,dtcol);
524:   }
525:   return(0);
526: }

528: /* ------------------------------------------------------------------------ */

530: EXTERN_C_BEGIN
531: int PCCreate_LU(PC pc)
532: {
533:   int   ierr,size;
534:   PC_LU *dir;

537:   PetscNew(PC_LU,&dir);
538:   PetscLogObjectMemory(pc,sizeof(PC_LU));

540:   dir->fact             = 0;
541:   dir->inplace          = 0;
542:   dir->info.fill        = 5.0;
543:   dir->info.dtcol       = 0.0; /* default to no pivoting; this is only thing PETSc LU supports */
544:   dir->info.damping     = 0.0;
545:   dir->info.damp        = 0.0;
546:   dir->info.zeropivot   = 1.e-12;
547:   dir->col              = 0;
548:   dir->row              = 0;
549:   MPI_Comm_size(pc->comm,&size);
550:   if (size == 1) {
551:     PetscStrallocpy(MATORDERING_ND,&dir->ordering);
552:   } else {
553:     PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
554:   }
555:   dir->reusefill        = PETSC_FALSE;
556:   dir->reuseordering    = PETSC_FALSE;
557:   pc->data              = (void*)dir;

559:   pc->ops->destroy           = PCDestroy_LU;
560:   pc->ops->apply             = PCApply_LU;
561:   pc->ops->applytranspose    = PCApplyTranspose_LU;
562:   pc->ops->setup             = PCSetUp_LU;
563:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
564:   pc->ops->view              = PCView_LU;
565:   pc->ops->applyrichardson   = 0;
566:   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;

568:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
569:                     PCLUSetFill_LU);
570:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
571:                     PCLUSetDamping_LU);
572:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
573:                     PCLUSetUseInPlace_LU);
574:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
575:                     PCLUSetMatOrdering_LU);
576:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
577:                     PCLUSetReuseOrdering_LU);
578:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
579:                     PCLUSetReuseFill_LU);
580:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
581:                     PCLUSetPivoting_LU);
582:   return(0);
583: }
584: EXTERN_C_END