Actual source code: factor.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <../src/ksp/pc/impls/factor/factor.h>  /*I "petscpc.h" I*/

  6: /*@
  7:     PCFactorSetUpMatSolverPackage - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
  8:        set the options for that particular factorization object.

 10:   Input Parameter:
 11: .  pc  - the preconditioner context

 13:   Notes: After you have called this function (which has to be after the KSPSetOperators() or PCSetOperators()) you can call PCFactorGetMatrix() and then set factor options on that matrix.

 15: .seealso: PCFactorSetMatSolverPackage(), PCFactorGetMatrix()

 17:   Level: intermediate

 19: @*/
 20: PetscErrorCode PCFactorSetUpMatSolverPackage(PC pc)
 21: {

 26:   PetscTryMethod(pc,"PCFactorSetUpMatSolverPackage_C",(PC),(pc));
 27:   return(0);
 28: }

 32: /*@
 33:    PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero

 35:    Logically Collective on PC

 37:    Input Parameters:
 38: +  pc - the preconditioner context
 39: -  zero - all pivots smaller than this will be considered zero

 41:    Options Database Key:
 42: .  -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot

 44:    Level: intermediate

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

 48: .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
 49: @*/
 50: PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
 51: {

 57:   PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));
 58:   return(0);
 59: }

 63: /*@
 64:    PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
 65:      numerical factorization, thus the matrix has nonzero pivots

 67:    Logically Collective on PC

 69:    Input Parameters:
 70: +  pc - the preconditioner context
 71: -  shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS

 73:    Options Database Key:
 74: .  -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types

 76:    Level: intermediate

 78: .keywords: PC, set, factorization,

 80: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
 81: @*/
 82: PetscErrorCode  PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
 83: {

 89:   PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));
 90:   return(0);
 91: }

 95: /*@
 96:    PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
 97:      numerical factorization, thus the matrix has nonzero pivots

 99:    Logically Collective on PC

101:    Input Parameters:
102: +  pc - the preconditioner context
103: -  shiftamount - amount of shift

105:    Options Database Key:
106: .  -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default

108:    Level: intermediate

110: .keywords: PC, set, factorization,

112: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
113: @*/
114: PetscErrorCode  PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
115: {

121:   PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));
122:   return(0);
123: }

127: /*
128:    PCFactorSetDropTolerance - The preconditioner will use an ILU
129:    based on a drop tolerance. (Under development)

131:    Logically Collective on PC

133:    Input Parameters:
134: +  pc - the preconditioner context
135: .  dt - the drop tolerance, try from 1.e-10 to .1
136: .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
137: -  maxrowcount - the max number of nonzeros allowed in a row, best value
138:                  depends on the number of nonzeros in row of original matrix

140:    Options Database Key:
141: .  -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance

143:    Level: intermediate

145:       There are NO default values for the 3 parameters, you must set them with reasonable values for your
146:       matrix. We don't know how to compute reasonable values.

148: .keywords: PC, levels, reordering, factorization, incomplete, ILU
149: */
150: PetscErrorCode  PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
151: {

158:   PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));
159:   return(0);
160: }

164: /*@
165:    PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot

167:    Not Collective

169:    Input Parameters:
170: .  pc - the preconditioner context

172:    Output Parameter:
173: .  pivot - the tolerance

175:    Level: intermediate


178: .seealso: PCFactorSetZeroPivot()
179: @*/
180: PetscErrorCode  PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
181: {

186:   PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));
187:   return(0);
188: }

192: /*@
193:    PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot

195:    Not Collective

197:    Input Parameters:
198: .  pc - the preconditioner context

200:    Output Parameter:
201: .  shift - how much to shift the diagonal entry

203:    Level: intermediate


206: .seealso: PCFactorSetShiftAmount(), PCFactorSetShiftType(), PCFactorGetShiftType()
207: @*/
208: PetscErrorCode  PCFactorGetShiftAmount(PC pc,PetscReal *shift)
209: {

214:   PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));
215:   return(0);
216: }

