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/sles/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: MatLUInfo info;
18: } PC_LU;
21: EXTERN_C_BEGIN
22: int PCLUSetReuseOrdering_LU(PC pc,PetscTruth flag)
23: {
24: PC_LU *lu;
27: lu = (PC_LU*)pc->data;
28: lu->reuseordering = flag;
29: return(0);
30: }
31: EXTERN_C_END
33: EXTERN_C_BEGIN
34: int PCLUSetReuseFill_LU(PC pc,PetscTruth flag)
35: {
36: PC_LU *lu;
39: lu = (PC_LU*)pc->data;
40: lu->reusefill = flag;
41: return(0);
42: }
43: EXTERN_C_END
45: static int PCSetFromOptions_LU(PC pc)
46: {
47: PC_LU *lu = (PC_LU*)pc->data;
48: int ierr;
49: PetscTruth flg,set;
50: char tname[256];
51: PetscFList ordlist;
52: PetscReal tol;
55: MatOrderingRegisterAll(PETSC_NULL);
56: PetscOptionsHead("LU options");
57: PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);
58: if (flg) {
59: PCLUSetUseInPlace(pc);
60: }
61: PetscOptionsReal("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);
63: PetscOptionsHasName(pc->prefix,"-pc_lu_damping",&flg);
64: if (flg) {
65: PCLUSetDamping(pc,0.0);
66: }
67: PetscOptionsReal("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);
68: PetscOptionsReal("-pc_lu_zeropivot","Pivot is considered zero if less than","PCLUSetSetZeroPivot",lu->info.zeropivot,&lu->info.zeropivot,0);
70: PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);
71: if (flg) {
72: PCLUSetReuseFill(pc,PETSC_TRUE);
73: }
74: PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);
75: if (flg) {
76: PCLUSetReuseOrdering(pc,PETSC_TRUE);
77: }
79: MatGetOrderingList(&ordlist);
80: PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
81: if (flg) {
82: PCLUSetMatOrdering(pc,tname);
83: }
84: PetscOptionsReal("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);
86: PetscOptionsReal("-pc_lu_pivoting","Pivoting tolerance (used only for some factorization)","PCLUSetPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);
88: flg = lu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
89: PetscOptionsLogical("-pc_lu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCLUSetPivotInBlocks",flg,&flg,&set);
90: if (set) {
91: PCLUSetPivotInBlocks(pc,flg);
92: }
93: PetscOptionsTail();
94: return(0);
95: }
97: static int PCView_LU(PC pc,PetscViewer viewer)
98: {
99: PC_LU *lu = (PC_LU*)pc->data;
100: int ierr;
101: PetscTruth isascii,isstring;
104: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
105: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
106: if (isascii) {
107: MatInfo info;
109: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," LU: in-place factorizationn");}
110: else {PetscViewerASCIIPrintf(viewer," LU: out-of-place factorizationn");}
111: PetscViewerASCIIPrintf(viewer," matrix ordering: %sn",lu->ordering);
112: PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %gn",lu->info.zeropivot);
113: if (lu->fact) {
114: MatGetInfo(lu->fact,MAT_LOCAL,&info);
115: PetscViewerASCIIPrintf(viewer," LU nonzeros %gn",info.nz_used);
116: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
117: MatView(lu->fact,viewer);
118: PetscViewerPopFormat(viewer);
119: }
120: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorizationn");}
121: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorizationn");}
122: } else if (isstring) {
123: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
124: } else {
125: SETERRQ1(1,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
126: }
127: return(0);
128: }
130: static int PCGetFactoredMatrix_LU(PC pc,Mat *mat)
131: {
132: PC_LU *dir = (PC_LU*)pc->data;
135: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
136: *mat = dir->fact;
137: return(0);
138: }
140: static int PCSetUp_LU(PC pc)
141: {
142: int ierr;
143: PetscTruth flg;
144: PC_LU *dir = (PC_LU*)pc->data;
147: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
149: if (dir->inplace) {
150: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
151: if (dir->col) {ISDestroy(dir->col);}
152: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
153: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
154: MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);
155: dir->fact = pc->pmat;
156: } else {
157: MatInfo info;
158: if (!pc->setupcalled) {
159: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
160: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
161: if (flg) {
162: PetscReal tol = 1.e-10;
163: PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
164: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
165: }
166: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
167: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
168: MatGetInfo(dir->fact,MAT_LOCAL,&info);
169: dir->actualfill = info.fill_ratio_needed;
170: PetscLogObjectParent(pc,dir->fact);
171: } else if (pc->flag != SAME_NONZERO_PATTERN) {
172: if (!dir->reuseordering) {
173: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
174: if (dir->col) {ISDestroy(dir->col);}
175: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
176: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
177: if (flg) {
178: PetscReal tol = 1.e-10;
179: PetscOptionsGetReal(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
180: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
181: }
182: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
183: }
184: MatDestroy(dir->fact);
185: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
186: MatGetInfo(dir->fact,MAT_LOCAL,&info);
187: dir->actualfill = info.fill_ratio_needed;
188: PetscLogObjectParent(pc,dir->fact);
189: }
190: MatLUFactorNumeric(pc->pmat,&dir->fact);
191: }
192: return(0);
193: }
195: static int PCDestroy_LU(PC pc)
196: {
197: PC_LU *dir = (PC_LU*)pc->data;
198: int ierr;
201: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
202: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
203: if (dir->col) {ISDestroy(dir->col);}
204: PetscStrfree(dir->ordering);
205: PetscFree(dir);
206: return(0);
207: }
209: static int PCApply_LU(PC pc,Vec x,Vec y)
210: {
211: PC_LU *dir = (PC_LU*)pc->data;
212: int ierr;
215: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
216: else {MatSolve(dir->fact,x,y);}
217: return(0);
218: }
220: static int PCApplyTranspose_LU(PC pc,Vec x,Vec y)
221: {
222: PC_LU *dir = (PC_LU*)pc->data;
223: int ierr;
226: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
227: else {MatSolveTranspose(dir->fact,x,y);}
228: return(0);
229: }
231: /* -----------------------------------------------------------------------------------*/
233: EXTERN_C_BEGIN
234: int PCLUSetFill_LU(PC pc,PetscReal fill)
235: {
236: PC_LU *dir;
239: dir = (PC_LU*)pc->data;
240: dir->info.fill = fill;
241: return(0);
242: }
243: EXTERN_C_END
245: EXTERN_C_BEGIN
246: int PCLUSetDamping_LU(PC pc,PetscReal damping)
247: {
248: PC_LU *dir;
251: dir = (PC_LU*)pc->data;
252: dir->info.damping = damping;
253: dir->info.damp = 1.0;
254: return(0);
255: }
256: EXTERN_C_END
258: EXTERN_C_BEGIN
259: int PCLUSetUseInPlace_LU(PC pc)
260: {
261: PC_LU *dir;
264: dir = (PC_LU*)pc->data;
265: dir->inplace = 1;
266: return(0);
267: }
268: EXTERN_C_END
270: EXTERN_C_BEGIN
271: int PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
272: {
273: PC_LU *dir = (PC_LU*)pc->data;
274: int ierr;
277: PetscStrfree(dir->ordering);
278: PetscStrallocpy(ordering,&dir->ordering);
279: return(0);
280: }
281: EXTERN_C_END
283: EXTERN_C_BEGIN
284: int PCLUSetPivoting_LU(PC pc,PetscReal dtcol)
285: {
286: PC_LU *dir = (PC_LU*)pc->data;
289: if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(1,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
290: dir->info.dtcol = dtcol;
291: return(0);
292: }
293: EXTERN_C_END
295: EXTERN_C_BEGIN
296: int PCLUSetPivotInBlocks_LU(PC pc,PetscTruth pivot)
297: {
298: PC_LU *dir = (PC_LU*)pc->data;
301: dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
302: return(0);
303: }
304: EXTERN_C_END
306: /* -----------------------------------------------------------------------------------*/
308: /*@
309: PCLUSetReuseOrdering - When similar matrices are factored, this
310: causes the ordering computed in the first factor to be used for all
311: following factors; applies to both fill and drop tolerance LUs.
313: Collective on PC
315: Input Parameters:
316: + pc - the preconditioner context
317: - flag - PETSC_TRUE to reuse else PETSC_FALSE
319: Options Database Key:
320: . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
322: Level: intermediate
324: .keywords: PC, levels, reordering, factorization, incomplete, LU
326: .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
327: @*/
328: int PCLUSetReuseOrdering(PC pc,PetscTruth flag)
329: {
330: int ierr,(*f)(PC,PetscTruth);
334: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)(void))&f);
335: if (f) {
336: (*f)(pc,flag);
337: }
338: return(0);
339: }
341: /*@
342: PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
343: this causes later ones to use the fill computed in the initial factorization.
345: Collective on PC
347: Input Parameters:
348: + pc - the preconditioner context
349: - flag - PETSC_TRUE to reuse else PETSC_FALSE
351: Options Database Key:
352: . -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
354: Level: intermediate
356: .keywords: PC, levels, reordering, factorization, incomplete, LU
358: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
359: @*/
360: int PCLUSetReuseFill(PC pc,PetscTruth flag)
361: {
362: int ierr,(*f)(PC,PetscTruth);
366: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)(void))&f);
367: if (f) {
368: (*f)(pc,flag);
369: }
370: return(0);
371: }
373: /*@
374: PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
375: fill = number nonzeros in factor/number nonzeros in original matrix.
377: Collective on PC
378:
379: Input Parameters:
380: + pc - the preconditioner context
381: - fill - amount of expected fill
383: Options Database Key:
384: . -pc_lu_fill <fill> - Sets fill amount
386: Level: intermediate
388: Note:
389: For sparse matrix factorizations it is difficult to predict how much
390: fill to expect. By running with the option -log_info PETSc will print the
391: actual amount of fill used; allowing you to set the value accurately for
392: future runs. Default PETSc uses a value of 5.0
394: .keywords: PC, set, factorization, direct, fill
396: .seealso: PCILUSetFill()
397: @*/
398: int PCLUSetFill(PC pc,PetscReal fill)
399: {
400: int ierr,(*f)(PC,PetscReal);
404: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
405: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)(void))&f);
406: if (f) {
407: (*f)(pc,fill);
408: }
409: return(0);
410: }
412: /*@
413: PCLUSetDamping - adds this quantity to the diagonal of the matrix during the
414: LU numerical factorization
416: Collective on PC
417:
418: Input Parameters:
419: + pc - the preconditioner context
420: - damping - amount of damping
422: Options Database Key:
423: . -pc_lu_damping <damping> - Sets damping amount
425: Level: intermediate
427: .keywords: PC, set, factorization, direct, fill
429: .seealso: PCILUSetFill(), PCILUSetDamp()
430: @*/
431: int PCLUSetDamping(PC pc,PetscReal damping)
432: {
433: int ierr,(*f)(PC,PetscReal);
437: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)(void))&f);
438: if (f) {
439: (*f)(pc,damping);
440: }
441: return(0);
442: }
444: /*@
445: PCLUSetUseInPlace - Tells the system to do an in-place factorization.
446: For dense matrices, this enables the solution of much larger problems.
447: For sparse matrices the factorization cannot be done truly in-place
448: so this does not save memory during the factorization, but after the matrix
449: is factored, the original unfactored matrix is freed, thus recovering that
450: space.
452: Collective on PC
454: Input Parameters:
455: . pc - the preconditioner context
457: Options Database Key:
458: . -pc_lu_in_place - Activates in-place factorization
460: Notes:
461: PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when
462: a different matrix is provided for the multiply and the preconditioner in
463: a call to SLESSetOperators().
464: This is because the Krylov space methods require an application of the
465: matrix multiplication, which is not possible here because the matrix has
466: been factored in-place, replacing the original matrix.
468: Level: intermediate
470: .keywords: PC, set, factorization, direct, inplace, in-place, LU
472: .seealso: PCILUSetUseInPlace()
473: @*/
474: int PCLUSetUseInPlace(PC pc)
475: {
476: int ierr,(*f)(PC);
480: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)(void))&f);
481: if (f) {
482: (*f)(pc);
483: }
484: return(0);
485: }
487: /*@C
488: PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to
489: be used in the LU factorization.
491: Collective on PC
493: Input Parameters:
494: + pc - the preconditioner context
495: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
497: Options Database Key:
498: . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
500: Level: intermediate
502: Notes: nested dissection is used by default
504: .seealso: PCILUSetMatOrdering()
505: @*/
506: int PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
507: {
508: int ierr,(*f)(PC,MatOrderingType);
511: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)(void))&f);
512: if (f) {
513: (*f)(pc,ordering);
514: }
515: return(0);
516: }
518: /*@
519: PCLUSetPivoting - Determines when pivoting is done during LU.
520: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
521: it is never done. For the Matlab and SuperLU factorization this is used.
523: Collective on PC
525: Input Parameters:
526: + pc - the preconditioner context
527: - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable)
529: Options Database Key:
530: . -pc_lu_pivoting - dttol
532: Level: intermediate
534: .seealso: PCILUSetMatOrdering(), PCLUSetPivotInBlocks()
535: @*/
536: int PCLUSetPivoting(PC pc,PetscReal dtcol)
537: {
538: int ierr,(*f)(PC,PetscReal);
541: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivoting_C",(void (**)(void))&f);
542: if (f) {
543: (*f)(pc,dtcol);
544: }
545: return(0);
546: }
548: /*@
549: PCLUSetPivotInBlocks - Determines if pivoting is done while factoring each block
550: with BAIJ or SBAIJ matrices
552: Collective on PC
554: Input Parameters:
555: + pc - the preconditioner context
556: - pivot - PETSC_TRUE or PETSC_FALSE
558: Options Database Key:
559: . -pc_lu_pivot_in_blocks <true,false>
561: Level: intermediate
563: .seealso: PCILUSetMatOrdering(), PCLUSetPivoting()
564: @*/
565: int PCLUSetPivotInBlocks(PC pc,PetscTruth pivot)
566: {
567: int ierr,(*f)(PC,PetscTruth);
570: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetPivotInBlocks_C",(void (**)(void))&f);
571: if (f) {
572: (*f)(pc,pivot);
573: }
574: return(0);
575: }
577: /* ------------------------------------------------------------------------ */
579: EXTERN_C_BEGIN
580: int PCCreate_LU(PC pc)
581: {
582: int ierr,size;
583: PC_LU *dir;
586: PetscNew(PC_LU,&dir);
587: PetscLogObjectMemory(pc,sizeof(PC_LU));
589: dir->fact = 0;
590: dir->inplace = 0;
591: dir->info.fill = 5.0;
592: dir->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */
593: dir->info.damping = 0.0;
594: dir->info.damp = 0.0;
595: dir->info.zeropivot = 1.e-12;
596: dir->info.pivotinblocks = 1.0;
597: dir->col = 0;
598: dir->row = 0;
599: MPI_Comm_size(pc->comm,&size);
600: if (size == 1) {
601: PetscStrallocpy(MATORDERING_ND,&dir->ordering);
602: } else {
603: PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
604: }
605: dir->reusefill = PETSC_FALSE;
606: dir->reuseordering = PETSC_FALSE;
607: pc->data = (void*)dir;
609: pc->ops->destroy = PCDestroy_LU;
610: pc->ops->apply = PCApply_LU;
611: pc->ops->applytranspose = PCApplyTranspose_LU;
612: pc->ops->setup = PCSetUp_LU;
613: pc->ops->setfromoptions = PCSetFromOptions_LU;
614: pc->ops->view = PCView_LU;
615: pc->ops->applyrichardson = 0;
616: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
618: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
619: PCLUSetFill_LU);
620: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
621: PCLUSetDamping_LU);
622: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
623: PCLUSetUseInPlace_LU);
624: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
625: PCLUSetMatOrdering_LU);
626: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
627: PCLUSetReuseOrdering_LU);
628: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
629: PCLUSetReuseFill_LU);
630: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivoting_C","PCLUSetPivoting_LU",
631: PCLUSetPivoting_LU);
632: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetPivotInBlocks_C","PCLUSetPivotInBlocks_LU",
633: PCLUSetPivotInBlocks_LU);
634: return(0);
635: }
636: EXTERN_C_END