Actual source code: lu.c

  1: /*
  2:    Defines a direct factorization preconditioner for any Mat implementation
  3:    Note: this need not be consided a preconditioner since it supplies
  4:          a direct solver.
  5: */
 6:  #include src/ksp/pc/pcimpl.h

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


 23: PetscErrorCode PCLUSetZeroPivot_LU(PC pc,PetscReal z)
 24: {
 25:   PC_LU *lu;

 28:   lu                 = (PC_LU*)pc->data;
 29:   lu->info.zeropivot = z;
 30:   return(0);
 31: }

 37: PetscErrorCode PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag)
 38: {
 39:   PC_LU *lu;

 42:   lu                = (PC_LU*)pc->data;
 43:   lu->reuseordering = flag;
 44:   return(0);
 45: }

 51: PetscErrorCode PCLUSetReuseFill_LU(PC pc,PetscTruth flag)
 52: {
 53:   PC_LU *lu;

 56:   lu = (PC_LU*)pc->data;
 57:   lu->reusefill = flag;
 58:   return(0);
 59: }

 65: PetscErrorCode PCLUSetShift_LU(PC pc,PetscTruth shift)
 66: {
 67:   PC_LU *dir;
 68: 
 70:   dir = (PC_LU*)pc->data;
 71:   dir->info.shift = shift;
 72:   if (shift) dir->info.shift_fraction = 0.0;
 73:   return(0);
 74: }

 79: static PetscErrorCode PCSetFromOptions_LU(PC pc)
 80: {
 81:   PC_LU          *lu = (PC_LU*)pc->data;
 83:   PetscTruth     flg,set;
 84:   char           tname[256];
 85:   PetscFList     ordlist;
 86:   PetscReal      tol;

 89:   MatOrderingRegisterAll(PETSC_NULL);
 90:   PetscOptionsHead("LU options");
 91:     PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);
 92:     if (flg) {
 93:       PCLUSetUseInPlace(pc);
 94:     }
 95:     PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);

 97:     PetscOptionsName("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",&flg);
 98:     if (flg) {
 99:         PCLUSetDamping(pc,(PetscReal) PETSC_DECIDE);
100:     }
101:     PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);
102:     PetscOptionsName("-pc_lu_shift","Manteuffel shift applied to diagonal","PCLUSetShift",&flg);
103:     if (flg) {
104:       PCLUSetShift(pc,PETSC_TRUE);
105:     }
106:     PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);

108:     PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);
109:     if (flg) {
110:       PCLUSetReuseFill(pc,PETSC_TRUE);
111:     }
112:     PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);
113:     if (flg) {
114:       PCLUSetReuseOrdering(pc,PETSC_TRUE);
115:     }

117:     MatGetOrderingList(&ordlist);
118:     PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
119:     if (flg) {
120:       PCLUSetMatOrdering(pc,tname);
121:     }
122:     PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);

124:     PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);

126:     flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
127:     PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);
128:     if (set) {
129:       PCLUSetPivotInBlocks(pc,flg);
130:     }
131:   PetscOptionsTail();
132:   return(0);
133: }

137: static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer)
138: {
139:   PC_LU          *lu = (PC_LU*)pc->data;
141:   PetscTruth     iascii,isstring;

144:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
145:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
146:   if (iascii) {
147:     MatInfo info;

149:     if (lu->inplace) {PetscViewerASCIIPrintf(viewer,"  LU: in-place factorization\n");}
150:     else             {PetscViewerASCIIPrintf(viewer,"  LU: out-of-place factorization\n");}
151:     PetscViewerASCIIPrintf(viewer,"    matrix ordering: %s\n",lu->ordering);
152:     PetscViewerASCIIPrintf(viewer,"  LU: tolerance for zero pivot %g\n",lu->info.zeropivot);
153:     if (lu->info.shift) {PetscViewerASCIIPrintf(viewer,"  LU: using Manteuffel shift\n");}
154:     if (lu->fact) {
155:       MatGetInfo(lu->fact,MAT_LOCAL,&info);
156:       PetscViewerASCIIPrintf(viewer,"    LU nonzeros %g\n",info.nz_used);
157:       PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
158:       MatView(lu->fact,viewer);
159:       PetscViewerPopFormat(viewer);
160:     }
161:     if (lu->reusefill)    {PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");}
162:     if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");}
163:   } else if (isstring) {
164:     PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
165:   } else {
166:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
167:   }
168:   return(0);
169: }

173: static PetscErrorCode PCGetFactoredMatrix_LU(PC pc,Mat *mat)
174: {
175:   PC_LU *dir = (PC_LU*)pc->data;

178:   if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
179:   *mat = dir->fact;
180:   return(0);
181: }

185: static PetscErrorCode PCSetUp_LU(PC pc)
186: {
188:   PetscTruth     flg;
189:   PC_LU          *dir = (PC_LU*)pc->data;

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

194:   if (dir->inplace) {
195:     if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
196:     if (dir->col) {ISDestroy(dir->col);}
197:     MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
198:     if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
199:     MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);
200:     dir->fact = pc->pmat;
201:   } else {
202:     MatInfo info;
203:     if (!pc->setupcalled) {
204:       MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
205:       PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
206:       if (flg) {
207:         PetscReal tol = 1.e-10;
208:         PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
209:         MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
210:       }
211:       if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
212:       MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
213:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
214:       dir->actualfill = info.fill_ratio_needed;
215:       PetscLogObjectParent(pc,dir->fact);
216:     } else if (pc->flag != SAME_NONZERO_PATTERN) {
217:       if (!dir->reuseordering) {
218:         if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
219:         if (dir->col) {ISDestroy(dir->col);}
220:         MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
221:         PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
222:         if (flg) {
223:           PetscReal tol = 1.e-10;
224:           PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
225:           MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
226:         }
227:         if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
228:       }
229:       MatDestroy(dir->fact);
230:       MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
231:       MatGetInfo(dir->fact,MAT_LOCAL,&info);
232:       dir->actualfill = info.fill_ratio_needed;
233:       PetscLogObjectParent(pc,dir->fact);
234:     }
235:     MatLUFactorNumeric(pc->pmat,&dir->fact);
236:   }
237:   return(0);
238: }

