Actual source code: factor.c

  1: #define PETSCKSP_DLL

 3:  #include ../src/ksp/pc/impls/factor/factor.h

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

 10:    Collective on PC
 11:    
 12:    Input Parameters:
 13: +  pc - the preconditioner context
 14: -  zero - all pivots smaller than this will be considered zero

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

 19:    Level: intermediate

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

 23: .seealso: PCFactorSetShiftNonzero(), PCFactorSetShiftPd()
 24: @*/
 25: PetscErrorCode  PCFactorSetZeroPivot(PC pc,PetscReal zero)
 26: {
 27:   PetscErrorCode ierr,(*f)(PC,PetscReal);

 31:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",(void (**)(void))&f);
 32:   if (f) {
 33:     (*f)(pc,zero);
 34:   }
 35:   return(0);
 36: }

 40: /*@
 41:    PCFactorSetShiftNonzero - adds this quantity to the diagonal of the matrix during 
 42:      numerical factorization, thus the matrix has nonzero pivots

 44:    Collective on PC
 45:    
 46:    Input Parameters:
 47: +  pc - the preconditioner context
 48: -  shift - amount of shift

 50:    Options Database Key:
 51: .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default

 53:    Note: If 0.0 is given, then no shift is used. If a diagonal element is classified as a zero
 54:          pivot, then the shift is doubled until this is alleviated.

 56:    Level: intermediate

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

 60: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftPd()
 61: @*/
 62: PetscErrorCode  PCFactorSetShiftNonzero(PC pc,PetscReal shift)
 63: {
 64:   PetscErrorCode ierr,(*f)(PC,PetscReal);

 68:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftNonzero_C",(void (**)(void))&f);
 69:   if (f) {
 70:     (*f)(pc,shift);
 71:   }
 72:   return(0);
 73: }

 77: /*@
 78:    PCFactorSetShiftPd - specify whether to use Manteuffel shifting.
 79:    If a matrix factorisation breaks down because of nonpositive pivots,
 80:    adding sufficient identity to the diagonal will remedy this.
 81:    Setting this causes a bisection method to find the minimum shift that
 82:    will lead to a well-defined matrix factor.

 84:    Collective on PC

 86:    Input parameters:
 87: +  pc - the preconditioner context
 88: -  shifting - PETSC_TRUE to set shift else PETSC_FALSE

 90:    Options Database Key:
 91: .  -pc_factor_shift_positive_definite true or false - Activate/Deactivate PCFactorSetShiftPd(); the value
 92:    is optional with false being the default

 94:    Level: intermediate

 96: .keywords: PC, indefinite, factorization

 98: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftNonzero()
 99: @*/
100: PetscErrorCode  PCFactorSetShiftPd(PC pc,PetscTruth shift)
101: {
102:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

106:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftPd_C",(void (**)(void))&f);
107:   if (f) {
108:     (*f)(pc,shift);
109:   }
110:   return(0);
111: }

115: /*@
116:    PCFactorSetUseDropTolerance - The preconditioner will use an ILU 
117:    based on a drop tolerance.

119:    Collective on PC

121:    Input Parameters:
122: +  pc - the preconditioner context
123: .  dt - the drop tolerance, try from 1.e-10 to .1
124: .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
125: -  maxrowcount - the max number of nonzeros allowed in a row, best value
126:                  depends on the number of nonzeros in row of original matrix

128:    Options Database Key:
129: .  -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance

131:    Level: intermediate

133:     Notes:
134:       This uses the iludt() code of Saad's SPARSKIT package

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

139: .keywords: PC, levels, reordering, factorization, incomplete, ILU
140: @*/
141: PetscErrorCode  PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
142: {
143:   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);

147:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);
148:   if (f) {
149:     (*f)(pc,dt,dtcol,maxrowcount);
150:   }
151:   return(0);
152: }

