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/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 Cholesky */
17: MatFactorInfo info;
18: } PC_Cholesky;
20: EXTERN_C_BEGIN
23: int PCCholeskySetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
24: {
25: PC_Cholesky *lu;
26:
28: lu = (PC_Cholesky*)pc->data;
29: lu->reuseordering = flag;
30: return(0);
31: }
32: EXTERN_C_END
34: EXTERN_C_BEGIN
37: int PCCholeskySetReuseFill_Cholesky(PC pc,PetscTruth flag)
38: {
39: PC_Cholesky *lu;
40:
42: lu = (PC_Cholesky*)pc->data;
43: lu->reusefill = flag;
44: return(0);
45: }
46: EXTERN_C_END
50: static int PCSetFromOptions_Cholesky(PC pc)
51: {
52: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
53: int ierr;
54: PetscTruth flg;
55: char tname[256];
56: PetscFList ordlist;
57:
59: MatOrderingRegisterAll(PETSC_NULL);
60: PetscOptionsHead("Cholesky options");
61: PetscOptionsName("-pc_cholesky_in_place","Form Cholesky in the same memory as the matrix","PCCholeskySetUseInPlace",&flg);
62: if (flg) {
63: PCCholeskySetUseInPlace(pc);
64: }
65: PetscOptionsReal("-pc_cholesky_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCCholeskySetFill",lu->info.fill,&lu->info.fill,0);
66:
67: PetscOptionsName("-pc_cholesky_reuse_fill","Use fill from previous factorization","PCCholeskySetReuseFill",&flg);
68: if (flg) {
69: PCCholeskySetReuseFill(pc,PETSC_TRUE);
70: }
71: PetscOptionsName("-pc_cholesky_reuse_ordering","Reuse ordering from previous factorization","PCCholeskySetReuseOrdering",&flg);
72: if (flg) {
73: PCCholeskySetReuseOrdering(pc,PETSC_TRUE);
74: }
75:
76: MatGetOrderingList(&ordlist);
77: PetscOptionsList("-pc_cholesky_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCCholeskySetMatOrdering",ordlist,lu->ordering,tname,256,&flg);
78: if (flg) {
79: PCCholeskySetMatOrdering(pc,tname);
80: }
81: PetscOptionsName("-pc_cholesky_damping","Damping added to diagonal","PCCholestkySetDamping",&flg);
82: if (flg) {
83: PCCholeskySetDamping(pc,(PetscReal) PETSC_DECIDE);
84: }
85: PetscOptionsReal("-pc_cholesky_damping","Damping added to diagonal","PCCholeskySetDamping",lu->info.damping,&lu->info.damping,0);
86: PetscOptionsName("-pc_cholesky_shift","Manteuffel shift applied to diagonal","PCCholeskySetShift",&flg);
87: if (flg) {
88: printf("doing cholesky shift\n");
89: PCCholeskySetShift(pc,PETSC_TRUE);
90: }
92: PetscOptionsTail();
93: return(0);
94: }
98: static int PCView_Cholesky(PC pc,PetscViewer viewer)
99: {
100: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
101: int ierr;
102: PetscTruth isascii,isstring;
103:
105: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
106: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
107: if (isascii) {
108: MatInfo info;
109:
110: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");}
111: else {PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");}
112: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);
113: if (lu->fact) {
114: MatGetInfo(lu->fact,MAT_LOCAL,&info);
115: PetscViewerASCIIPrintf(viewer," Cholesky nonzeros %g\n",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 factorization\n");}
121: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
122: } else if (isstring) {
123: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
124: } else {
125: SETERRQ1(1,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
126: }
127: return(0);
128: }
132: static int PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
133: {
134: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
135:
137: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
138: *mat = dir->fact;
139: return(0);
140: }
144: static int PCSetUp_Cholesky(PC pc)
145: {
146: int ierr;
147: PetscTruth flg;
148: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
151: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
152:
153: if (dir->inplace) {
154: if (dir->row && dir->col && (dir->row != dir->col)) {
155: ISDestroy(dir->row);
156: dir->row = 0;
157: }
158: if (dir->col) {
159: ISDestroy(dir->col);
160: dir->col = 0;
161: }
162: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
163: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
164: ISDestroy(dir->col);
165: dir->col=0;
166: }
167: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
168: MatCholeskyFactor(pc->pmat,dir->row,&dir->info);
169: dir->fact = pc->pmat;
170: } else {
171: MatInfo info;
172: if (!pc->setupcalled) {
173: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
174: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
175: ISDestroy(dir->col);
176: dir->col=0;
177: }
178: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
179: if (flg) {
180: PetscReal tol = 1.e-10;
181: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
182: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
183: }
184: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
185: MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
186: MatGetInfo(dir->fact,MAT_LOCAL,&info);
187: dir->actualfill = info.fill_ratio_needed;
188: PetscLogObjectParent(pc,dir->fact);
189: } else if (pc->flag != SAME_NONZERO_PATTERN) {
190: if (!dir->reuseordering) {
191: if (dir->row && dir->col && (dir->row != dir->col)) {
192: ISDestroy(dir->row);
193: dir->row = 0;
194: }
195: if (dir->col) {
196: ISDestroy(dir->col);
197: dir->col =0;
198: }
199: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
200: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
201: ISDestroy(dir->col);
202: dir->col=0;
203: }
204: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
205: if (flg) {
206: PetscReal tol = 1.e-10;
207: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
208: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
209: }
210: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
211: }
212: MatDestroy(dir->fact);
213: MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
214: MatGetInfo(dir->fact,MAT_LOCAL,&info);
215: dir->actualfill = info.fill_ratio_needed;
216: PetscLogObjectParent(pc,dir->fact);
217: }
218: MatCholeskyFactorNumeric(pc->pmat,&dir->fact);
219: }
220: return(0);
221: }
225: static int PCDestroy_Cholesky(PC pc)
226: {
227: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
228: int ierr;
231: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
232: if (dir->row) {ISDestroy(dir->row);}
233: if (dir->col) {ISDestroy(dir->col);}
234: PetscStrfree(dir->ordering);
235: PetscFree(dir);
236: return(0);
237: }
241: static int PCApply_Cholesky(PC pc,Vec x,Vec y)
242: {
243: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
244: int ierr;
245:
247: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
248: else {MatSolve(dir->fact,x,y);}
249: return(0);
250: }
254: static int PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
255: {
256: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
257: int ierr;
260: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
261: else {MatSolveTranspose(dir->fact,x,y);}
262: return(0);
263: }
265: /* -----------------------------------------------------------------------------------*/
267: EXTERN_C_BEGIN
270: int PCCholeskySetFill_Cholesky(PC pc,PetscReal fill)
271: {
272: PC_Cholesky *dir;
273:
275: dir = (PC_Cholesky*)pc->data;
276: dir->info.fill = fill;
277: return(0);
278: }
279: EXTERN_C_END
281: EXTERN_C_BEGIN
284: int PCCholeskySetDamping_Cholesky(PC pc,PetscReal damping)
285: {
286: PC_Cholesky *dir;
287:
289: dir = (PC_Cholesky*)pc->data;
290: if (damping == (PetscReal) PETSC_DECIDE) {
291: dir->info.damping = 1.e-12;
292: } else {
293: dir->info.damping = damping;
294: }
295: return(0);
296: }
297: EXTERN_C_END
299: EXTERN_C_BEGIN
302: int PCCholeskySetShift_Cholesky(PC pc,PetscTruth shift)
303: {
304: PC_Cholesky *dir;
305:
307: dir = (PC_Cholesky*)pc->data;
308: dir->info.shift = shift;
309: if (shift) dir->info.shift_fraction = 0.0;
310: return(0);
311: }
312: EXTERN_C_END
314: EXTERN_C_BEGIN
317: int PCCholeskySetUseInPlace_Cholesky(PC pc)
318: {
319: PC_Cholesky *dir;
322: dir = (PC_Cholesky*)pc->data;
323: dir->inplace = 1;
324: return(0);
325: }
326: EXTERN_C_END
328: EXTERN_C_BEGIN
331: int PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
332: {
333: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
334: int ierr;
337: PetscStrfree(dir->ordering);
338: PetscStrallocpy(ordering,&dir->ordering);
339: return(0);
340: }
341: EXTERN_C_END
343: /* -----------------------------------------------------------------------------------*/
347: /*@
348: PCCholeskySetReuseOrdering - When similar matrices are factored, this
349: causes the ordering computed in the first factor to be used for all
350: following factors.
352: Collective on PC
354: Input Parameters:
355: + pc - the preconditioner context
356: - flag - PETSC_TRUE to reuse else PETSC_FALSE
358: Options Database Key:
359: . -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
361: Level: intermediate
363: .keywords: PC, levels, reordering, factorization, incomplete, LU
365: .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
366: @*/
367: int PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
368: {
369: int ierr,(*f)(PC,PetscTruth);
373: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)(void))&f);
374: if (f) {
375: (*f)(pc,flag);
376: }
377: return(0);
378: }
382: /*@
383: PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
384: this causes later ones to use the fill computed in the initial factorization.
386: Collective on PC
388: Input Parameters:
389: + pc - the preconditioner context
390: - flag - PETSC_TRUE to reuse else PETSC_FALSE
392: Options Database Key:
393: . -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
395: Level: intermediate
397: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
399: .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
400: @*/
401: int PCCholeskySetReuseFill(PC pc,PetscTruth flag)
402: {
403: int ierr,(*f)(PC,PetscTruth);
407: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)(void))&f);
408: if (f) {
409: (*f)(pc,flag);
410: }
411: return(0);
412: }
416: /*@
417: PCCholeskySetFill - Indicates the amount of fill you expect in the factored matrix,
418: fill = number nonzeros in factor/number nonzeros in original matrix.
420: Collective on PC
421:
422: Input Parameters:
423: + pc - the preconditioner context
424: - fill - amount of expected fill
426: Options Database Key:
427: . -pc_cholesky_fill <fill> - Sets fill amount
429: Level: intermediate
431: Note:
432: For sparse matrix factorizations it is difficult to predict how much
433: fill to expect. By running with the option -log_info PETSc will print the
434: actual amount of fill used; allowing you to set the value accurately for
435: future runs. Default PETSc uses a value of 5.0
437: .keywords: PC, set, factorization, direct, fill
439: .seealso: PCILUSetFill()
440: @*/
441: int PCCholeskySetFill(PC pc,PetscReal fill)
442: {
443: int ierr,(*f)(PC,PetscReal);
447: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
448: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetFill_C",(void (**)(void))&f);
449: if (f) {
450: (*f)(pc,fill);
451: }
452: return(0);
453: }
457: /*@
458: PCCholeskySetDamping - Adds this quantity to the diagonal of the matrix during the
459: Cholesky numerical factorization.
461: Collective on PC
462:
463: Input Parameters:
464: + pc - the preconditioner context
465: - damping - amount of damping
467: Options Database Key:
468: . -pc_cholesky_damping <damping> - Sets damping amount
470: Level: intermediate
472: .keywords: PC, set, factorization, direct, fill
474: .seealso: PCCholeskySetFill(), PCILUSetDamping()
475: @*/
476: int PCCholeskySetDamping(PC pc,PetscReal damping)
477: {
478: int ierr,(*f)(PC,PetscReal);
482: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetDamping_C",(void (**)(void))&f);
483: if (f) {
484: (*f)(pc,damping);
485: }
486: return(0);
487: }
491: /*@
492: PCCholeskySetShift - specify whether to use Manteuffel shifting of Cholesky.
493: If an Cholesky factorisation breaks down because of nonpositive pivots,
494: adding sufficient identity to the diagonal will remedy this.
495: Setting this causes a bisection method to find the minimum shift that
496: will lead to a well-defined Cholesky.
498: Input parameters:
499: + pc - the preconditioner context
500: - shifting - PETSC_TRUE to set shift else PETSC_FALSE
502: Options Database Key:
503: . -pc_ilu_shift - Activate PCCholeskySetShift()
505: Level: intermediate
507: .keywords: PC, indefinite, factorization, incomplete, Cholesky
509: .seealso: PCILUSetShift()
510: @*/
511: int PCCholeskySetShift(PC pc,PetscTruth shift)
512: {
513: int ierr,(*f)(PC,PetscTruth);
517: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetShift_C",(void (**)(void))&f);
518: if (f) {
519: (*f)(pc,shift);
520: }
521: return(0);
522: }
526: /*@
527: PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
528: For dense matrices, this enables the solution of much larger problems.
529: For sparse matrices the factorization cannot be done truly in-place
530: so this does not save memory during the factorization, but after the matrix
531: is factored, the original unfactored matrix is freed, thus recovering that
532: space.
534: Collective on PC
536: Input Parameters:
537: . pc - the preconditioner context
539: Options Database Key:
540: . -pc_cholesky_in_place - Activates in-place factorization
542: Notes:
543: PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
544: a different matrix is provided for the multiply and the preconditioner in
545: a call to KSPSetOperators().
546: This is because the Krylov space methods require an application of the
547: matrix multiplication, which is not possible here because the matrix has
548: been factored in-place, replacing the original matrix.
550: Level: intermediate
552: .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
554: .seealso: PCICholeskySetUseInPlace()
555: @*/
556: int PCCholeskySetUseInPlace(PC pc)
557: {
558: int ierr,(*f)(PC);
562: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)(void))&f);
563: if (f) {
564: (*f)(pc);
565: }
566: return(0);
567: }
571: /*@
572: PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
573: be used it the Cholesky factorization.
575: Collective on PC
577: Input Parameters:
578: + pc - the preconditioner context
579: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
581: Options Database Key:
582: . -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
584: Level: intermediate
586: .seealso: PCICholeskySetMatOrdering()
587: @*/
588: int PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
589: {
590: int ierr,(*f)(PC,MatOrderingType);
593: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)(void))&f);
594: if (f) {
595: (*f)(pc,ordering);
596: }
597: return(0);
598: }
600: /*MC
601: PCCholesky - Uses a direct solver, based on Cholesky factorization, as a preconditioner
603: Options Database Keys:
604: + -pc_cholesky_reuse_ordering - Activate PCLUSetReuseOrdering()
605: . -pc_cholesky_reuse_fill - Activates PCLUSetReuseFill()
606: . -pc_cholesky_fill <fill> - Sets fill amount
607: . -pc_cholesky_damping <damping> - Sets damping amount
608: . -pc_cholesky_shift - Activates Manteuffel shift
609: . -pc_cholesky_in_place - Activates in-place factorization
610: - -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
612: Notes: Not all options work for all matrix formats
614: Level: beginner
616: Concepts: Cholesky factorization, direct solver
618: Notes: Usually this will compute an "exact" solution in one iteration and does
619: not need a Krylov method (i.e. you can use -ksp_type preonly, or
620: KSPSetType(ksp,KSPPREONLY) for the Krylov method
622: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
623: PCILU, PCLU, PCICC, PCCholeskySetReuseOrdering(), PCCholeskySetReuseFill(), PCGetFactoredMatrix(),
624: PCCholeskySetFill(), PCCholeskySetDamping(), PCCholeskySetShift(),
625: PCCholeskySetUseInPlace(), PCCholeskySetMatOrdering()
627: M*/
629: EXTERN_C_BEGIN
632: int PCCreate_Cholesky(PC pc)
633: {
634: int ierr;
635: PC_Cholesky *dir;
638: PetscNew(PC_Cholesky,&dir);
639: PetscLogObjectMemory(pc,sizeof(PC_Cholesky));
641: dir->fact = 0;
642: dir->inplace = 0;
643: dir->info.fill = 5.0;
644: dir->info.damping = 0.0;
645: dir->info.shift = PETSC_FALSE;
646: dir->info.shift_fraction = 0.0;
647: dir->info.pivotinblocks = 1.0;
648: dir->col = 0;
649: dir->row = 0;
650: PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
651: dir->reusefill = PETSC_FALSE;
652: dir->reuseordering = PETSC_FALSE;
653: pc->data = (void*)dir;
655: pc->ops->destroy = PCDestroy_Cholesky;
656: pc->ops->apply = PCApply_Cholesky;
657: pc->ops->applytranspose = PCApplyTranspose_Cholesky;
658: pc->ops->setup = PCSetUp_Cholesky;
659: pc->ops->setfromoptions = PCSetFromOptions_Cholesky;
660: pc->ops->view = PCView_Cholesky;
661: pc->ops->applyrichardson = 0;
662: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;
664: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetFill_C","PCCholeskySetFill_Cholesky",
665: PCCholeskySetFill_Cholesky);
666: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetDamping_C","PCCholeskySetDamping_Cholesky",
667: PCCholeskySetDamping_Cholesky);
668: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetShift_C","PCCholeskySetShift_Cholesky",
669: PCCholeskySetShift_Cholesky);
670: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetUseInPlace_C","PCCholeskySetUseInPlace_Cholesky",
671: PCCholeskySetUseInPlace_Cholesky);
672: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetMatOrdering_C","PCCholeskySetMatOrdering_Cholesky",
673: PCCholeskySetMatOrdering_Cholesky);
674: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseOrdering_C","PCCholeskySetReuseOrdering_Cholesky",
675: PCCholeskySetReuseOrdering_Cholesky);
676: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseFill_C","PCCholeskySetReuseFill_Cholesky",
677: PCCholeskySetReuseFill_Cholesky);
678: return(0);
679: }
680: EXTERN_C_END