242: static PetscErrorCode PCDestroy_LU(PC pc)
243: {
244:   PC_LU          *dir = (PC_LU*)pc->data;

248:   if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
249:   if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
250:   if (dir->col) {ISDestroy(dir->col);}
251:   PetscStrfree(dir->ordering);
252:   PetscFree(dir);
253:   return(0);
254: }

258: static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y)
259: {
260:   PC_LU          *dir = (PC_LU*)pc->data;

264:   if (dir->inplace) {MatSolve(pc->pmat,x,y);}
265:   else              {MatSolve(dir->fact,x,y);}
266:   return(0);
267: }

271: static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y)
272: {
273:   PC_LU          *dir = (PC_LU*)pc->data;

277:   if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
278:   else              {MatSolveTranspose(dir->fact,x,y);}
279:   return(0);
280: }

282: /* -----------------------------------------------------------------------------------*/

287: PetscErrorCode PCLUSetFill_LU(PC pc,PetscReal fill)
288: {
289:   PC_LU *dir;

292:   dir = (PC_LU*)pc->data;
293:   dir->info.fill = fill;
294:   return(0);
295: }

301: PetscErrorCode PCLUSetDamping_LU(PC pc,PetscReal damping)
302: {
303:   PC_LU *dir;

306:   dir = (PC_LU*)pc->data;
307:   if (damping == (PetscReal) PETSC_DECIDE) {
308:     dir->info.damping = 1.e-12;
309:   } else {
310:     dir->info.damping = damping;
311:   }
312:   return(0);
313: }

319: PetscErrorCode PCLUSetUseInPlace_LU(PC pc)
320: {
321:   PC_LU *dir;

324:   dir = (PC_LU*)pc->data;
325:   dir->inplace = PETSC_TRUE;
326:   return(0);
327: }

333: PetscErrorCode PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
334: {
335:   PC_LU          *dir = (PC_LU*)pc->data;

339:   PetscStrfree(dir->ordering);
340:   PetscStrallocpy(ordering,&dir->ordering);
341:   return(0);
342: }

348: PetscErrorCode PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
349: {
350:   PC_LU *dir = (PC_LU*)pc->data;

353:   if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
354:   dir->info.dtcol = dtcol;
355:   return(0);
356: }

362: PetscErrorCode PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
363: {
364:   PC_LU *dir = (PC_LU*)pc->data;

367:   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
368:   return(0);
369: }

372: /* -----------------------------------------------------------------------------------*/

376: /*@
377:    PCLUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero

379:    Collective on PC
380:    
381:    Input Parameters:
382: +  pc - the preconditioner context
383: -  zero - all pivots smaller than this will be considered zero

385:    Options Database Key:
386: .  -pc_lu_zeropivot <zero> - Sets the zero pivot size

388:    Level: intermediate

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

392: .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot()
393: @*/
394: PetscErrorCode PCLUSetZeroPivot(PC pc,PetscReal zero)
395: {
396:   PetscErrorCode ierr,(*f)(PC,PetscReal);

400:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetZeroPivot_C",(void (**)(void))&f);
401:   if (f) {
402:     (*f)(pc,zero);
403:   }
404:   return(0);
405: }

