Actual source code: factor.c
petsc-3.6.4 2016-04-12
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: PCFactorGetLevels - Gets the number of levels of fill to use.
167: Logically Collective on PC
169: Input Parameters:
170: . pc - the preconditioner context
172: Output Parameter:
173: . levels - number of levels of fill
175: Level: intermediate
177: .keywords: PC, levels, fill, factorization, incomplete, ILU
178: @*/
179: PetscErrorCode PCFactorGetLevels(PC pc,PetscInt *levels)
180: {
185: PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
186: return(0);
187: }
191: /*@
192: PCFactorSetLevels - Sets the number of levels of fill to use.
194: Logically Collective on PC
196: Input Parameters:
197: + pc - the preconditioner context
198: - levels - number of levels of fill
200: Options Database Key:
201: . -pc_factor_levels <levels> - Sets fill level
203: Level: intermediate
205: .keywords: PC, levels, fill, factorization, incomplete, ILU
206: @*/
207: PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels)
208: {
213: if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
215: PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
216: return(0);
217: }
221: /*@
222: PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
223: treated as level 0 fill even if there is no non-zero location.
225: Logically Collective on PC
227: Input Parameters:
228: + pc - the preconditioner context
229: - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
231: Options Database Key:
232: . -pc_factor_diagonal_fill
234: Notes:
235: Does not apply with 0 fill.
237: Level: intermediate
239: .keywords: PC, levels, fill, factorization, incomplete, ILU
241: .seealso: PCFactorGetAllowDiagonalFill()
243: @*/
244: PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
245: {
250: PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
251: return(0);
252: }
256: /*@
257: PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
258: treated as level 0 fill even if there is no non-zero location.
260: Logically Collective on PC
262: Input Parameter:
263: . pc - the preconditioner context
265: Output Parameter:
266: . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
268: Options Database Key:
269: . -pc_factor_diagonal_fill
271: Notes:
272: Does not apply with 0 fill.
274: Level: intermediate
276: .keywords: PC, levels, fill, factorization, incomplete, ILU
278: .seealso: PCFactorSetAllowDiagonalFill()
280: @*/
281: PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
282: {
287: PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
288: return(0);
289: }
293: /*@
294: PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
296: Logically Collective on PC
298: Input Parameters:
299: + pc - the preconditioner context
300: - tol - diagonal entries smaller than this in absolute value are considered zero
302: Options Database Key:
303: . -pc_factor_nonzeros_along_diagonal <tol>
305: Level: intermediate
307: .keywords: PC, set, factorization, direct, fill
309: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
310: @*/
311: PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
312: {
318: PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
319: return(0);
320: }
324: /*@C
325: PCFactorSetMatSolverPackage - sets the software that is used to perform the factorization
327: Logically Collective on PC
329: Input Parameters:
330: + pc - the preconditioner context
331: - stype - for example, superlu, superlu_dist
333: Options Database Key:
334: . -pc_factor_mat_solver_package <stype> - petsc, superlu, superlu_dist, mumps, cusparse
336: Level: intermediate
338: Note:
339: By default this will use the PETSc factorization if it exists
342: .keywords: PC, set, factorization, direct, fill
344: .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
346: @*/
347: PetscErrorCode PCFactorSetMatSolverPackage(PC pc,const MatSolverPackage stype)
348: {
353: PetscTryMethod(pc,"PCFactorSetMatSolverPackage_C",(PC,const MatSolverPackage),(pc,stype));
354: return(0);
355: }
359: /*@C
360: PCFactorGetMatSolverPackage - gets the software that is used to perform the factorization
362: Not Collective
364: Input Parameter:
365: . pc - the preconditioner context
367: Output Parameter:
368: . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
370: Level: intermediate
373: .keywords: PC, set, factorization, direct, fill
375: .seealso: MatGetFactor(), MatSolverPackage, PCFactorGetMatSolverPackage()
377: @*/
378: PetscErrorCode PCFactorGetMatSolverPackage(PC pc,const MatSolverPackage *stype)
379: {
380: PetscErrorCode ierr,(*f)(PC,const MatSolverPackage*);
384: PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",&f);
385: if (f) {
386: (*f)(pc,stype);
387: } else {
388: *stype = NULL;
389: }
390: return(0);
391: }
395: /*@
396: PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
397: fill = number nonzeros in factor/number nonzeros in original matrix.
399: Not Collective, each process can expect a different amount of fill
401: Input Parameters:
402: + pc - the preconditioner context
403: - fill - amount of expected fill
405: Options Database Key:
406: . -pc_factor_fill <fill> - Sets fill amount
408: Level: intermediate
410: Note:
411: For sparse matrix factorizations it is difficult to predict how much
412: fill to expect. By running with the option -info PETSc will print the
413: actual amount of fill used; allowing you to set the value accurately for
414: future runs. Default PETSc uses a value of 5.0
416: This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
419: .keywords: PC, set, factorization, direct, fill
421: @*/
422: PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill)
423: {
428: if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
429: PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
430: return(0);
431: }
435: /*@
436: PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
437: For dense matrices, this enables the solution of much larger problems.
438: For sparse matrices the factorization cannot be done truly in-place
439: so this does not save memory during the factorization, but after the matrix
440: is factored, the original unfactored matrix is freed, thus recovering that
441: space.
443: Logically Collective on PC
445: Input Parameters:
446: + pc - the preconditioner context
447: - flg - PETSC_TRUE to enable, PETSC_FALSE to disable
449: Options Database Key:
450: . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
452: Notes:
453: PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
454: a different matrix is provided for the multiply and the preconditioner in
455: a call to KSPSetOperators().
456: This is because the Krylov space methods require an application of the
457: matrix multiplication, which is not possible here because the matrix has
458: been factored in-place, replacing the original matrix.
460: Level: intermediate
462: .keywords: PC, set, factorization, direct, inplace, in-place, LU
464: .seealso: PCFactorGetUseInPlace()
465: @*/
466: PetscErrorCode PCFactorSetUseInPlace(PC pc,PetscBool flg)
467: {
472: PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
473: return(0);
474: }
478: /*@
479: PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
481: Logically Collective on PC
483: Input Parameter:
484: . pc - the preconditioner context
486: Output Parameter:
487: . flg - PETSC_TRUE to enable, PETSC_FALSE to disable
489: Level: intermediate
491: .keywords: PC, set, factorization, direct, inplace, in-place, LU
493: .seealso: PCFactorSetUseInPlace()
494: @*/
495: PetscErrorCode PCFactorGetUseInPlace(PC pc,PetscBool *flg)
496: {
501: PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
502: return(0);
503: }
507: /*@C
508: PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
509: be used in the LU factorization.
511: Logically Collective on PC
513: Input Parameters:
514: + pc - the preconditioner context
515: - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
517: Options Database Key:
518: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
520: Level: intermediate
522: Notes: nested dissection is used by default
524: For Cholesky and ICC and the SBAIJ format reorderings are not available,
525: since only the upper triangular part of the matrix is stored. You can use the
526: SeqAIJ format in this case to get reorderings.
528: @*/
529: PetscErrorCode PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
530: {
535: PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));
536: return(0);
537: }
541: /*@
542: PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
543: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
544: it is never done. For the MATLAB and SuperLU factorization this is used.
546: Logically Collective on PC
548: Input Parameters:
549: + pc - the preconditioner context
550: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
552: Options Database Key:
553: . -pc_factor_pivoting <dtcol>
555: Level: intermediate
557: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
558: @*/
559: PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
560: {
566: PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
567: return(0);
568: }
572: /*@
573: PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
574: with BAIJ or SBAIJ matrices
576: Logically Collective on PC
578: Input Parameters:
579: + pc - the preconditioner context
580: - pivot - PETSC_TRUE or PETSC_FALSE
582: Options Database Key:
583: . -pc_factor_pivot_in_blocks <true,false>
585: Level: intermediate
587: .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
588: @*/
589: PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
590: {
596: PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
597: return(0);
598: }
602: /*@
603: PCFactorSetReuseFill - When matrices with same different nonzero structure are factored,
604: this causes later ones to use the fill ratio computed in the initial factorization.
606: Logically Collective on PC
608: Input Parameters:
609: + pc - the preconditioner context
610: - flag - PETSC_TRUE to reuse else PETSC_FALSE
612: Options Database Key:
613: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
615: Level: intermediate
617: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
619: .seealso: PCFactorSetReuseOrdering()
620: @*/
621: PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag)
622: {
628: PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
629: return(0);
630: }