Actual source code: lu.c
1: /*$Id: lu.c,v 1.144 2001/04/10 19:36:09 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" /*I "petscpc.h" I*/
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;
50: char tname[256];
51: PetscFList ordlist;
54: MatOrderingRegisterAll(PETSC_NULL);
55: PetscOptionsHead("LU options");
56: PetscOptionsName("-pc_lu_in_place","Form LU in the same memory as the matrix","PCLUSetUseInPlace",&flg);
57: if (flg) {
58: PCLUSetUseInPlace(pc);
59: }
60: PetscOptionsDouble("-pc_lu_fill","Expected non-zeros in LU/non-zeros in matrix","PCLUSetFill",lu->info.fill,&lu->info.fill,0);
62: PetscOptionsHasName(pc->prefix,"-pc_lu_damping",&flg);
63: if (flg) {
64: PCLUSetDamping(pc,0.0);
65: }
66: PetscOptionsDouble("-pc_lu_damping","Damping added to diagonal","PCLUSetDamping",lu->info.damping,&lu->info.damping,0);
68: PetscOptionsName("-pc_lu_reuse_fill","Use fill from previous factorization","PCLUSetReuseFill",&flg);
69: if (flg) {
70: PCLUSetReuseFill(pc,PETSC_TRUE);
71: }
72: PetscOptionsName("-pc_lu_reuse_ordering","Reuse ordering from previous factorization","PCLUSetReuseOrdering",&flg);
73: if (flg) {
74: PCLUSetReuseOrdering(pc,PETSC_TRUE);
75: }
77: MatGetOrderingList(&ordlist);
78: PetscOptionsList("-pc_lu_mat_ordering_type","Reordering to reduce nonzeros in LU","PCLUSetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
79: if (flg) {
80: PCLUSetMatOrdering(pc,tname);
81: }
82: PetscOptionsDouble("-pc_lu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,0,0);
84: PetscOptionsDouble("-pc_lu_column_pivoting","Column pivoting tolerance (not used)","PCLUSetColumnPivoting",lu->info.dtcol,&lu->info.dtcol,&flg);
85: PetscOptionsTail();
86: return(0);
87: }
89: static int PCView_LU(PC pc,PetscViewer viewer)
90: {
91: PC_LU *lu = (PC_LU*)pc->data;
92: int ierr;
93: PetscTruth isascii,isstring;
96: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
97: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
98: if (isascii) {
99: MatInfo info;
101: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," LU: in-place factorizationn");}
102: else {PetscViewerASCIIPrintf(viewer," LU: out-of-place factorizationn");}
103: PetscViewerASCIIPrintf(viewer," matrix ordering: %sn",lu->ordering);
104: if (lu->fact) {
105: MatGetInfo(lu->fact,MAT_LOCAL,&info);
106: PetscViewerASCIIPrintf(viewer," LU nonzeros %gn",info.nz_used);
107: }
108: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorizationn");}
109: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorizationn");}
110: } else if (isstring) {
111: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
112: } else {
113: SETERRQ1(1,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name);
114: }
115: return(0);
116: }
118: static int PCGetFactoredMatrix_LU(PC pc,Mat *mat)
119: {
120: PC_LU *dir = (PC_LU*)pc->data;
123: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
124: *mat = dir->fact;
125: return(0);
126: }
128: static int PCSetUp_LU(PC pc)
129: {
130: int ierr;
131: PetscTruth flg;
132: PC_LU *dir = (PC_LU*)pc->data;
135: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
137: if (dir->inplace) {
138: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
139: if (dir->col) {ISDestroy(dir->col);}
140: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
141: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
142: MatLUFactor(pc->pmat,dir->row,dir->col,&dir->info);
143: dir->fact = pc->pmat;
144: } else {
145: MatInfo info;
146: if (!pc->setupcalled) {
147: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
148: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
149: if (flg) {
150: PetscReal tol = 1.e-10;
151: PetscOptionsGetDouble(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
152: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
153: }
154: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
155: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
156: MatGetInfo(dir->fact,MAT_LOCAL,&info);
157: dir->actualfill = info.fill_ratio_needed;
158: PetscLogObjectParent(pc,dir->fact);
159: } else if (pc->flag != SAME_NONZERO_PATTERN) {
160: if (!dir->reuseordering) {
161: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
162: if (dir->col) {ISDestroy(dir->col);}
163: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
164: PetscOptionsHasName(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&flg);
165: if (flg) {
166: PetscReal tol = 1.e-10;
167: PetscOptionsGetDouble(pc->prefix,"-pc_lu_nonzeros_along_diagonal",&tol,PETSC_NULL);
168: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->col);
169: }
170: if (dir->row) {PetscLogObjectParent(pc,dir->row); PetscLogObjectParent(pc,dir->col);}
171: }
172: MatDestroy(dir->fact);
173: MatLUFactorSymbolic(pc->pmat,dir->row,dir->col,&dir->info,&dir->fact);
174: MatGetInfo(dir->fact,MAT_LOCAL,&info);
175: dir->actualfill = info.fill_ratio_needed;
176: PetscLogObjectParent(pc,dir->fact);
177: }
178: MatLUFactorNumeric(pc->pmat,&dir->fact);
179: }
180: return(0);
181: }
183: static int PCDestroy_LU(PC pc)
184: {
185: PC_LU *dir = (PC_LU*)pc->data;
186: int ierr;
189: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
190: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
191: if (dir->col) {ISDestroy(dir->col);}
192: PetscStrfree(dir->ordering);
193: PetscFree(dir);
194: return(0);
195: }
197: static int PCApply_LU(PC pc,Vec x,Vec y)
198: {
199: PC_LU *dir = (PC_LU*)pc->data;
200: int ierr;
203: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
204: else {MatSolve(dir->fact,x,y);}
205: return(0);
206: }
208: static int PCApplyTranspose_LU(PC pc,Vec x,Vec y)
209: {
210: PC_LU *dir = (PC_LU*)pc->data;
211: int ierr;
214: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
215: else {MatSolveTranspose(dir->fact,x,y);}
216: return(0);
217: }
219: /* -----------------------------------------------------------------------------------*/
221: EXTERN_C_BEGIN
222: int PCLUSetFill_LU(PC pc,PetscReal fill)
223: {
224: PC_LU *dir;
227: dir = (PC_LU*)pc->data;
228: dir->info.fill = fill;
229: return(0);
230: }
231: EXTERN_C_END
233: EXTERN_C_BEGIN
234: int PCLUSetDamping_LU(PC pc,PetscReal damping)
235: {
236: PC_LU *dir;
239: dir = (PC_LU*)pc->data;
240: dir->info.damping = damping;
241: dir->info.damp = 1.0;
242: return(0);
243: }
244: EXTERN_C_END
246: EXTERN_C_BEGIN
247: int PCLUSetUseInPlace_LU(PC pc)
248: {
249: PC_LU *dir;
252: dir = (PC_LU*)pc->data;
253: dir->inplace = 1;
254: return(0);
255: }
256: EXTERN_C_END
258: EXTERN_C_BEGIN
259: int PCLUSetMatOrdering_LU(PC pc,MatOrderingType ordering)
260: {
261: PC_LU *dir = (PC_LU*)pc->data;
262: int ierr;
265: PetscStrfree(dir->ordering);
266: PetscStrallocpy(ordering,&dir->ordering);
267: return(0);
268: }
269: EXTERN_C_END
271: EXTERN_C_BEGIN
272: int PCLUSetColumnPivoting_LU(PC pc,PetscReal dtcol)
273: {
274: PC_LU *dir = (PC_LU*)pc->data;
277: if (dtcol < 0.0 || dtcol > 1.0) SETERRQ1(1,"Column pivot tolerance is %g must be between 0 and 1",dtcol);
278: dir->info.dtcol = dtcol;
279: return(0);
280: }
281: EXTERN_C_END
283: /* -----------------------------------------------------------------------------------*/
285: /*@
286: PCLUSetReuseOrdering - When similar matrices are factored, this
287: causes the ordering computed in the first factor to be used for all
288: following factors; applies to both fill and drop tolerance LUs.
290: Collective on PC
292: Input Parameters:
293: + pc - the preconditioner context
294: - flag - PETSC_TRUE to reuse else PETSC_FALSE
296: Options Database Key:
297: . -pc_lu_reuse_ordering - Activate PCLUSetReuseOrdering()
299: Level: intermediate
301: .keywords: PC, levels, reordering, factorization, incomplete, LU
303: .seealso: PCLUSetReuseFill(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill()
304: @*/
305: int PCLUSetReuseOrdering(PC pc,PetscTruth flag)
306: {
307: int ierr,(*f)(PC,PetscTruth);
311: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseOrdering_C",(void (**)())&f);
312: if (f) {
313: (*f)(pc,flag);
314: }
315: return(0);
316: }
318: /*@
319: PCLUSetReuseFill - When matrices with same nonzero structure are LU factored,
320: this causes later ones to use the fill computed in the initial factorization.
322: Collective on PC
324: Input Parameters:
325: + pc - the preconditioner context
326: - flag - PETSC_TRUE to reuse else PETSC_FALSE
328: Options Database Key:
329: . -pc_lu_reuse_fill - Activates PCLUSetReuseFill()
331: Level: intermediate
333: .keywords: PC, levels, reordering, factorization, incomplete, LU
335: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCILUDTSetReuseFill()
336: @*/
337: int PCLUSetReuseFill(PC pc,PetscTruth flag)
338: {
339: int ierr,(*f)(PC,PetscTruth);
343: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetReuseFill_C",(void (**)())&f);
344: if (f) {
345: (*f)(pc,flag);
346: }
347: return(0);
348: }
350: /*@
351: PCLUSetFill - Indicate the amount of fill you expect in the factored matrix,
352: fill = number nonzeros in factor/number nonzeros in original matrix.
354: Collective on PC
355:
356: Input Parameters:
357: + pc - the preconditioner context
358: - fill - amount of expected fill
360: Options Database Key:
361: . -pc_lu_fill <fill> - Sets fill amount
363: Level: intermediate
365: Note:
366: For sparse matrix factorizations it is difficult to predict how much
367: fill to expect. By running with the option -log_info PETSc will print the
368: actual amount of fill used; allowing you to set the value accurately for
369: future runs. Default PETSc uses a value of 5.0
371: .keywords: PC, set, factorization, direct, fill
373: .seealso: PCILUSetFill()
374: @*/
375: int PCLUSetFill(PC pc,PetscReal fill)
376: {
377: int ierr,(*f)(PC,PetscReal);
381: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
382: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetFill_C",(void (**)())&f);
383: if (f) {
384: (*f)(pc,fill);
385: }
386: return(0);
387: }
389: /*@
390: PCLUSetDamping - adds this quantity to the diagonal of the matrix during the
391: LU numerical factorization
393: Collective on PC
394:
395: Input Parameters:
396: + pc - the preconditioner context
397: - damping - amount of damping
399: Options Database Key:
400: . -pc_lu_damping <damping> - Sets damping amount
402: Level: intermediate
404: .keywords: PC, set, factorization, direct, fill
406: .seealso: PCILUSetFill(), PCILUSetDamp()
407: @*/
408: int PCLUSetDamping(PC pc,PetscReal damping)
409: {
410: int ierr,(*f)(PC,PetscReal);
414: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetDamping_C",(void (**)())&f);
415: if (f) {
416: (*f)(pc,damping);
417: }
418: return(0);
419: }
421: /*@
422: PCLUSetUseInPlace - Tells the system to do an in-place factorization.
423: For dense matrices, this enables the solution of much larger problems.
424: For sparse matrices the factorization cannot be done truly in-place
425: so this does not save memory during the factorization, but after the matrix
426: is factored, the original unfactored matrix is freed, thus recovering that
427: space.
429: Collective on PC
431: Input Parameters:
432: . pc - the preconditioner context
434: Options Database Key:
435: . -pc_lu_in_place - Activates in-place factorization
437: Notes:
438: PCLUSetUseInplace() can only be used with the KSP method KSPPREONLY or when
439: a different matrix is provided for the multiply and the preconditioner in
440: a call to SLESSetOperators().
441: This is because the Krylov space methods require an application of the
442: matrix multiplication, which is not possible here because the matrix has
443: been factored in-place, replacing the original matrix.
445: Level: intermediate
447: .keywords: PC, set, factorization, direct, inplace, in-place, LU
449: .seealso: PCILUSetUseInPlace()
450: @*/
451: int PCLUSetUseInPlace(PC pc)
452: {
453: int ierr,(*f)(PC);
457: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetUseInPlace_C",(void (**)())&f);
458: if (f) {
459: (*f)(pc);
460: }
461: return(0);
462: }
464: /*@
465: PCLUSetMatOrdering - Sets the ordering routine (to reduce fill) to
466: be used in the LU factorization.
468: Collective on PC
470: Input Parameters:
471: + pc - the preconditioner context
472: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
474: Options Database Key:
475: . -pc_lu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
477: Level: intermediate
479: Notes: nested dissection is used by default
481: .seealso: PCILUSetMatOrdering()
482: @*/
483: int PCLUSetMatOrdering(PC pc,MatOrderingType ordering)
484: {
485: int ierr,(*f)(PC,MatOrderingType);
488: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetMatOrdering_C",(void (**)())&f);
489: if (f) {
490: (*f)(pc,ordering);
491: }
492: return(0);
493: }
495: /*@
496: PCLUSetColumnPivoting - Determines when column pivoting is done during LU.
497: For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices
498: it is never done. For the Matlab factorization this is used.
500: Collective on PC
502: Input Parameters:
503: + pc - the preconditioner context
504: - dtcol - 0.0 implies no pivoting, 1.0 complete column pivoting (slower, requires more memory but more stable)
506: Options Database Key:
507: . -pc_lu_column_pivoting - dttol
509: Level: intermediate
511: .seealso: PCILUSetMatOrdering()
512: @*/
513: int PCLUSetColumnPivoting(PC pc,PetscReal dtcol)
514: {
515: int ierr,(*f)(PC,PetscReal);
518: PetscObjectQueryFunction((PetscObject)pc,"PCLUSetColumnPivoting_C",(void (**)())&f);
519: if (f) {
520: (*f)(pc,dtcol);
521: }
522: return(0);
523: }
525: /* ------------------------------------------------------------------------ */
527: EXTERN_C_BEGIN
528: int PCCreate_LU(PC pc)
529: {
530: int ierr;
531: PC_LU *dir;
534: PetscNew(PC_LU,&dir);
535: PetscLogObjectMemory(pc,sizeof(PC_LU));
537: dir->fact = 0;
538: dir->inplace = 0;
539: dir->info.fill = 5.0;
540: dir->info.dtcol = 0.0; /* default to no pivoting; this is only thing PETSc LU supports */
541: dir->info.damping = 0.0;
542: dir->info.damp = 0.0;
543: dir->col = 0;
544: dir->row = 0;
545: PetscStrallocpy(MATORDERING_ND,&dir->ordering);
546: dir->reusefill = PETSC_FALSE;
547: dir->reuseordering = PETSC_FALSE;
548: pc->data = (void*)dir;
550: pc->ops->destroy = PCDestroy_LU;
551: pc->ops->apply = PCApply_LU;
552: pc->ops->applytranspose = PCApplyTranspose_LU;
553: pc->ops->setup = PCSetUp_LU;
554: pc->ops->setfromoptions = PCSetFromOptions_LU;
555: pc->ops->view = PCView_LU;
556: pc->ops->applyrichardson = 0;
557: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_LU;
559: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetFill_C","PCLUSetFill_LU",
560: PCLUSetFill_LU);
561: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetDamping_C","PCLUSetDamping_LU",
562: PCLUSetDamping_LU);
563: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetUseInPlace_C","PCLUSetUseInPlace_LU",
564: PCLUSetUseInPlace_LU);
565: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetMatOrdering_C","PCLUSetMatOrdering_LU",
566: PCLUSetMatOrdering_LU);
567: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseOrdering_C","PCLUSetReuseOrdering_LU",
568: PCLUSetReuseOrdering_LU);
569: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetReuseFill_C","PCLUSetReuseFill_LU",
570: PCLUSetReuseFill_LU);
571: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCLUSetColumnPivoting_C","PCLUSetColumnPivoting_LU",
572: PCLUSetColumnPivoting_LU);
573: return(0);
574: }
575: EXTERN_C_END