156: /*@
157:    PCFactorSetLevels - Sets the number of levels of fill to use.

159:    Collective on PC

161:    Input Parameters:
162: +  pc - the preconditioner context
163: -  levels - number of levels of fill

165:    Options Database Key:
166: .  -pc_factor_levels <levels> - Sets fill level

168:    Level: intermediate

170: .keywords: PC, levels, fill, factorization, incomplete, ILU
171: @*/
172: PetscErrorCode  PCFactorSetLevels(PC pc,PetscInt levels)
173: {
174:   PetscErrorCode ierr,(*f)(PC,PetscInt);

178:   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
179:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);
180:   if (f) {
181:     (*f)(pc,levels);
182:   }
183:   return(0);
184: }

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

192:    Collective on PC

194:    Input Parameters:
195: +  pc - the preconditioner context

197:    Options Database Key:
198: .  -pc_factor_diagonal_fill

200:    Notes:
201:    Does not apply with 0 fill.

203:    Level: intermediate

205: .keywords: PC, levels, fill, factorization, incomplete, ILU
206: @*/
207: PetscErrorCode  PCFactorSetAllowDiagonalFill(PC pc)
208: {
209:   PetscErrorCode ierr,(*f)(PC);

213:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);
214:   if (f) {
215:     (*f)(pc);
216:   }
217:   return(0);
218: }

222: /*@
223:     PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot.

225:     Collective on PC

227:     Input Parameters:
228: +   pc - the preconditioner context
229: -   shift - amount to shift or PETSC_DECIDE

231:     Options Database Key:
232: .   -pc_factor_shift_in_blocks <shift>

234:     Level: intermediate

236: .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks()
237: @*/
238: PetscErrorCode  PCFactorSetShiftInBlocks(PC pc,PetscReal shift)
239: {
240:   PetscErrorCode ierr,(*f)(PC,PetscReal);

243:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);
244:   if (f) {
245:     (*f)(pc,shift);
246:   }
247:   return(0);
248: }

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

255:    Collective on PC
256:    
257:    Input Parameters:
258: +  pc - the preconditioner context
259: -  tol - diagonal entries smaller than this in absolute value are considered zero

261:    Options Database Key:
262: .  -pc_factor_nonzeros_along_diagonal

264:    Level: intermediate

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

268: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
269: @*/
270: PetscErrorCode  PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
271: {
272:   PetscErrorCode ierr,(*f)(PC,PetscReal);

276:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",(void (**)(void))&f);
277:   if (f) {
278:     (*f)(pc,rtol);
279:   }
280:   return(0);
281: }

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

288:    Collective on PC
289:    
290:    Input Parameters:
291: +  pc - the preconditioner context
292: -  stype - for example, spooles, superlu, superlu_dist

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

297:    Level: intermediate

299:    Note:
300:      By default this will use the PETSc factorization if it exists
301:      

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

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

307: @*/
308: PetscErrorCode  PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
309: {
310:   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage);

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

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

326:    Collective on PC
327:    
328:    Input Parameter:
329: .  pc - the preconditioner context

331:    Output Parameter:
332: .   stype - for example, spooles, superlu, superlu_dist

334:    Level: intermediate


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

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

341: @*/
342: PetscErrorCode  PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
343: {
344:   PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);

348:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",(void (**)(void))&f);
349:   if (f) {
350:     (*f)(pc,stype);
351:   }
352:   return(0);
353: }