220: /*@
221:    PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected

223:    Not Collective

225:    Input Parameters:
226: .  pc - the preconditioner context

228:    Output Parameter:
229: .  type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO,  MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS

231:    Level: intermediate


234: .seealso: PCFactorSetShiftType(), MatFactorShiftType, PCFactorSetShiftAmount(), PCFactorGetShiftAmount()
235: @*/
236: PetscErrorCode  PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
237: {

242:   PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));
243:   return(0);
244: }

248: /*@
249:    PCFactorGetLevels - Gets the number of levels of fill to use.

251:    Logically Collective on PC

253:    Input Parameters:
254: .  pc - the preconditioner context

256:    Output Parameter:
257: .  levels - number of levels of fill

259:    Level: intermediate

261: .keywords: PC, levels, fill, factorization, incomplete, ILU
262: @*/
263: PetscErrorCode  PCFactorGetLevels(PC pc,PetscInt *levels)
264: {

269:   PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
270:   return(0);
271: }

275: /*@
276:    PCFactorSetLevels - Sets the number of levels of fill to use.

278:    Logically Collective on PC

280:    Input Parameters:
281: +  pc - the preconditioner context
282: -  levels - number of levels of fill

284:    Options Database Key:
285: .  -pc_factor_levels <levels> - Sets fill level

287:    Level: intermediate

289: .keywords: PC, levels, fill, factorization, incomplete, ILU
290: @*/
291: PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
292: {

297:   if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
299:   PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
300:   return(0);
301: }

305: /*@
306:    PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
307:    treated as level 0 fill even if there is no non-zero location.

309:    Logically Collective on PC

311:    Input Parameters:
312: +  pc - the preconditioner context
313: -  flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off

315:    Options Database Key:
316: .  -pc_factor_diagonal_fill

318:    Notes:
319:    Does not apply with 0 fill.

321:    Level: intermediate

323: .keywords: PC, levels, fill, factorization, incomplete, ILU

325: .seealso: PCFactorGetAllowDiagonalFill()

327: @*/
328: PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
329: {

334:   PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
335:   return(0);
336: }

340: /*@
341:    PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
342:        treated as level 0 fill even if there is no non-zero location.

344:    Logically Collective on PC

346:    Input Parameter:
347: .  pc - the preconditioner context

349:    Output Parameter:
350: .   flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off

352:    Options Database Key:
353: .  -pc_factor_diagonal_fill

355:    Notes:
356:    Does not apply with 0 fill.

358:    Level: intermediate

360: .keywords: PC, levels, fill, factorization, incomplete, ILU

362: .seealso: PCFactorSetAllowDiagonalFill()

364: @*/
365: PetscErrorCode  PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
366: {

371:   PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
372:   return(0);
373: }

377: /*@
378:    PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal

380:    Logically Collective on PC

382:    Input Parameters:
383: +  pc - the preconditioner context
384: -  tol - diagonal entries smaller than this in absolute value are considered zero

386:    Options Database Key:
387: .  -pc_factor_nonzeros_along_diagonal <tol>

389:    Level: intermediate

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

393: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
394: @*/
395: PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
396: {

402:   PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
403:   return(0);
404: }

408: /*@C
409:    PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization

411:    Logically Collective on PC

413:    Input Parameters:
414: +  pc - the preconditioner context
415: -  stype - for example, superlu, superlu_dist

417:    Options Database Key:
418: .  -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse

420:    Level: intermediate

422:    Note:
423:      By default this will use the PETSc factorization if it exists


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

428: .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()

430: @*/
431: PetscErrorCode  PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
432: {

437:   PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));
438:   return(0);
439: }

443: /*@C
444:    PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization

446:    Not Collective

448:    Input Parameter:
449: .  pc - the preconditioner context

451:    Output Parameter:
452: .   stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)

454:    Level: intermediate


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

459: .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()

461: @*/
462: PetscErrorCode  PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
463: {
464:   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);

468:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);
469:   if (f) {
470:     (*f)(pc,stype);
471:   } else {
472:     *stype = NULL;
473:   }
474:   return(0);
475: }