409: /*@
410:    PCLUSetShift - specify whether to use Manteuffel shifting of LU.
411:    If an LU factorisation breaks down because of nonpositive pivots,
412:    adding sufficient identity to the diagonal will remedy this.
413:    Setting this causes a bisection method to find the minimum shift that
414:    will lead to a well-defined LU.

416:    Input parameters:
417: +  pc - the preconditioner context
418: -  shifting - PETSC_TRUE to set shift else PETSC_FALSE

420:    Options Database Key:
421: .  -pc_lu_shift - Activate PCLUSetShift()

423:    Level: intermediate

425: .keywords: PC, indefinite, factorization

427: .seealso: PCLUSetDamping(), PCILUSetShift()
428: @*/
429: PetscErrorCode PCLUSetShift(PC pc,PetscTruth shifting)
430: {
431:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

435:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetShift_C",(void (**)(void))&f);
436:   if (f) {
437:     (*f)(pc,shifting);
438:   }
439:   return(0);
440: }


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

450:    Collective on PC

452:    Input Parameters:
453: +  pc - the preconditioner context
454: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

456:    Options Database Key:
457: .  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()

459:    Level: intermediate

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

463: .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
464: @*/
465: PetscErrorCode PCLUSetReuseOrdering(PC pc,PetscTruth flag)
466: {
467:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

471:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);
472:   if (f) {
473:     (*f)(pc,flag);
474:   }
475:   return(0);
476: }

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

484:    Collective on PC

486:    Input Parameters:
487: +  pc - the preconditioner context
488: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

490:    Options Database Key:
491: .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()

493:    Level: intermediate

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

497: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
498: @*/
499: PetscErrorCode PCLUSetReuseFill(PC pc,PetscTruth flag)
500: {
501:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

505:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);
506:   if (f) {
507:     (*f)(pc,flag);
508:   }
509:   return(0);
510: }

514: /*@
515:    PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
516:    fill = number nonzeros in factor/number nonzeros in original matrix.

518:    Collective on PC
519:    
520:    Input Parameters:
521: +  pc - the preconditioner context
522: -  fill - amount of expected fill

524:    Options Database Key:
525: .  -pc_lu_fill <fill> - Sets fill amount

527:    Level: intermediate

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

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

537: .seealso: PCILUSetFill()
538: @*/
539: PetscErrorCode PCLUSetFill(PC pc,PetscReal fill)
540: {
541:   PetscErrorCode ierr,(*f)(PC,PetscReal);

545:   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
546:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);
547:   if (f) {
548:     (*f)(pc,fill);
549:   }
550:   return(0);
551: }

555: /*@
556:    PCLUSetDamping - adds this quantity to the diagonal of the matrix during the 
557:      LU numerical factorization

559:    Collective on PC
560:    
561:    Input Parameters:
562: +  pc - the preconditioner context
563: -  damping - amount of damping  (use PETSC_DECIDE for default of 1.e-12)

565:    Options Database Key:
566: .  -pc_lu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default

568:    Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
569:          pivot, then the damping is doubled until this is alleviated.

571:    Level: intermediate

573: .keywords: PC, set, factorization, direct, fill
574: .seealso: PCILUSetFill(), PCILUSetDamp()
575: @*/
576: PetscErrorCode PCLUSetDamping(PC pc,PetscReal damping)
577: {
578:   PetscErrorCode ierr,(*f)(PC,PetscReal);

582:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);
583:   if (f) {
584:     (*f)(pc,damping);
585:   }
586:   return(0);
587: }

591: /*@
592:    PCLUSetUseInPlace - Tells the system to do an in-place factorization.
593:    For dense matrices, this enables the solution of much larger problems. 
594:    For sparse matrices the factorization cannot be done truly in-place 
595:    so this does not save memory during the factorization, but after the matrix
596:    is factored, the original unfactored matrix is freed, thus recovering that
597:    space.

599:    Collective on PC

601:    Input Parameters:
602: .  pc - the preconditioner context

604:    Options Database Key:
605: .  -pc_lu_in_place - Activates in-place factorization

607:    Notes:
608:    PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when 
609:    a different matrix is provided for the multiply and the preconditioner in 
610:    a call to KSPSetOperators().
611:    This is because the Krylov space methods require an application of the 
612:    matrix multiplication, which is not possible here because the matrix has 
613:    been factored in-place, replacing the original matrix.

615:    Level: intermediate

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

619: .seealso: PCILUSetUseInPlace()
620: @*/
621: PetscErrorCode PCLUSetUseInPlace(PC pc)
622: {
623:   PetscErrorCode ierr,(*f)(PC);

627:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);
628:   if (f) {
629:     (*f)(pc);
630:   }
631:   return(0);
632: }

636: /*@C
637:     PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to 
638:     be used in the LU factorization.

640:     Collective on PC

642:     Input Parameters:
643: +   pc - the preconditioner context
644: -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM

646:     Options Database Key:
647: .   -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine

649:     Level: intermediate

651:     Notes: nested dissection is used by default

653: .seealso: PCILUSetMatOrdering()
654: @*/
655: PetscErrorCode PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
656: {
657:   PetscErrorCode ierr,(*f)(PC,MatOrderingType);

660:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);
661:   if (f) {
662:     (*f)(pc,ordering);
663:   }
664:   return(0);
665: }

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

