Actual source code: factor.c
petsc-3.13.6 2020-09-29
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: Level: intermediate
52: .seealso: PCFactorSetMatSolverType(), PCFactorGetMatrix()
53: @*/
54: PetscErrorCode PCFactorSetUpMatSolverType(PC pc)
55: {
60: PetscTryMethod(pc,"PCFactorSetUpMatSolverType_C",(PC),(pc));
61: return(0);
62: }
64: /*@
65: PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
67: Logically Collective on PC
69: Input Parameters:
70: + pc - the preconditioner context
71: - zero - all pivots smaller than this will be considered zero
73: Options Database Key:
74: . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot
76: Level: intermediate
78: .seealso: PCFactorSetShiftType(), PCFactorSetShiftAmount()
79: @*/
80: PetscErrorCode PCFactorSetZeroPivot(PC pc,PetscReal zero)
81: {
87: PetscTryMethod(pc,"PCFactorSetZeroPivot_C",(PC,PetscReal),(pc,zero));
88: return(0);
89: }
91: /*@
92: PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during
93: numerical factorization, thus the matrix has nonzero pivots
95: Logically Collective on PC
97: Input Parameters:
98: + pc - the preconditioner context
99: - shifttype - type of shift; one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, MAT_SHIFT_INBLOCKS
101: Options Database Key:
102: . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types
104: Level: intermediate
106: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftAmount()
107: @*/
108: PetscErrorCode PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype)
109: {
115: PetscTryMethod(pc,"PCFactorSetShiftType_C",(PC,MatFactorShiftType),(pc,shifttype));
116: return(0);
117: }
119: /*@
120: PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during
121: numerical factorization, thus the matrix has nonzero pivots
123: Logically Collective on PC
125: Input Parameters:
126: + pc - the preconditioner context
127: - shiftamount - amount of shift
129: Options Database Key:
130: . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default
132: Level: intermediate
134: .seealso: PCFactorSetZeroPivot(), PCFactorSetShiftType()
135: @*/
136: PetscErrorCode PCFactorSetShiftAmount(PC pc,PetscReal shiftamount)
137: {
143: PetscTryMethod(pc,"PCFactorSetShiftAmount_C",(PC,PetscReal),(pc,shiftamount));
144: return(0);
145: }
147: /*@
148: PCFactorSetDropTolerance - The preconditioner will use an ILU
149: based on a drop tolerance. (Under development)
151: Logically Collective on PC
153: Input Parameters:
154: + pc - the preconditioner context
155: . dt - the drop tolerance, try from 1.e-10 to .1
156: . dtcol - tolerance for column pivot, good values [0.1 to 0.01]
157: - maxrowcount - the max number of nonzeros allowed in a row, best value
158: depends on the number of nonzeros in row of original matrix
160: Options Database Key:
161: . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
163: Level: intermediate
165: There are NO default values for the 3 parameters, you must set them with reasonable values for your
166: matrix. We don't know how to compute reasonable values.
168: @*/
169: PetscErrorCode PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
170: {
177: PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount));
178: return(0);
179: }
181: /*@
182: PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot
184: Not Collective
186: Input Parameters:
187: . pc - the preconditioner context
189: Output Parameter:
190: . pivot - the tolerance
192: Level: intermediate
194: .seealso: PCFactorSetZeroPivot()
195: @*/
196: PetscErrorCode PCFactorGetZeroPivot(PC pc,PetscReal *pivot)
197: {
202: PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot));
203: return(0);
204: }
206: /*@
207: PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot
209: Not Collective
211: Input Parameters:
212: . pc - the preconditioner context
214: Output Parameter:
215: . shift - how much to shift the diagonal entry
217: Level: intermediate
219: .seealso: PCFactorSetShiftAmount(), PCFactorSetShiftType(), PCFactorGetShiftType()
220: @*/
221: PetscErrorCode PCFactorGetShiftAmount(PC pc,PetscReal *shift)
222: {
227: PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift));
228: return(0);
229: }
231: /*@
232: PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected
234: Not Collective
236: Input Parameters:
237: . pc - the preconditioner context
239: Output Parameter:
240: . type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS
242: Level: intermediate
244: .seealso: PCFactorSetShiftType(), MatFactorShiftType, PCFactorSetShiftAmount(), PCFactorGetShiftAmount()
245: @*/
246: PetscErrorCode PCFactorGetShiftType(PC pc,MatFactorShiftType *type)
247: {
252: PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type));
253: return(0);
254: }
256: /*@
257: PCFactorGetLevels - Gets the number of levels of fill to use.
259: Logically Collective on PC
261: Input Parameters:
262: . pc - the preconditioner context
264: Output Parameter:
265: . levels - number of levels of fill
267: Level: intermediate
269: @*/
270: PetscErrorCode PCFactorGetLevels(PC pc,PetscInt *levels)
271: {
276: PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels));
277: return(0);
278: }
280: /*@
281: PCFactorSetLevels - Sets the number of levels of fill to use.
283: Logically Collective on PC
285: Input Parameters:
286: + pc - the preconditioner context
287: - levels - number of levels of fill
289: Options Database Key:
290: . -pc_factor_levels <levels> - Sets fill level
292: Level: intermediate
294: @*/
295: PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels)
296: {
301: if (levels < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
303: PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels));
304: return(0);
305: }
307: /*@
308: PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be
309: treated as level 0 fill even if there is no non-zero location.
311: Logically Collective on PC
313: Input Parameters:
314: + pc - the preconditioner context
315: - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
317: Options Database Key:
318: . -pc_factor_diagonal_fill
320: Notes:
321: Does not apply with 0 fill.
323: Level: intermediate
325: .seealso: PCFactorGetAllowDiagonalFill()
326: @*/
327: PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg)
328: {
333: PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg));
334: return(0);
335: }
337: /*@
338: PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are
339: treated as level 0 fill even if there is no non-zero location.
341: Logically Collective on PC
343: Input Parameter:
344: . pc - the preconditioner context
346: Output Parameter:
347: . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off
349: Options Database Key:
350: . -pc_factor_diagonal_fill
352: Notes:
353: Does not apply with 0 fill.
355: Level: intermediate
357: .seealso: PCFactorSetAllowDiagonalFill()
358: @*/
359: PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg)
360: {
365: PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg));
366: return(0);
367: }
369: /*@
370: PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
372: Logically Collective on PC
374: Input Parameters:
375: + pc - the preconditioner context
376: - tol - diagonal entries smaller than this in absolute value are considered zero
378: Options Database Key:
379: . -pc_factor_nonzeros_along_diagonal <tol>
381: Level: intermediate
383: .seealso: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetZeroPivot(), MatReorderForNonzeroDiagonal()
384: @*/
385: PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
386: {
392: PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol));
393: return(0);
394: }
396: /*@C
397: PCFactorSetMatSolverType - sets the software that is used to perform the factorization
399: Logically Collective on PC
401: Input Parameters:
402: + pc - the preconditioner context
403: - stype - for example, superlu, superlu_dist
405: Options Database Key:
406: . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse
408: Level: intermediate
410: Note:
411: By default this will use the PETSc factorization if it exists
413: .seealso: MatGetFactor(), MatSolverType, PCFactorGetMatSolverType()
414: @*/
415: PetscErrorCode PCFactorSetMatSolverType(PC pc,MatSolverType stype)
416: {
421: PetscTryMethod(pc,"PCFactorSetMatSolverType_C",(PC,MatSolverType),(pc,stype));
422: return(0);
423: }
425: /*@C
426: PCFactorGetMatSolverType - gets the software that is used to perform the factorization
428: Not Collective
430: Input Parameter:
431: . pc - the preconditioner context
433: Output Parameter:
434: . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package)
436: Level: intermediate
438: .seealso: MatGetFactor(), MatSolverType, PCFactorGetMatSolverType()
439: @*/
440: PetscErrorCode PCFactorGetMatSolverType(PC pc,MatSolverType *stype)
441: {
442: PetscErrorCode ierr,(*f)(PC,MatSolverType*);
446: PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",&f);
447: if (f) {
448: (*f)(pc,stype);
449: } else {
450: *stype = NULL;
451: }
452: return(0);
453: }
455: /*@
456: PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix,
457: fill = number nonzeros in factor/number nonzeros in original matrix.
459: Not Collective, each process can expect a different amount of fill
461: Input Parameters:
462: + pc - the preconditioner context
463: - fill - amount of expected fill
465: Options Database Key:
466: . -pc_factor_fill <fill> - Sets fill amount
468: Level: intermediate
470: Note:
471: For sparse matrix factorizations it is difficult to predict how much
472: fill to expect. By running with the option -info PETSc will print the
473: actual amount of fill used; allowing you to set the value accurately for
474: future runs. Default PETSc uses a value of 5.0
476: This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels.
478: @*/
479: PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill)
480: {
485: if (fill < 1.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
486: PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill));
487: return(0);
488: }
490: /*@
491: PCFactorSetUseInPlace - Tells the system to do an in-place factorization.
492: For dense matrices, this enables the solution of much larger problems.
493: For sparse matrices the factorization cannot be done truly in-place
494: so this does not save memory during the factorization, but after the matrix
495: is factored, the original unfactored matrix is freed, thus recovering that
496: space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place.
498: Logically Collective on PC
500: Input Parameters:
501: + pc - the preconditioner context
502: - flg - PETSC_TRUE to enable, PETSC_FALSE to disable
504: Options Database Key:
505: . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization
507: Notes:
508: PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when
509: a different matrix is provided for the multiply and the preconditioner in
510: a call to KSPSetOperators().
511: This is because the Krylov space methods require an Section 1.5 Writing Application Codes with PETSc of the
512: matrix multiplication, which is not possible here because the matrix has
513: been factored in-place, replacing the original matrix.
515: Level: intermediate
517: .seealso: PCFactorGetUseInPlace()
518: @*/
519: PetscErrorCode PCFactorSetUseInPlace(PC pc,PetscBool flg)
520: {
525: PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg));
526: return(0);
527: }
529: /*@
530: PCFactorGetUseInPlace - Determines if an in-place factorization is being used.
532: Logically Collective on PC
534: Input Parameter:
535: . pc - the preconditioner context
537: Output Parameter:
538: . flg - PETSC_TRUE to enable, PETSC_FALSE to disable
540: Level: intermediate
542: .seealso: PCFactorSetUseInPlace()
543: @*/
544: PetscErrorCode PCFactorGetUseInPlace(PC pc,PetscBool *flg)
545: {
550: PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg));
551: return(0);
552: }
554: /*@C
555: PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to
556: be used in the LU factorization.
558: Logically Collective on PC
560: Input Parameters:
561: + pc - the preconditioner context
562: - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM
564: Options Database Key:
565: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
567: Level: intermediate
569: Notes:
570: nested dissection is used by default
572: For Cholesky and ICC and the SBAIJ format the only reordering available is natural since only the upper half of the matrix is stored
573: and reordering this matrix is very expensive.
575: You can use SeqAIJ matrix with Cholesky and ICC and use any ordering
577: .seealso: MatOrderingType
579: @*/
580: PetscErrorCode PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering)
581: {
586: PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering));
587: return(0);
588: }
590: /*@
591: PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization.
592: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
593: it is never done. For the MATLAB and SuperLU factorization this is used.
595: Logically Collective on PC
597: Input Parameters:
598: + pc - the preconditioner context
599: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
601: Options Database Key:
602: . -pc_factor_pivoting <dtcol>
604: Level: intermediate
606: .seealso: PCILUSetMatOrdering(), PCFactorSetPivotInBlocks()
607: @*/
608: PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol)
609: {
615: PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol));
616: return(0);
617: }
619: /*@
620: PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block
621: with BAIJ or SBAIJ matrices
623: Logically Collective on PC
625: Input Parameters:
626: + pc - the preconditioner context
627: - pivot - PETSC_TRUE or PETSC_FALSE
629: Options Database Key:
630: . -pc_factor_pivot_in_blocks <true,false>
632: Level: intermediate
634: .seealso: PCILUSetMatOrdering(), PCFactorSetColumnPivot()
635: @*/
636: PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot)
637: {
643: PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot));
644: return(0);
645: }
647: /*@
648: PCFactorSetReuseFill - When matrices with different nonzero structure are factored,
649: this causes later ones to use the fill ratio computed in the initial factorization.
651: Logically Collective on PC
653: Input Parameters:
654: + pc - the preconditioner context
655: - flag - PETSC_TRUE to reuse else PETSC_FALSE
657: Options Database Key:
658: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
660: Level: intermediate
662: .seealso: PCFactorSetReuseOrdering()
663: @*/
664: PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag)
665: {
671: PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag));
672: return(0);
673: }
675: PetscErrorCode PCFactorInitialize(PC pc)
676: {
678: PC_Factor *fact = (PC_Factor*)pc->data;
681: MatFactorInfoInitialize(&fact->info);
682: fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE;
683: fact->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
684: fact->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON;
685: fact->info.pivotinblocks = 1.0;
686: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
688: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);
689: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);
690: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);
691: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);
692: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);
693: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);
694: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",PCFactorGetMatSolverType_Factor);
695: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",PCFactorSetMatSolverType_Factor);
696: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",PCFactorSetUpMatSolverType_Factor);
697: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);
698: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);
699: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);
700: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);
701: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);
702: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor);
703: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);
704: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor);
705: PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor);
706: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor);
707: PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor);
708: return(0);
709: }