479: /*@
480:    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
481:    fill = number nonzeros in factor/number nonzeros in original matrix.

483:    Not Collective, each process can expect a different amount of fill

485:    Input Parameters:
486: +  pc - the preconditioner context
487: -  fill - amount of expected fill

489:    Options Database Key:
490: .  -pc_factor_fill <fill> - Sets fill amount

492:    Level: intermediate

494:    Note:
495:    For sparse matrix factorizations it is difficult to predict how much
496:    fill to expect. By running with the option -info PETSc will print the
497:    actual amount of fill used; allowing you to set the value accurately for
498:    future runs. Default PETSc uses a value of 5.0

500:    This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.


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

505: @*/
506: PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
507: {

512:   if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
513:   PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
514:   return(0);
515: }

519: /*@
520:    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
521:    For dense matrices, this enables the solution of much larger problems.
522:    For sparse matrices the factorization cannot be done truly in-place
523:    so this does not save memory during the factorization, but after the matrix
524:    is factored, the original unfactored matrix is freed, thus recovering that
525:    space.

527:    Logically Collective on PC

529:    Input Parameters:
530: +  pc - the preconditioner context
531: -  flg - PETSC_TRUE to enable, PETSC_FALSE to disable

533:    Options Database Key:
534: .  -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization

536:    Notes:
537:    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
538:    a different matrix is provided for the multiply and the preconditioner in
539:    a call to KSPSetOperators().
540:    This is because the Krylov space methods require an application of the
541:    matrix multiplication, which is not possible here because the matrix has
542:    been factored in-place, replacing the original matrix.

544:    Level: intermediate

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

548: .seealso: PCFactorGetUseInPlace()
549: @*/
550: PetscErrorCode  PCFactorSetUseInPlace(PC pc,PetscBool flg)
551: {

556:   PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
557:   return(0);
558: }

562: /*@
563:    PCFactorGetUseInPlace - Determines if an in-place factorization is being used.

565:    Logically Collective on PC

567:    Input Parameter:
568: .  pc - the preconditioner context

570:    Output Parameter:
571: .  flg - PETSC_TRUE to enable, PETSC_FALSE to disable

573:    Level: intermediate

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

577: .seealso: PCFactorSetUseInPlace()
578: @*/
579: PetscErrorCode  PCFactorGetUseInPlace(PC pc,PetscBool *flg)
580: {

585:   PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
586:   return(0);
587: }

591: /*@C
592:     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
593:     be used in the LU factorization.

595:     Logically Collective on PC

597:     Input Parameters:
598: +   pc - the preconditioner context
599: -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM

601:     Options Database Key:
602: .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

604:     Level: intermediate

606:     Notes: nested dissection is used by default

608:     For Cholesky and ICC and the SBAIJ format reorderings are not available,
609:     since only the upper triangular part of the matrix is stored. You can use the
610:     SeqAIJ format in this case to get reorderings.

612: @*/
613: PetscErrorCode  PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
614: {

619:   PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));
620:   return(0);
621: }

625: /*@
626:     PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
627:       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
628:       it is never done. For the MATLAB and SuperLU factorization this is used.

630:     Logically Collective on PC

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

636:     Options Database Key:
637: .   -pc_factor_pivoting <dtcol>

639:     Level: intermediate

641: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
642: @*/
643: PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
644: {

650:   PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
651:   return(0);
652: }

656: /*@
657:     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
658:       with BAIJ or SBAIJ matrices

660:     Logically Collective on PC

662:     Input Parameters:
663: +   pc - the preconditioner context
664: -   pivot - PETSC_TRUE or PETSC_FALSE

666:     Options Database Key:
667: .   -pc_factor_pivot_in_blocks <true,false>

669:     Level: intermediate

671: .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
672: @*/
673: PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
674: {

680:   PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
681:   return(0);
682: }

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

690:    Logically Collective on PC

692:    Input Parameters:
693: +  pc - the preconditioner context
694: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

696:    Options Database Key:
697: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()

699:    Level: intermediate

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

703: .seealso: PCFactorSetReuseOrdering()
704: @*/
705: PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool flag)
706: {

712:   PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
713:   return(0);
714: }