674:     Collective on PC

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

680:     Options Database Key:
681: .   -pc_lu_pivoting <dtcol>

683:     Level: intermediate

685: .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks()
686: @*/
687: PetscErrorCode PCLUSetPivoting(PC pc,PetscReal dtcol)
688: {
689:   PetscErrorCode ierr,(*f)(PC,PetscReal);

692:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);
693:   if (f) {
694:     (*f)(pc,dtcol);
695:   }
696:   return(0);
697: }

701: /*@
702:     PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block
703:       with BAIJ or SBAIJ matrices

705:     Collective on PC

707:     Input Parameters:
708: +   pc - the preconditioner context
709: -   pivot - PETSC_TRUE or PETSC_FALSE

711:     Options Database Key:
712: .   -pc_lu_pivot_in_blocks <true,false>

714:     Level: intermediate

716: .seealso: PCILUSetMatOrdering(), PCLUSetPivoting()
717: @*/
718: PetscErrorCode PCLUSetPivotInBlocks(PC pc,PetscTruth pivot)
719: {
720:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

723:   PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);
724:   if (f) {
725:     (*f)(pc,pivot);
726:   }
727:   return(0);
728: }

730: /* ------------------------------------------------------------------------ */

732: /*MC
733:    PCLU - Uses a direct solver, based on LU factorization, as a preconditioner

735:    Options Database Keys:
736: +  -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
737: .  -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
738: .  -pc_lu_fill <fill> - Sets fill amount
739: .  -pc_lu_damping <damping> - Sets damping amount
740: .  -pc_lu_in_place - Activates in-place factorization
741: .  -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
742: .  -pc_lu_pivoting <dtcol>
743: -  -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
744:                                          stability of factorization.

746:    Notes: Not all options work for all matrix formats
747:           Run with -help to see additional options for particular matrix formats or factorization
748:           algorithms

750:    Level: beginner

752:    Concepts: LU factorization, direct solver

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

758: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
759:            PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(),
760:            PCLUSetFill(), PCLUSetDamping(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCLUSetPivoting(),
761:            PCLUSetPivotingInBlocks()
762: M*/

767: PetscErrorCode PCCreate_LU(PC pc)
768: {
770:   PetscMPIInt    size;
771:   PC_LU          *dir;

774:   PetscNew(PC_LU,&dir);
775:   PetscLogObjectMemory(pc,sizeof(PC_LU));

777:   MatFactorInfoInitialize(&dir->info);
778:   dir->fact                = 0;
779:   dir->inplace             = PETSC_FALSE;
780:   dir->info.fill           = 5.0;
781:   dir->info.dtcol          = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
782:   dir->info.damping        = 0.0;
783:   dir->info.zeropivot      = 1.e-12;
784:   dir->info.pivotinblocks  = 1.0;
785:   dir->info.shift          = PETSC_FALSE;
786:   dir->info.shift_fraction = 0.0;
787:   dir->col                 = 0;
788:   dir->row                 = 0;
789:   MPI_Comm_size(pc->comm,&size);
790:   if (size == 1) {
791:     PetscStrallocpy(MATORDERING_ND,&dir->ordering);
792:   } else {
793:     PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
794:   }
795:   dir->reusefill        = PETSC_FALSE;
796:   dir->reuseordering    = PETSC_FALSE;
797:   pc->data              = (void*)dir;

799:   pc->ops->destroy           = PCDestroy_LU;
800:   pc->ops->apply             = PCApply_LU;
801:   pc->ops->applytranspose    = PCApplyTranspose_LU;
802:   pc->ops->setup             = PCSetUp_LU;
803:   pc->ops->setfromoptions    = PCSetFromOptions_LU;
804:   pc->ops->view              = PCView_LU;
805:   pc->ops->applyrichardson   = 0;
806:   pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;

808:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
809:                     PCLUSetFill_LU);
810:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
811:                     PCLUSetDamping_LU);
812:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetShift_C","PCLUSetShift_LU",
813:                     PCLUSetShift_LU);
814:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
815:                     PCLUSetUseInPlace_LU);
816:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
817:                     PCLUSetMatOrdering_LU);
818:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
819:                     PCLUSetReuseOrdering_LU);
820:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
821:                     PCLUSetReuseFill_LU);
822:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
823:                     PCLUSetPivoting_LU);
824:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU",
825:                     PCLUSetPivotInBlocks_LU);
826:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetZeroPivot_C","PCLUSetZeroPivot_LU",
827:                     PCLUSetZeroPivot_LU);
828:   return(0);
829: }