Actual source code: factor.c

petsc-3.3-p7 2013-05-11
  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
 36:    
 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
 68:    
 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
100:    
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:    PCFactorSetLevels - Sets the number of levels of fill to use.

167:    Logically Collective on PC

169:    Input Parameters:
170: +  pc - the preconditioner context
171: -  levels - number of levels of fill

173:    Options Database Key:
174: .  -pc_factor_levels <levels> - Sets fill level

176:    Level: intermediate

178: .keywords: PC, levels, fill, factorization, incomplete, ILU
179: @*/
180: PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
181: {

186:   if (levels < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
188:   PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
189:   return(0);
190: }

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

198:    Logically Collective on PC

200:    Input Parameters:
201: +  pc - the preconditioner context

203:    Options Database Key:
204: .  -pc_factor_diagonal_fill

206:    Notes:
207:    Does not apply with 0 fill.

209:    Level: intermediate

211: .keywords: PC, levels, fill, factorization, incomplete, ILU
212: @*/
213: PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc)
214: {

219:   PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC),(pc));
220:   return(0);
221: }

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

228:    Logically Collective on PC
229:    
230:    Input Parameters:
231: +  pc - the preconditioner context
232: -  tol - diagonal entries smaller than this in absolute value are considered zero

234:    Options Database Key:
235: .  -pc_factor_nonzeros_along_diagonal

237:    Level: intermediate

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

241: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
242: @*/
243: PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
244: {

250:   PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
251:   return(0);
252: }

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

259:    Logically Collective on PC
260:    
261:    Input Parameters:
262: +  pc - the preconditioner context
263: -  stype - for example, spooles, superlu, superlu_dist

265:    Options Database Key:
266: .  -pc_factor_mat_solver_package <stype> - spooles, petsc, superlu, superlu_dist, mumps

268:    Level: intermediate

270:    Note:
271:      By default this will use the PETSc factorization if it exists
272:      

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

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

278: @*/
279: PetscErrorCode  PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
280: {

285:   PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));
286:   return(0);
287: }

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

294:    Not Collective
295:    
296:    Input Parameter:
297: .  pc - the preconditioner context

299:    Output Parameter:
300: .   stype - for example, spooles, superlu, superlu_dist (PETSC_NULL if the PC does not have a solver package)

302:    Level: intermediate


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

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

309: @*/
310: PetscErrorCode  PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
311: {
312:   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);

316:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void(**)(void))&f);
317:   if (f) {
318:     (*f)(pc,stype);
319:   } else {
320:     *stype = PETSC_NULL;
321:   }
322:   return(0);
323: }

327: /*@
328:    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
329:    fill = number nonzeros in factor/number nonzeros in original matrix.

331:    Not Collective, each process can expect a different amount of fill
332:    
333:    Input Parameters:
334: +  pc - the preconditioner context
335: -  fill - amount of expected fill

337:    Options Database Key:
338: .  -pc_factor_fill <fill> - Sets fill amount

340:    Level: intermediate

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

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

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

353: @*/
354: PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
355: {

360:   if (fill < 1.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
361:   PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
362:   return(0);
363: }

367: /*@
368:    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
369:    For dense matrices, this enables the solution of much larger problems. 
370:    For sparse matrices the factorization cannot be done truly in-place 
371:    so this does not save memory during the factorization, but after the matrix
372:    is factored, the original unfactored matrix is freed, thus recovering that
373:    space.

375:    Logically Collective on PC

377:    Input Parameters:
378: .  pc - the preconditioner context

380:    Options Database Key:
381: .  -pc_factor_in_place - Activates in-place factorization

383:    Notes:
384:    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 
385:    a different matrix is provided for the multiply and the preconditioner in 
386:    a call to KSPSetOperators().
387:    This is because the Krylov space methods require an application of the 
388:    matrix multiplication, which is not possible here because the matrix has 
389:    been factored in-place, replacing the original matrix.

391:    Level: intermediate

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

395: .seealso: PCILUSetUseInPlace()
396: @*/
397: PetscErrorCode  PCFactorSetUseInPlace(PC pc)
398: {

403:   PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC),(pc));
404:   return(0);
405: }

409: /*@C
410:     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 
411:     be used in the LU factorization.

413:     Logically Collective on PC

415:     Input Parameters:
416: +   pc - the preconditioner context
417: -   ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM

419:     Options Database Key:
420: .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

422:     Level: intermediate

424:     Notes: nested dissection is used by default

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

430: @*/
431: PetscErrorCode  PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
432: {

437:   PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,const MatOrderingType),(pc,ordering));
438:   return(0);
439: }

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

448:     Logically Collective on PC

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

454:     Options Database Key:
455: .   -pc_factor_pivoting <dtcol>

457:     Level: intermediate

459: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
460: @*/
461: PetscErrorCode  PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
462: {

468:   PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
469:   return(0);
470: }

474: /*@
475:     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
476:       with BAIJ or SBAIJ matrices

478:     Logically Collective on PC

480:     Input Parameters:
481: +   pc - the preconditioner context
482: -   pivot - PETSC_TRUE or PETSC_FALSE

484:     Options Database Key:
485: .   -pc_factor_pivot_in_blocks <true,false>

487:     Level: intermediate

489: .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
490: @*/
491: PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscBool  pivot)
492: {

498:   PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
499:   return(0);
500: }

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

508:    Logically Collective on PC

510:    Input Parameters:
511: +  pc - the preconditioner context
512: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

514:    Options Database Key:
515: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()

517:    Level: intermediate

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

521: .seealso: PCFactorSetReuseOrdering()
522: @*/
523: PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscBool  flag)
524: {

530:   PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
531:   return(0);
532: }