Actual source code: factor.c
petsc-3.11.4 2019-09-28
2: #include <../src/ksp/pc/impls/factor/factor.h>
4: static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag)
5: {
6: PC_Factor *lu = (PC_Factor*)pc->data;
9: lu->reuseordering = flag;
10: return(0);
11: }
13: static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc,PetscBool flag)
14: {
15: PC_Factor *lu = (PC_Factor*)pc->data;
18: lu->reusefill = flag;
19: return(0);
20: }
22: static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg)
23: {
24: PC_Factor *dir = (PC_Factor*)pc->data;
27: dir->inplace = flg;
28: return(0);
29: }
31: static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc,PetscBool *flg)
32: {
33: PC_Factor *dir = (PC_Factor*)pc->data;
36: *flg = dir->inplace;
37: return(0);
38: }
40: /*@
41: PCFactorSetUpMatSolverType - Can be called after KSPSetOperators() or PCSetOperators(), causes MatGetFactor() to be called so then one may
42: set the options for that particular factorization object.
44: Input Parameter:
45: . pc - the preconditioner context
47: Notes:
48: 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.
50: .seealso: PCFactorSetMatSolverType(), PCFactorGetMatrix()
52: Level: intermediate
54: @*/
55: PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
56: {
61: PetscTryMethod(pc,"PCFactorSetUpMatSolverType_C",(PC),(pc));
62: return(0);
63: }
65: /*@
66: PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
68: Logically Collective on PC
70: Input Parameters:
71: + pc - the preconditioner context
72: - zero - all pivots smaller than this will be considered zero
74: Options Database Key:
75: . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
77: Level: intermediate
79: .keywords: PC, set, factorization, direct, fill
81: .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
82: @*/
83: PetscErrorCode PCFactorSetZeroPivot(PC pc,PetscReal zero)
84: {
90: PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));
91: return(0);
92: }
94: /*@
95: PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
96: numerical factorization, thus the matrix has nonzero pivots
98: Logically Collective on PC
100: Input Parameters:
101: + pc - the preconditioner context
102: - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
104: Options Database Key:
105: . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types
107: Level: intermediate
109: .keywords: PC, set, factorization,
111: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
112: @*/
113: PetscErrorCode PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
114: {
120: PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));
121: return(0);
122: }
124: /*@
125: PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
126: numerical factorization, thus the matrix has nonzero pivots
128: Logically Collective on PC
130: Input Parameters:
131: + pc - the preconditioner context
132: - shiftamount - amount of shift
134: Options Database Key:
135: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
137: Level: intermediate
139: .keywords: PC, set, factorization,
141: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
142: @*/
143: PetscErrorCode PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
144: {
150: PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));
151: return(0);
152: }
154: /*
155: PCFactorSetDropTolerance - The preconditioner will use an ILU
156: based on a drop tolerance. (Under development)
158: Logically Collective on PC
160: Input Parameters:
161: + pc - the preconditioner context
162: . dt - the drop tolerance, try from 1.e-10 to .1
163: . dtcol - tolerance for column pivot, good values [0.1 to 0.01]
164: - maxrowcount - the max number of nonzeros allowed in a row, best value
165: depends on the number of nonzeros in row of original matrix
167: Options Database Key:
168: . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
170: Level: intermediate
172: There are NO default values for the 3 parameters, you must set them with reasonable values for your
173: matrix. We don't know how to compute reasonable values.
175: .keywords: PC, levels, reordering, factorization, incomplete, ILU
176: */
177: PetscErrorCode PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
178: {
185: PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));
186: return(0);
187: }
189: /*@
190: PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
192: Not Collective
194: Input Parameters:
195: . pc - the preconditioner context
197: Output Parameter:
198: . pivot - the tolerance
200: Level: intermediate
203: .seealso: PCFactorSetZeroPivot()
204: @*/
205: PetscErrorCode PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
206: {
211: PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));
212: return(0);
213: }
215: /*@
216: PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
218: Not Collective
220: Input Parameters:
221: . pc - the preconditioner context
223: Output Parameter:
224: . shift - how much to shift the diagonal entry
226: Level: intermediate
229: .seealso: PCFactorSetShiftAmount(), PCFactorSetShiftType(), PCFactorGetShiftType()
230: @*/
231: PetscErrorCode PCFactorGetShiftAmount(PC pc,PetscReal *shift)
232: {
237: PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));
238: return(0);
239: }
241: /*@
242: PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
244: Not Collective
246: Input Parameters:
247: . pc - the preconditioner context
249: Output Parameter:
250: . type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS
252: Level: intermediate
255: .seealso: PCFactorSetShiftType(), MatFactorShiftType, PCFactorSetShiftAmount(), PCFactorGetShiftAmount()
256: @*/
257: PetscErrorCode PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
258: {
263: PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));
264: return(0);
265: }
267: /*@
268: PCFactorGetLevels - Gets the number of levels of fill to use.
270: Logically Collective on PC
272: Input Parameters:
273: . pc - the preconditioner context
275: Output Parameter:
276: . levels - number of levels of fill
278: Level: intermediate
280: .keywords: PC, levels, fill, factorization, incomplete, ILU
281: @*/
282: PetscErrorCode PCFactorGetLevels(PC pc,PetscInt *levels)
283: {
288: PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
289: return(0);
290: }
292: /*@
293: PCFactorSetLevels - Sets the number of levels of fill to use.
295: Logically Collective on PC
297: Input Parameters:
298: + pc - the preconditioner context
299: - levels - number of levels of fill
301: Options Database Key:
302: . -pc_factor_levels <levels> - Sets fill level
304: Level: intermediate
306: .keywords: PC, levels, fill, factorization, incomplete, ILU
307: @*/
308: PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels)
309: {
314: if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
316: PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
317: return(0);
318: }
320: /*@
321: PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
322: treated as level 0 fill even if there is no non-zero location.
324: Logically Collective on PC
326: Input Parameters:
327: + pc - the preconditioner context
328: - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
330: Options Database Key:
331: . -pc_factor_diagonal_fill
333: Notes:
334: Does not apply with 0 fill.
336: Level: intermediate
338: .keywords: PC, levels, fill, factorization, incomplete, ILU
340: .seealso: PCFactorGetAllowDiagonalFill()
342: @*/
343: PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
344: {
349: PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
350: return(0);
351: }
353: /*@
354: PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
355: treated as level 0 fill even if there is no non-zero location.
357: Logically Collective on PC
359: Input Parameter:
360: . pc - the preconditioner context
362: Output Parameter:
363: . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
365: Options Database Key:
366: . -pc_factor_diagonal_fill
368: Notes:
369: Does not apply with 0 fill.
371: Level: intermediate
373: .keywords: PC, levels, fill, factorization, incomplete, ILU
375: .seealso: PCFactorSetAllowDiagonalFill()
377: @*/
378: PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
379: {
384: PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
385: return(0);
386: }
388: /*@
389: PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
391: Logically Collective on PC
393: Input Parameters:
394: + pc - the preconditioner context
395: - tol - diagonal entries smaller than this in absolute value are considered zero
397: Options Database Key:
398: . -pc_factor_nonzeros_along_diagonal <tol>
400: Level: intermediate
402: .keywords: PC, set, factorization, direct, fill
404: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
405: @*/
406: PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
407: {
413: PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
414: return(0);
415: }
417: /*@C
418: PCFactorSetMatSolverType - sets the software that is used to perform the factorization
420: Logically Collective on PC
422: Input Parameters:
423: + pc - the preconditioner context
424: - stype - for example, superlu, superlu_dist
426: Options Database Key:
427: . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
429: Level: intermediate
431: Note:
432: By default this will use the PETSc factorization if it exists
435: .keywords: PC, set, factorization, direct, fill
437: .seealso: MatGetFactor(), MatSolverType, PCFactorGetMatSolverType()
439: @*/
440: PetscErrorCode PCFactorSetMatSolverType(PC pc,MatSolverType stype)
441: {
446: PetscTryMethod(pc,"PCFactorSetMatSolverType_C",(PC,MatSolverType),(pc,stype));
447: return(0);
448: }
450: /*@C
451: PCFactorGetMatSolverType - gets the software that is used to perform the factorization
453: Not Collective
455: Input Parameter:
456: . pc - the preconditioner context
458: Output Parameter:
459: . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
461: Level: intermediate
464: .keywords: PC, set, factorization, direct, fill
466: .seealso: MatGetFactor(), MatSolverType, PCFactorGetMatSolverType()
468: @*/
469: PetscErrorCode PCFactorGetMatSolverType(PC pc,MatSolverType *stype)
470: {
471: PetscErrorCode ierr,(*f)(PC,MatSolverType*);
475: PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",&f);
476: if (f) {
477: (*f)(pc,stype);
478: } else {
479: *stype = NULL;
480: }
481: return(0);
482: }
484: /*@
485: PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
486: fill = number nonzeros in factor/number nonzeros in original matrix.
488: Not Collective, each process can expect a different amount of fill
490: Input Parameters:
491: + pc - the preconditioner context
492: - fill - amount of expected fill
494: Options Database Key:
495: . -pc_factor_fill <fill> - Sets fill amount
497: Level: intermediate
499: Note:
500: For sparse matrix factorizations it is difficult to predict how much
501: fill to expect. By running with the option -info PETSc will print the
502: actual amount of fill used; allowing you to set the value accurately for
503: future runs. Default PETSc uses a value of 5.0
505: This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
508: .keywords: PC, set, factorization, direct, fill
510: @*/
511: PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill)
512: {
517: if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
518: PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
519: return(0);
520: }
522: /*@
523: PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
524: For dense matrices, this enables the solution of much larger problems.
525: For sparse matrices the factorization cannot be done truly in-place
526: so this does not save memory during the factorization, but after the matrix
527: is factored, the original unfactored matrix is freed, thus recovering that
528: space.
530: Logically Collective on PC
532: Input Parameters:
533: + pc - the preconditioner context
534: - flg - PETSC_TRUE to enable, PETSC_FALSE to disable
536: Options Database Key:
537: . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
539: Notes:
540: PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
541: a different matrix is provided for the multiply and the preconditioner in
542: a call to KSPSetOperators().
543: This is because the Krylov space methods require an application of the
544: matrix multiplication, which is not possible here because the matrix has
545: been factored in-place, replacing the original matrix.
547: Level: intermediate
549: .keywords: PC, set, factorization, direct, inplace, in-place, LU
551: .seealso: PCFactorGetUseInPlace()
552: @*/
553: PetscErrorCode PCFactorSetUseInPlace(PC pc,PetscBool flg)
554: {
559: PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
560: return(0);
561: }
563: /*@
564: PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
566: Logically Collective on PC
568: Input Parameter:
569: . pc - the preconditioner context
571: Output Parameter:
572: . flg - PETSC_TRUE to enable, PETSC_FALSE to disable
574: Level: intermediate
576: .keywords: PC, set, factorization, direct, inplace, in-place, LU
578: .seealso: PCFactorSetUseInPlace()
579: @*/
580: PetscErrorCode PCFactorGetUseInPlace(PC pc,PetscBool *flg)
581: {
586: PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
587: return(0);
588: }
590: /*@C
591: PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
592: be used in the LU factorization.
594: Logically Collective on PC
596: Input Parameters:
597: + pc - the preconditioner context
598: - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
600: Options Database Key:
601: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
603: Level: intermediate
605: Notes:
606: 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: }
623: /*@
624: PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
625: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
626: it is never done. For the MATLAB and SuperLU factorization this is used.
628: Logically Collective on PC
630: Input Parameters:
631: + pc - the preconditioner context
632: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
634: Options Database Key:
635: . -pc_factor_pivoting <dtcol>
637: Level: intermediate
639: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
640: @*/
641: PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
642: {
648: PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
649: return(0);
650: }
652: /*@
653: PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
654: with BAIJ or SBAIJ matrices
656: Logically Collective on PC
658: Input Parameters:
659: + pc - the preconditioner context
660: - pivot - PETSC_TRUE or PETSC_FALSE
662: Options Database Key:
663: . -pc_factor_pivot_in_blocks <true,false>
665: Level: intermediate
667: .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
668: @*/
669: PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
670: {
676: PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
677: return(0);
678: }
680: /*@
681: PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
682: this causes later ones to use the fill ratio computed in the initial factorization.
684: Logically Collective on PC
686: Input Parameters:
687: + pc - the preconditioner context
688: - flag - PETSC_TRUE to reuse else PETSC_FALSE
690: Options Database Key:
691: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
693: Level: intermediate
695: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
697: .seealso: PCFactorSetReuseOrdering()
698: @*/
699: PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag)
700: {
706: PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
707: return(0);
708: }
710: PetscErrorCode PCFactorInitialize(PC pc)
711: {
713: PC_Factor *fact = (PC_Factor*)pc->data;
716: MatFactorInfoInitialize(&fact->info);
717: fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
718: fact->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
719: fact->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON;
720: fact->info.pivotinblocks = 1.0;
721: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
723: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
724: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);
725: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
726: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);
727: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
728: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);
729: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",PCFactorGetMatSolverType_Factor);
730: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",PCFactorSetMatSolverType_Factor);
731: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",PCFactorSetUpMatSolverType_Factor);
732: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);
733: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);
734: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);
735: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);
736: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);
737: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);
738: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);
739: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor);
740: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor);
741: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor);
742: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor);
743: return(0);
744: }