Actual source code: lu.c
1: /*$Id: lu.c,v 1.149 2001/08/07 21:30:26 bsmith Exp $*/
2: /*
3: Defines a direct factorization preconditioner for any Mat implementation
4: Note: this need not be consided a preconditioner since it supplies
5: a direct solver.
6: */
7: #include src/ksp/pc/pcimpl.h
9: typedef struct {
10: Mat fact; /* factored matrix */
11: PetscReal actualfill; /* actual fill in factor */
12: int inplace; /* flag indicating in-place factorization */
13: IS row,col; /* index sets used for reordering */
14: MatOrderingType ordering; /* matrix ordering */
15: PetscTruth reuseordering; /* reuses previous reordering computed */
16: PetscTruth reusefill; /* reuse fill from previous LU */
17: MatFactorInfo info;
18: } PC_LU;
21: EXTERN_C_BEGIN
24: int PCLUSetZeroPivot_LU(PC pc,PetscReal z)
25: {
26: PC_LU *lu;
29: lu = (PC_LU*)pc->data;
30: lu->info.zeropivot = z;
31: return(0);
32: }
33: EXTERN_C_END
35: EXTERN_C_BEGIN
38: int PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag)
39: {
40: PC_LU *lu;
43: lu = (PC_LU*)pc->data;
44: lu->reuseordering = flag;
45: return(0);
46: }
47: EXTERN_C_END
49: EXTERN_C_BEGIN
52: int PCLUSetReuseFill_LU(PC pc,PetscTruth flag)
53: {
54: PC_LU *lu;
57: lu = (PC_LU*)pc->data;
58: lu->reusefill = flag;
59: return(0);
60: }
61: EXTERN_C_END
63: EXTERN_C_BEGIN
66: int PCLUSetShift_LU(PC pc,PetscTruth shift)
67: {
68: PC_LU *dir;
69:
71: dir = (PC_LU*)pc->data;
72: dir->info.shift = shift;
73: if (shift) dir->info.shift_fraction = 0.0;
74: return(0);
75: }
76: EXTERN_C_END
80: static int PCSetFromOptions_LU(PC pc)
81: {
82: PC_LU *lu = (PC_LU*)pc->data;
83: int ierr;
84: PetscTruth flg,set;
85: char tname[256];
86: PetscFList ordlist;
87: PetscReal tol;
90: MatOrderingRegisterAll(PETSC_NULL);
91: PetscOptionsHead("LU options");
92: PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);
93: if (flg) {
94: PCLUSetUseInPlace(pc);
95: }
96: PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);
98: PetscOptionsName("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",&flg);
99: if (flg) {
100: PCLUSetDamping(pc,(PetscReal) PETSC_DECIDE);
101: }
102: PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);
103: PetscOptionsName("-pc_lu_shift","Manteuffel shift applied to diagonal","PCLUSetShift",&flg);
104: if (flg) {
105: PCLUSetShift(pc,PETSC_TRUE);
106: }
107: PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);
109: PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);
110: if (flg) {
111: PCLUSetReuseFill(pc,PETSC_TRUE);
112: }
113: PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);
114: if (flg) {
115: PCLUSetReuseOrdering(pc,PETSC_TRUE);
116: }
118: MatGetOrderingList(&ordlist);
119: PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
120: if (flg) {
121: PCLUSetMatOrdering(pc,tname);
122: }
123: PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);
125: PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);
127: flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
128: PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);
129: if (set) {
130: PCLUSetPivotInBlocks(pc,flg);
131: }
132: PetscOptionsTail();
133: return(0);
134: }
138: static int PCView_LU(PC pc,PetscViewer viewer)
139: {
140: PC_LU *lu = (PC_LU*)pc->data;
141: int ierr;
142: PetscTruth isascii,isstring;
145: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
146: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
147: if (isascii) {
148: MatInfo info;
150: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");}
151: else {PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");}
152: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);
153: PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %g\n",lu->info.zeropivot);
154: if (lu->info.shift) {PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");}
155: if (lu->fact) {
156: MatGetInfo(lu->fact,MAT_LOCAL,&info);
157: PetscViewerASCIIPrintf(viewer," LU nonzeros %g\n",info.nz_used);
158: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
159: MatView(lu->fact,viewer);
160: PetscViewerPopFormat(viewer);
161: }
162: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
163: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
164: } else if (isstring) {
165: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
166: } else {
167: SETERRQ1(1,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
168: }
169: return(0);
170: }
174: static int PCGetFactoredMatrix_LU(PC pc,Mat *mat)
175: {
176: PC_LU *dir = (PC_LU*)pc->data;
179: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
180: *mat = dir->fact;
181: return(0);
182: }
186: static int PCSetUp_LU(PC pc)
187: {
188: int ierr;
189: PetscTruth flg;
190: PC_LU *dir = (PC_LU*)pc->data;
193: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
195: if (dir->inplace) {
196: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
197: if (dir->col) {ISDestroy(dir->col);}
198: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
199: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
200: MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);
201: dir->fact = pc->pmat;
202: } else {
203: MatInfo info;
204: if (!pc->setupcalled) {
205: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
206: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
207: if (flg) {
208: PetscReal tol = 1.e-10;
209: PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
210: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
211: }
212: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
213: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
214: MatGetInfo(dir->fact,MAT_LOCAL,&info);
215: dir->actualfill = info.fill_ratio_needed;
216: PetscLogObjectParent(pc,dir->fact);
217: } else if (pc->flag != SAME_NONZERO_PATTERN) {
218: if (!dir->reuseordering) {
219: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
220: if (dir->col) {ISDestroy(dir->col);}
221: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
222: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
223: if (flg) {
224: PetscReal tol = 1.e-10;
225: PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
226: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
227: }
228: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
229: }
230: MatDestroy(dir->fact);
231: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
232: MatGetInfo(dir->fact,MAT_LOCAL,&info);
233: dir->actualfill = info.fill_ratio_needed;
234: PetscLogObjectParent(pc,dir->fact);
235: }
236: MatLUFactorNumeric(pc->pmat,&dir->fact);
237: }
238: return(0);
239: }
243: static int PCDestroy_LU(PC pc)
244: {
245: PC_LU *dir = (PC_LU*)pc->data;
246: int ierr;
249: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
250: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
251: if (dir->col) {ISDestroy(dir->col);}
252: PetscStrfree(dir->ordering);
253: PetscFree(dir);
254: return(0);
255: }
259: static int PCApply_LU(PC pc,Vec x,Vec y)
260: {
261: PC_LU *dir = (PC_LU*)pc->data;
262: int ierr;
265: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
266: else {MatSolve(dir->fact,x,y);}
267: return(0);
268: }
272: static int PCApplyTranspose_LU(PC pc,Vec x,Vec y)
273: {
274: PC_LU *dir = (PC_LU*)pc->data;
275: int ierr;
278: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
279: else {MatSolveTranspose(dir->fact,x,y);}
280: return(0);
281: }
283: /* -----------------------------------------------------------------------------------*/
285: EXTERN_C_BEGIN
288: int PCLUSetFill_LU(PC pc,PetscReal fill)
289: {
290: PC_LU *dir;
293: dir = (PC_LU*)pc->data;
294: dir->info.fill = fill;
295: return(0);
296: }
297: EXTERN_C_END
299: EXTERN_C_BEGIN
302: int PCLUSetDamping_LU(PC pc,PetscReal damping)
303: {
304: PC_LU *dir;
307: dir = (PC_LU*)pc->data;
308: if (damping == (PetscReal) PETSC_DECIDE) {
309: dir->info.damping = 1.e-12;
310: } else {
311: dir->info.damping = damping;
312: }
313: return(0);
314: }
315: EXTERN_C_END
317: EXTERN_C_BEGIN
320: int PCLUSetUseInPlace_LU(PC pc)
321: {
322: PC_LU *dir;
325: dir = (PC_LU*)pc->data;
326: dir->inplace = 1;
327: return(0);
328: }
329: EXTERN_C_END
331: EXTERN_C_BEGIN
334: int PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
335: {
336: PC_LU *dir = (PC_LU*)pc->data;
337: int ierr;
340: PetscStrfree(dir->ordering);
341: PetscStrallocpy(ordering,&dir->ordering);
342: return(0);
343: }
344: EXTERN_C_END
346: EXTERN_C_BEGIN
349: int PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
350: {
351: PC_LU *dir = (PC_LU*)pc->data;
354: if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(1,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
355: dir->info.dtcol = dtcol;
356: return(0);
357: }
358: EXTERN_C_END
360: EXTERN_C_BEGIN
363: int PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
364: {
365: PC_LU *dir = (PC_LU*)pc->data;
368: dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
369: return(0);
370: }
371: EXTERN_C_END
373: /* -----------------------------------------------------------------------------------*/
377: /*@
378: PCLUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
380: Collective on PC
381:
382: Input Parameters:
383: + pc - the preconditioner context
384: - zero - all pivots smaller than this will be considered zero
386: Options Database Key:
387: . -pc_lu_zeropivot <zero> - Sets the zero pivot size
389: Level: intermediate
391: .keywords: PC, set, factorization, direct, fill
393: .seealso: PCLUSetFill(), PCLUSetDamp(), PCILUSetZeroPivot()
394: @*/
395: int PCLUSetZeroPivot(PC pc,PetscReal zero)
396: {
397: int ierr,(*f)(PC,PetscReal);
401: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetZeroPivot_C",(void (**)(void))&f);
402: if (f) {
403: (*f)(pc,zero);
404: }
405: return(0);
406: }
410: /*@
411: PCLUSetShift - specify whether to use Manteuffel shifting of LU.
412: If an LU factorisation breaks down because of nonpositive pivots,
413: adding sufficient identity to the diagonal will remedy this.
414: Setting this causes a bisection method to find the minimum shift that
415: will lead to a well-defined LU.
417: Input parameters:
418: + pc - the preconditioner context
419: - shifting - PETSC_TRUE to set shift else PETSC_FALSE
421: Options Database Key:
422: . -pc_lu_shift - Activate PCLUSetShift()
424: Level: intermediate
426: .keywords: PC, indefinite, factorization
428: .seealso: PCLUSetDamping(), PCILUSetShift()
429: @*/
430: int PCLUSetShift(PC pc,PetscTruth shifting)
431: {
432: int ierr,(*f)(PC,PetscTruth);
436: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetShift_C",(void (**)(void))&f);
437: if (f) {
438: (*f)(pc,shifting);
439: }
440: return(0);
441: }
446: /*@
447: PCLUSetReuseOrdering - When similar matrices are factored, this
448: causes the ordering computed in the first factor to be used for all
449: following factors; applies to both fill and drop tolerance LUs.
451: Collective on PC
453: Input Parameters:
454: + pc - the preconditioner context
455: - flag - PETSC_TRUE to reuse else PETSC_FALSE
457: Options Database Key:
458: . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
460: Level: intermediate
462: .keywords: PC, levels, reordering, factorization, incomplete, LU
464: .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
465: @*/
466: int PCLUSetReuseOrdering(PC pc,PetscTruth flag)
467: {
468: int ierr,(*f)(PC,PetscTruth);
472: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);
473: if (f) {
474: (*f)(pc,flag);
475: }
476: return(0);
477: }
481: /*@
482: PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
483: this causes later ones to use the fill computed in the initial factorization.
485: Collective on PC
487: Input Parameters:
488: + pc - the preconditioner context
489: - flag - PETSC_TRUE to reuse else PETSC_FALSE
491: Options Database Key:
492: . -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
494: Level: intermediate
496: .keywords: PC, levels, reordering, factorization, incomplete, LU
498: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
499: @*/
500: int PCLUSetReuseFill(PC pc,PetscTruth flag)
501: {
502: int ierr,(*f)(PC,PetscTruth);
506: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);
507: if (f) {
508: (*f)(pc,flag);
509: }
510: return(0);
511: }
515: /*@
516: PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
517: fill = number nonzeros in factor/number nonzeros in original matrix.
519: Collective on PC
520:
521: Input Parameters:
522: + pc - the preconditioner context
523: - fill - amount of expected fill
525: Options Database Key:
526: . -pc_lu_fill <fill> - Sets fill amount
528: Level: intermediate
530: Note:
531: For sparse matrix factorizations it is difficult to predict how much
532: fill to expect. By running with the option -log_info PETSc will print the
533: actual amount of fill used; allowing you to set the value accurately for
534: future runs. Default PETSc uses a value of 5.0
536: .keywords: PC, set, factorization, direct, fill
538: .seealso: PCILUSetFill()
539: @*/
540: int PCLUSetFill(PC pc,PetscReal fill)
541: {
542: int ierr,(*f)(PC,PetscReal);
546: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
547: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);
548: if (f) {
549: (*f)(pc,fill);
550: }
551: return(0);
552: }
556: /*@
557: PCLUSetDamping - adds this quantity to the diagonal of the matrix during the
558: LU numerical factorization
560: Collective on PC
561:
562: Input Parameters:
563: + pc - the preconditioner context
564: - damping - amount of damping (use PETSC_DECIDE for default of 1.e-12)
566: Options Database Key:
567: . -pc_lu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default
569: Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
570: pivot, then the damping is doubled until this is alleviated.
572: Level: intermediate
574: .keywords: PC, set, factorization, direct, fill
575: .seealso: PCILUSetFill(), PCILUSetDamp()
576: @*/
577: int PCLUSetDamping(PC pc,PetscReal damping)
578: {
579: int ierr,(*f)(PC,PetscReal);
583: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);
584: if (f) {
585: (*f)(pc,damping);
586: }
587: return(0);
588: }
592: /*@
593: PCLUSetUseInPlace - Tells the system to do an in-place factorization.
594: For dense matrices, this enables the solution of much larger problems.
595: For sparse matrices the factorization cannot be done truly in-place
596: so this does not save memory during the factorization, but after the matrix
597: is factored, the original unfactored matrix is freed, thus recovering that
598: space.
600: Collective on PC
602: Input Parameters:
603: . pc - the preconditioner context
605: Options Database Key:
606: . -pc_lu_in_place - Activates in-place factorization
608: Notes:
609: PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when
610: a different matrix is provided for the multiply and the preconditioner in
611: a call to KSPSetOperators().
612: This is because the Krylov space methods require an application of the
613: matrix multiplication, which is not possible here because the matrix has
614: been factored in-place, replacing the original matrix.
616: Level: intermediate
618: .keywords: PC, set, factorization, direct, inplace, in-place, LU
620: .seealso: PCILUSetUseInPlace()
621: @*/
622: int PCLUSetUseInPlace(PC pc)
623: {
624: int ierr,(*f)(PC);
628: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);
629: if (f) {
630: (*f)(pc);
631: }
632: return(0);
633: }
637: /*@C
638: PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to
639: be used in the LU factorization.
641: Collective on PC
643: Input Parameters:
644: + pc - the preconditioner context
645: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
647: Options Database Key:
648: . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
650: Level: intermediate
652: Notes: nested dissection is used by default
654: .seealso: PCILUSetMatOrdering()
655: @*/
656: int PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
657: {
658: int ierr,(*f)(PC,MatOrderingType);
661: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);
662: if (f) {
663: (*f)(pc,ordering);
664: }
665: return(0);
666: }
670: /*@
671: PCLUSetPivoting - Determines when pivoting is done during LU.
672: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
673: it is never done. For the Matlab and SuperLU factorization this is used.
675: Collective on PC
677: Input Parameters:
678: + pc - the preconditioner context
679: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
681: Options Database Key:
682: . -pc_lu_pivoting <dtcol>
684: Level: intermediate
686: .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks()
687: @*/
688: int PCLUSetPivoting(PC pc,PetscReal dtcol)
689: {
690: int ierr,(*f)(PC,PetscReal);
693: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);
694: if (f) {
695: (*f)(pc,dtcol);
696: }
697: return(0);
698: }
702: /*@
703: PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block
704: with BAIJ or SBAIJ matrices
706: Collective on PC
708: Input Parameters:
709: + pc - the preconditioner context
710: - pivot - PETSC_TRUE or PETSC_FALSE
712: Options Database Key:
713: . -pc_lu_pivot_in_blocks <true,false>
715: Level: intermediate
717: .seealso: PCILUSetMatOrdering(), PCLUSetPivoting()
718: @*/
719: int PCLUSetPivotInBlocks(PC pc,PetscTruth pivot)
720: {
721: int ierr,(*f)(PC,PetscTruth);
724: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);
725: if (f) {
726: (*f)(pc,pivot);
727: }
728: return(0);
729: }
731: /* ------------------------------------------------------------------------ */
733: /*MC
734: PCLU - Uses a direct solver, based on LU factorization, as a preconditioner
736: Options Database Keys:
737: + -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
738: . -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
739: . -pc_lu_fill <fill> - Sets fill amount
740: . -pc_lu_damping <damping> - Sets damping amount
741: . -pc_lu_in_place - Activates in-place factorization
742: . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
743: . -pc_lu_pivoting <dtcol>
744: - -pc_lu_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase
745: stability of factorization.
747: Notes: Not all options work for all matrix formats
748: Run with -help to see additional options for particular matrix formats or factorization
749: algorithms
751: Level: beginner
753: Concepts: LU factorization, direct solver
755: Notes: Usually this will compute an "exact" solution in one iteration and does
756: not need a Krylov method (i.e. you can use -ksp_type preonly, or
757: KSPSetType(ksp,KSPPREONLY) for the Krylov method
759: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
760: PCILU, PCCHOLESKY, PCICC, PCLUSetReuseOrdering(), PCLUSetReuseFill(), PCGetFactoredMatrix(),
761: PCLUSetFill(), PCLUSetDamping(), PCLUSetUseInPlace(), PCLUSetMatOrdering(), PCLUSetPivoting(),
762: PCLUSetPivotingInBlocks()
763: M*/
765: EXTERN_C_BEGIN
768: int PCCreate_LU(PC pc)
769: {
770: int ierr,size;
771: PC_LU *dir;
774: PetscNew(PC_LU,&dir);
775: PetscLogObjectMemory(pc,sizeof(PC_LU));
777: dir->fact = 0;
778: dir->inplace = 0;
779: dir->info.fill = 5.0;
780: dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
781: dir->info.damping = 0.0;
782: dir->info.zeropivot = 1.e-12;
783: dir->info.pivotinblocks = 1.0;
784: dir->info.shift = PETSC_FALSE;
785: dir->info.shift_fraction = 0.0;
786: dir->col = 0;
787: dir->row = 0;
788: MPI_Comm_size(pc->comm,&size);
789: if (size == 1) {
790: PetscStrallocpy(MATORDERING_ND,&dir->ordering);
791: } else {
792: PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
793: }
794: dir->reusefill = PETSC_FALSE;
795: dir->reuseordering = PETSC_FALSE;
796: pc->data = (void*)dir;
798: pc->ops->destroy = PCDestroy_LU;
799: pc->ops->apply = PCApply_LU;
800: pc->ops->applytranspose = PCApplyTranspose_LU;
801: pc->ops->setup = PCSetUp_LU;
802: pc->ops->setfromoptions = PCSetFromOptions_LU;
803: pc->ops->view = PCView_LU;
804: pc->ops->applyrichardson = 0;
805: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
807: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
808: PCLUSetFill_LU);
809: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
810: PCLUSetDamping_LU);
811: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetShift_C","PCLUSetShift_LU",
812: PCLUSetShift_LU);
813: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
814: PCLUSetUseInPlace_LU);
815: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
816: PCLUSetMatOrdering_LU);
817: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
818: PCLUSetReuseOrdering_LU);
819: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
820: PCLUSetReuseFill_LU);
821: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
822: PCLUSetPivoting_LU);
823: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU",
824: PCLUSetPivotInBlocks_LU);
825: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetZeroPivot_C","PCLUSetZeroPivot_LU",
826: PCLUSetZeroPivot_LU);
827: return(0);
828: }
829: EXTERN_C_END