Actual source code: cholesky.c
1: /*$Id: cholesky.c,v 1.12 2001/08/06 21:16:36 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 Cholesky */
17: MatCholeskyInfo info;
18: } PC_Cholesky;
20: EXTERN_C_BEGIN
21: int PCCholeskySetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
22: {
23: PC_Cholesky *lu;
24:
26: lu = (PC_Cholesky*)pc->data;
27: lu->reuseordering = flag;
28: return(0);
29: }
30: EXTERN_C_END
32: EXTERN_C_BEGIN
33: int PCCholeskySetReuseFill_Cholesky(PC pc,PetscTruth flag)
34: {
35: PC_Cholesky *lu;
36:
38: lu = (PC_Cholesky*)pc->data;
39: lu->reusefill = flag;
40: return(0);
41: }
42: EXTERN_C_END
44: static int PCSetFromOptions_Cholesky(PC pc)
45: {
46: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
47: int ierr;
48: PetscTruth flg;
49: char tname[256];
50: PetscFList ordlist;
51:
53: MatOrderingRegisterAll(PETSC_NULL);
54: PetscOptionsHead("Cholesky options");
55: PetscOptionsName("-pc_cholesky_in_place","Form Cholesky in the same memory as the matrix","PCCholeskySetUseInPlace",&flg);
56: if (flg) {
57: PCCholeskySetUseInPlace(pc);
58: }
59: PetscOptionsReal("-pc_cholesky_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCCholeskySetFill",lu->info.fill,&lu->info.fill,0);
60:
61: PetscOptionsName("-pc_cholesky_reuse_fill","Use fill from previous factorization","PCCholeskySetReuseFill",&flg);
62: if (flg) {
63: PCCholeskySetReuseFill(pc,PETSC_TRUE);
64: }
65: PetscOptionsName("-pc_cholesky_reuse_ordering","Reuse ordering from previous factorization","PCCholeskySetReuseOrdering",&flg);
66: if (flg) {
67: PCCholeskySetReuseOrdering(pc,PETSC_TRUE);
68: }
69:
70: MatGetOrderingList(&ordlist);
71: PetscOptionsList("-pc_cholesky_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCCholeskySetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
72: if (flg) {
73: PCCholeskySetMatOrdering(pc,tname);
74: }
75:
76: PetscOptionsTail();
77: return(0);
78: }
80: static int PCView_Cholesky(PC pc,PetscViewer viewer)
81: {
82: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
83: int ierr;
84: PetscTruth isascii,isstring;
85:
87: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
88: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
89: if (isascii) {
90: MatInfo info;
91:
92: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorizationn");}
93: else {PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorizationn");}
94: PetscViewerASCIIPrintf(viewer," matrix ordering: %sn",lu->ordering);
95: if (lu->fact) {
96: MatGetInfo(lu->fact,MAT_LOCAL,&info);
97: PetscViewerASCIIPrintf(viewer," Cholesky nonzeros %gn",info.nz_used);
98: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
99: MatView(lu->fact,viewer);
100: PetscViewerPopFormat(viewer);
101: }
102: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorizationn");}
103: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorizationn");}
104: } else if (isstring) {
105: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
106: } else {
107: SETERRQ1(1,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
108: }
109: return(0);
110: }
112: static int PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
113: {
114: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
115:
117: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
118: *mat = dir->fact;
119: return(0);
120: }
122: static int PCSetUp_Cholesky(PC pc)
123: {
124: int ierr;
125: PetscTruth flg;
126: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
129: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
130:
131: if (dir->inplace) {
132: if (dir->row && dir->col && (dir->row != dir->col)) {
133: ISDestroy(dir->row);
134: dir->row = 0;
135: }
136: if (dir->col) {
137: ISDestroy(dir->col);
138: dir->col = 0;
139: }
140: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
141: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
142: ISDestroy(dir->col);
143: dir->col=0;
144: }
145: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
146: MatCholeskyFactor(pc->pmat,dir->row,dir->info.fill);
147: dir->fact = pc->pmat;
148: } else {
149: MatInfo info;
150: if (!pc->setupcalled) {
151: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
152: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
153: ISDestroy(dir->col);
154: dir->col=0;
155: }
156: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
157: if (flg) {
158: PetscReal tol = 1.e-10;
159: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
160: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
161: }
162: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
163: MatCholeskyFactorSymbolic(pc->pmat,dir->row,dir->info.fill,&dir->fact);
164: MatGetInfo(dir->fact,MAT_LOCAL,&info);
165: dir->actualfill = info.fill_ratio_needed;
166: PetscLogObjectParent(pc,dir->fact);
167: } else if (pc->flag != SAME_NONZERO_PATTERN) {
168: if (!dir->reuseordering) {
169: if (dir->row && dir->col && (dir->row != dir->col)) {
170: ISDestroy(dir->row);
171: dir->row = 0;
172: }
173: if (dir->col) {
174: ISDestroy(dir->col);
175: dir->col =0;
176: }
177: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
178: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
179: ISDestroy(dir->col);
180: dir->col=0;
181: }
182: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
183: if (flg) {
184: PetscReal tol = 1.e-10;
185: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
186: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
187: }
188: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
189: }
190: MatDestroy(dir->fact);
191: MatCholeskyFactorSymbolic(pc->pmat,dir->row,dir->info.fill,&dir->fact);
192: MatGetInfo(dir->fact,MAT_LOCAL,&info);
193: dir->actualfill = info.fill_ratio_needed;
194: PetscLogObjectParent(pc,dir->fact);
195: }
196: MatCholeskyFactorNumeric(pc->pmat,&dir->fact);
197: }
198: return(0);
199: }
201: static int PCDestroy_Cholesky(PC pc)
202: {
203: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
204: int ierr;
207: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
208: if (dir->row) {ISDestroy(dir->row);}
209: if (dir->col) {ISDestroy(dir->col);}
210: PetscStrfree(dir->ordering);
211: PetscFree(dir);
212: return(0);
213: }
215: static int PCApply_Cholesky(PC pc,Vec x,Vec y)
216: {
217: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
218: int ierr;
219:
221: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
222: else {MatSolve(dir->fact,x,y);}
223: return(0);
224: }
226: static int PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
227: {
228: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
229: int ierr;
232: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
233: else {MatSolveTranspose(dir->fact,x,y);}
234: return(0);
235: }
237: /* -----------------------------------------------------------------------------------*/
239: EXTERN_C_BEGIN
240: int PCCholeskySetFill_Cholesky(PC pc,PetscReal fill)
241: {
242: PC_Cholesky *dir;
243:
245: dir = (PC_Cholesky*)pc->data;
246: dir->info.fill = fill;
247: return(0);
248: }
249: EXTERN_C_END
251: EXTERN_C_BEGIN
252: int PCCholeskySetDamping_Cholesky(PC pc,PetscReal damping)
253: {
254: PC_Cholesky *dir;
255:
257: dir = (PC_Cholesky*)pc->data;
258: dir->info.damping = damping;
259: dir->info.damp = 1.0;
260: return(0);
261: }
262: EXTERN_C_END
264: EXTERN_C_BEGIN
265: int PCCholeskySetUseInPlace_Cholesky(PC pc)
266: {
267: PC_Cholesky *dir;
270: dir = (PC_Cholesky*)pc->data;
271: dir->inplace = 1;
272: return(0);
273: }
274: EXTERN_C_END
276: EXTERN_C_BEGIN
277: int PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
278: {
279: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
280: int ierr;
283: PetscStrfree(dir->ordering);
284: PetscStrallocpy(ordering,&dir->ordering);
285: return(0);
286: }
287: EXTERN_C_END
289: /* -----------------------------------------------------------------------------------*/
291: /*@
292: PCCholeskySetReuseOrdering - When similar matrices are factored, this
293: causes the ordering computed in the first factor to be used for all
294: following factors.
296: Collective on PC
298: Input Parameters:
299: + pc - the preconditioner context
300: - flag - PETSC_TRUE to reuse else PETSC_FALSE
302: Options Database Key:
303: . -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
305: Level: intermediate
307: .keywords: PC, levels, reordering, factorization, incomplete, LU
309: .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
310: @*/
311: int PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
312: {
313: int ierr,(*f)(PC,PetscTruth);
317: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)(void))&f);
318: if (f) {
319: (*f)(pc,flag);
320: }
321: return(0);
322: }
324: /*@
325: PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
326: this causes later ones to use the fill computed in the initial factorization.
328: Collective on PC
330: Input Parameters:
331: + pc - the preconditioner context
332: - flag - PETSC_TRUE to reuse else PETSC_FALSE
334: Options Database Key:
335: . -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
337: Level: intermediate
339: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
341: .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
342: @*/
343: int PCCholeskySetReuseFill(PC pc,PetscTruth flag)
344: {
345: int ierr,(*f)(PC,PetscTruth);
349: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)(void))&f);
350: if (f) {
351: (*f)(pc,flag);
352: }
353: return(0);
354: }
356: /*@
357: PCCholeskySetFill - Indicates the amount of fill you expect in the factored matrix,
358: fill = number nonzeros in factor/number nonzeros in original matrix.
360: Collective on PC
361:
362: Input Parameters:
363: + pc - the preconditioner context
364: - fill - amount of expected fill
366: Options Database Key:
367: . -pc_cholesky_fill <fill> - Sets fill amount
369: Level: intermediate
371: Note:
372: For sparse matrix factorizations it is difficult to predict how much
373: fill to expect. By running with the option -log_info PETSc will print the
374: actual amount of fill used; allowing you to set the value accurately for
375: future runs. Default PETSc uses a value of 5.0
377: .keywords: PC, set, factorization, direct, fill
379: .seealso: PCILUSetFill()
380: @*/
381: int PCCholeskySetFill(PC pc,PetscReal fill)
382: {
383: int ierr,(*f)(PC,PetscReal);
387: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
388: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetFill_C",(void (**)(void))&f);
389: if (f) {
390: (*f)(pc,fill);
391: }
392: return(0);
393: }
395: /*@
396: PCCholeskySetDamping - Adds this quantity to the diagonal of the matrix during the
397: Cholesky numerical factorization.
399: Collective on PC
400:
401: Input Parameters:
402: + pc - the preconditioner context
403: - damping - amount of damping
405: Options Database Key:
406: . -pc_cholesky_damping <damping> - Sets damping amount
408: Level: intermediate
410: .keywords: PC, set, factorization, direct, fill
412: .seealso: PCICholeskySetFill(), PCILUSetDamp()
413: @*/
414: int PCCholeskySetDamping(PC pc,PetscReal damping)
415: {
416: int ierr,(*f)(PC,PetscReal);
420: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetDamping_C",(void (**)(void))&f);
421: if (f) {
422: (*f)(pc,damping);
423: }
424: return(0);
425: }
427: /*@
428: PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
429: For dense matrices, this enables the solution of much larger problems.
430: For sparse matrices the factorization cannot be done truly in-place
431: so this does not save memory during the factorization, but after the matrix
432: is factored, the original unfactored matrix is freed, thus recovering that
433: space.
435: Collective on PC
437: Input Parameters:
438: . pc - the preconditioner context
440: Options Database Key:
441: . -pc_cholesky_in_place - Activates in-place factorization
443: Notes:
444: PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
445: a different matrix is provided for the multiply and the preconditioner in
446: a call to SLESSetOperators().
447: This is because the Krylov space methods require an application of the
448: matrix multiplication, which is not possible here because the matrix has
449: been factored in-place, replacing the original matrix.
451: Level: intermediate
453: .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
455: .seealso: PCICholeskySetUseInPlace()
456: @*/
457: int PCCholeskySetUseInPlace(PC pc)
458: {
459: int ierr,(*f)(PC);
463: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)(void))&f);
464: if (f) {
465: (*f)(pc);
466: }
467: return(0);
468: }
470: /*@
471: PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
472: be used it the Cholesky factorization.
474: Collective on PC
476: Input Parameters:
477: + pc - the preconditioner context
478: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
480: Options Database Key:
481: . -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
483: Level: intermediate
485: .seealso: PCICholeskySetMatOrdering()
486: @*/
487: int PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
488: {
489: int ierr,(*f)(PC,MatOrderingType);
492: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)(void))&f);
493: if (f) {
494: (*f)(pc,ordering);
495: }
496: return(0);
497: }
499: /* ------------------------------------------------------------------------ */
501: EXTERN_C_BEGIN
502: int PCCreate_Cholesky(PC pc)
503: {
504: int ierr;
505: PC_Cholesky *dir;
508: PetscNew(PC_Cholesky,&dir);
509: PetscLogObjectMemory(pc,sizeof(PC_Cholesky));
511: dir->fact = 0;
512: dir->inplace = 0;
513: dir->info.fill = 5.0;
514: dir->info.damping = 0.0;
515: dir->info.damp = 0.0;
516: dir->info.pivotinblocks = 1.0;
517: dir->col = 0;
518: dir->row = 0;
519: PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
520: dir->reusefill = PETSC_FALSE;
521: dir->reuseordering = PETSC_FALSE;
522: pc->data = (void*)dir;
524: pc->ops->destroy = PCDestroy_Cholesky;
525: pc->ops->apply = PCApply_Cholesky;
526: pc->ops->applytranspose = PCApplyTranspose_Cholesky;
527: pc->ops->setup = PCSetUp_Cholesky;
528: pc->ops->setfromoptions = PCSetFromOptions_Cholesky;
529: pc->ops->view = PCView_Cholesky;
530: pc->ops->applyrichardson = 0;
531: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;
533: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetFill_C","PCCholeskySetFill_Cholesky",
534: PCCholeskySetFill_Cholesky);
535: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetDamping_C","PCCholeskySetDamping_Cholesky",
536: PCCholeskySetDamping_Cholesky);
537: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetUseInPlace_C","PCCholeskySetUseInPlace_Cholesky",
538: PCCholeskySetUseInPlace_Cholesky);
539: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetMatOrdering_C","PCCholeskySetMatOrdering_Cholesky",
540: PCCholeskySetMatOrdering_Cholesky);
541: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseOrdering_C","PCCholeskySetReuseOrdering_Cholesky",
542: PCCholeskySetReuseOrdering_Cholesky);
543: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseFill_C","PCCholeskySetReuseFill_Cholesky",
544: PCCholeskySetReuseFill_Cholesky);
545: return(0);
546: }
547: EXTERN_C_END