357: /*@
358:    PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
359:    fill = number nonzeros in factor/number nonzeros in original matrix.

361:    Collective on PC
362:    
363:    Input Parameters:
364: +  pc - the preconditioner context
365: -  fill - amount of expected fill

367:    Options Database Key:
368: .  -pc_factor_fill <fill> - Sets fill amount

370:    Level: intermediate

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

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

380: @*/
381: PetscErrorCode  PCFactorSetFill(PC pc,PetscReal fill)
382: {
383:   PetscErrorCode ierr,(*f)(PC,PetscReal);

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

397: /*@
398:    PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
399:    For dense matrices, this enables the solution of much larger problems. 
400:    For sparse matrices the factorization cannot be done truly in-place 
401:    so this does not save memory during the factorization, but after the matrix
402:    is factored, the original unfactored matrix is freed, thus recovering that
403:    space.

405:    Collective on PC

407:    Input Parameters:
408: .  pc - the preconditioner context

410:    Options Database Key:
411: .  -pc_factor_in_place - Activates in-place factorization

413:    Notes:
414:    PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 
415:    a different matrix is provided for the multiply and the preconditioner in 
416:    a call to KSPSetOperators().
417:    This is because the Krylov space methods require an application of the 
418:    matrix multiplication, which is not possible here because the matrix has 
419:    been factored in-place, replacing the original matrix.

421:    Level: intermediate

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

425: .seealso: PCILUSetUseInPlace()
426: @*/
427: PetscErrorCode  PCFactorSetUseInPlace(PC pc)
428: {
429:   PetscErrorCode ierr,(*f)(PC);

433:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",(void (**)(void))&f);
434:   if (f) {
435:     (*f)(pc);
436:   }
437:   return(0);
438: }

442: /*@C
443:     PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 
444:     be used in the LU factorization.

446:     Collective on PC

448:     Input Parameters:
449: +   pc - the preconditioner context
450: -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM

452:     Options Database Key:
453: .   -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine

455:     Level: intermediate

457:     Notes: nested dissection is used by default

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

463: @*/
464: PetscErrorCode  PCFactorSetMatOrderingType(PC pc,const MatOrderingType ordering)
465: {
466:   PetscErrorCode ierr,(*f)(PC,const MatOrderingType);

469:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",(void (**)(void))&f);
470:   if (f) {
471:     (*f)(pc,ordering);
472:   }
473:   return(0);
474: }

478: /*@
479:     PCFactorSetPivoting - Determines when pivoting is done during LU. 
480:       For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
481:       it is never done. For the Matlab and SuperLU factorization this is used.

483:     Collective on PC

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

489:     Options Database Key:
490: .   -pc_factor_pivoting <dtcol>

492:     Level: intermediate

494: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
495: @*/
496: PetscErrorCode  PCFactorSetPivoting(PC pc,PetscReal dtcol)
497: {
498:   PetscErrorCode ierr,(*f)(PC,PetscReal);

501:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivoting_C",(void (**)(void))&f);
502:   if (f) {
503:     (*f)(pc,dtcol);
504:   }
505:   return(0);
506: }

510: /*@
511:     PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
512:       with BAIJ or SBAIJ matrices

514:     Collective on PC

516:     Input Parameters:
517: +   pc - the preconditioner context
518: -   pivot - PETSC_TRUE or PETSC_FALSE

520:     Options Database Key:
521: .   -pc_factor_pivot_in_blocks <true,false>

523:     Level: intermediate

525: .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting()
526: @*/
527: PetscErrorCode  PCFactorSetPivotInBlocks(PC pc,PetscTruth pivot)
528: {
529:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

532:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",(void (**)(void))&f);
533:   if (f) {
534:     (*f)(pc,pivot);
535:   }
536:   return(0);
537: }

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

545:    Collective on PC

547:    Input Parameters:
548: +  pc - the preconditioner context
549: -  flag - PETSC_TRUE to reuse else PETSC_FALSE

551:    Options Database Key:
552: .  -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()

554:    Level: intermediate

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

558: .seealso: PCFactorSetReuseOrdering()
559: @*/
560: PetscErrorCode  PCFactorSetReuseFill(PC pc,PetscTruth flag)
561: {
562:   PetscErrorCode ierr,(*f)(PC,PetscTruth);

566:   PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseFill_C",(void (**)(void))&f);
567:   if (f) {
568:     (*f)(pc,flag);
569:   }
570:   return(0);
571: }