Actual source code: cholesky.c
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: PetscTruth 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;
23: PetscErrorCode 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: }
37: PetscErrorCode 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: }
50: static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
51: {
52: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
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: PCCholeskySetShift(pc,PETSC_TRUE);
89: }
91: PetscOptionsTail();
92: return(0);
93: }
97: static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
98: {
99: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
101: PetscTruth iascii,isstring;
102:
104: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
105: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
106: if (iascii) {
107: MatInfo info;
108:
109: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");}
110: else {PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");}
111: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",lu->ordering);
112: if (lu->fact) {
113: MatGetInfo(lu->fact,MAT_LOCAL,&info);
114: PetscViewerASCIIPrintf(viewer," Cholesky nonzeros %g\n",info.nz_used);
115: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_FACTOR_INFO);
116: MatView(lu->fact,viewer);
117: PetscViewerPopFormat(viewer);
118: }
119: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
120: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
121: } else if (isstring) {
122: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
123: } else {
124: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
125: }
126: return(0);
127: }
131: static PetscErrorCode PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
132: {
133: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
134:
136: if (!dir->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
137: *mat = dir->fact;
138: return(0);
139: }
143: static PetscErrorCode PCSetUp_Cholesky(PC pc)
144: {
146: PetscTruth flg;
147: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
150: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
151:
152: if (dir->inplace) {
153: if (dir->row && dir->col && (dir->row != dir->col)) {
154: ISDestroy(dir->row);
155: dir->row = 0;
156: }
157: if (dir->col) {
158: ISDestroy(dir->col);
159: dir->col = 0;
160: }
161: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
162: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
163: ISDestroy(dir->col);
164: dir->col=0;
165: }
166: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
167: MatCholeskyFactor(pc->pmat,dir->row,&dir->info);
168: dir->fact = pc->pmat;
169: } else {
170: MatInfo info;
171: if (!pc->setupcalled) {
172: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
173: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
174: ISDestroy(dir->col);
175: dir->col=0;
176: }
177: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
178: if (flg) {
179: PetscReal tol = 1.e-10;
180: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
181: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
182: }
183: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
184: MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
185: MatGetInfo(dir->fact,MAT_LOCAL,&info);
186: dir->actualfill = info.fill_ratio_needed;
187: PetscLogObjectParent(pc,dir->fact);
188: } else if (pc->flag != SAME_NONZERO_PATTERN) {
189: if (!dir->reuseordering) {
190: if (dir->row && dir->col && (dir->row != dir->col)) {
191: ISDestroy(dir->row);
192: dir->row = 0;
193: }
194: if (dir->col) {
195: ISDestroy(dir->col);
196: dir->col =0;
197: }
198: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
199: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
200: ISDestroy(dir->col);
201: dir->col=0;
202: }
203: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
204: if (flg) {
205: PetscReal tol = 1.e-10;
206: PetscOptionsGetReal(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
207: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
208: }
209: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
210: }
211: MatDestroy(dir->fact);
212: MatCholeskyFactorSymbolic(pc->pmat,dir->row,&dir->info,&dir->fact);
213: MatGetInfo(dir->fact,MAT_LOCAL,&info);
214: dir->actualfill = info.fill_ratio_needed;
215: PetscLogObjectParent(pc,dir->fact);
216: }
217: MatCholeskyFactorNumeric(pc->pmat,&dir->fact);
218: }
219: return(0);
220: }
224: static PetscErrorCode PCDestroy_Cholesky(PC pc)
225: {
226: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
230: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
231: if (dir->row) {ISDestroy(dir->row);}
232: if (dir->col) {ISDestroy(dir->col);}
233: PetscStrfree(dir->ordering);
234: PetscFree(dir);
235: return(0);
236: }
240: static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
241: {
242: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
244:
246: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
247: else {MatSolve(dir->fact,x,y);}
248: return(0);
249: }
253: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
254: {
255: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
259: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
260: else {MatSolveTranspose(dir->fact,x,y);}
261: return(0);
262: }
264: /* -----------------------------------------------------------------------------------*/
269: PetscErrorCode PCCholeskySetFill_Cholesky(PC pc,PetscReal fill)
270: {
271: PC_Cholesky *dir;
272:
274: dir = (PC_Cholesky*)pc->data;
275: dir->info.fill = fill;
276: return(0);
277: }
283: PetscErrorCode PCCholeskySetDamping_Cholesky(PC pc,PetscReal damping)
284: {
285: PC_Cholesky *dir;
286:
288: dir = (PC_Cholesky*)pc->data;
289: if (damping == (PetscReal) PETSC_DECIDE) {
290: dir->info.damping = 1.e-12;
291: } else {
292: dir->info.damping = damping;
293: }
294: return(0);
295: }
301: PetscErrorCode PCCholeskySetShift_Cholesky(PC pc,PetscTruth shift)
302: {
303: PC_Cholesky *dir;
304:
306: dir = (PC_Cholesky*)pc->data;
307: dir->info.shift = shift;
308: if (shift) dir->info.shift_fraction = 0.0;
309: return(0);
310: }
316: PetscErrorCode PCCholeskySetUseInPlace_Cholesky(PC pc)
317: {
318: PC_Cholesky *dir;
321: dir = (PC_Cholesky*)pc->data;
322: dir->inplace = PETSC_TRUE;
323: return(0);
324: }
330: PetscErrorCode PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
331: {
332: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
336: PetscStrfree(dir->ordering);
337: PetscStrallocpy(ordering,&dir->ordering);
338: return(0);
339: }
342: /* -----------------------------------------------------------------------------------*/
346: /*@
347: PCCholeskySetReuseOrdering - When similar matrices are factored, this
348: causes the ordering computed in the first factor to be used for all
349: following factors.
351: Collective on PC
353: Input Parameters:
354: + pc - the preconditioner context
355: - flag - PETSC_TRUE to reuse else PETSC_FALSE
357: Options Database Key:
358: . -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
360: Level: intermediate
362: .keywords: PC, levels, reordering, factorization, incomplete, LU
364: .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
365: @*/
366: PetscErrorCode PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
367: {
368: PetscErrorCode ierr,(*f)(PC,PetscTruth);
372: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)(void))&f);
373: if (f) {
374: (*f)(pc,flag);
375: }
376: return(0);
377: }
381: /*@
382: PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
383: this causes later ones to use the fill computed in the initial factorization.
385: Collective on PC
387: Input Parameters:
388: + pc - the preconditioner context
389: - flag - PETSC_TRUE to reuse else PETSC_FALSE
391: Options Database Key:
392: . -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
394: Level: intermediate
396: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
398: .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
399: @*/
400: PetscErrorCode PCCholeskySetReuseFill(PC pc,PetscTruth flag)
401: {
402: PetscErrorCode ierr,(*f)(PC,PetscTruth);
406: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)(void))&f);
407: if (f) {
408: (*f)(pc,flag);
409: }
410: return(0);
411: }
415: /*@
416: PCCholeskySetFill - Indicates the amount of fill you expect in the factored matrix,
417: fill = number nonzeros in factor/number nonzeros in original matrix.
419: Collective on PC
420:
421: Input Parameters:
422: + pc - the preconditioner context
423: - fill - amount of expected fill
425: Options Database Key:
426: . -pc_cholesky_fill <fill> - Sets fill amount
428: Level: intermediate
430: Note:
431: For sparse matrix factorizations it is difficult to predict how much
432: fill to expect. By running with the option -log_info PETSc will print the
433: actual amount of fill used; allowing you to set the value accurately for
434: future runs. Default PETSc uses a value of 5.0
436: .keywords: PC, set, factorization, direct, fill
438: .seealso: PCILUSetFill()
439: @*/
440: PetscErrorCode PCCholeskySetFill(PC pc,PetscReal fill)
441: {
442: PetscErrorCode ierr,(*f)(PC,PetscReal);
446: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
447: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetFill_C",(void (**)(void))&f);
448: if (f) {
449: (*f)(pc,fill);
450: }
451: return(0);
452: }
456: /*@
457: PCCholeskySetDamping - Adds this quantity to the diagonal of the matrix during the
458: Cholesky numerical factorization.
460: Collective on PC
461:
462: Input Parameters:
463: + pc - the preconditioner context
464: - damping - amount of damping
466: Options Database Key:
467: . -pc_cholesky_damping <damping> - Sets damping amount
469: Level: intermediate
471: .keywords: PC, set, factorization, direct, fill
473: .seealso: PCCholeskySetFill(), PCILUSetDamping()
474: @*/
475: PetscErrorCode PCCholeskySetDamping(PC pc,PetscReal damping)
476: {
477: PetscErrorCode ierr,(*f)(PC,PetscReal);
481: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetDamping_C",(void (**)(void))&f);
482: if (f) {
483: (*f)(pc,damping);
484: }
485: return(0);
486: }
490: /*@
491: PCCholeskySetShift - specify whether to use Manteuffel shifting of Cholesky.
492: If an Cholesky factorisation breaks down because of nonpositive pivots,
493: adding sufficient identity to the diagonal will remedy this.
494: Setting this causes a bisection method to find the minimum shift that
495: will lead to a well-defined Cholesky.
497: Input parameters:
498: + pc - the preconditioner context
499: - shifting - PETSC_TRUE to set shift else PETSC_FALSE
501: Options Database Key:
502: . -pc_ilu_shift - Activate PCCholeskySetShift()
504: Level: intermediate
506: .keywords: PC, indefinite, factorization, incomplete, Cholesky
508: .seealso: PCILUSetShift()
509: @*/
510: PetscErrorCode PCCholeskySetShift(PC pc,PetscTruth shift)
511: {
512: PetscErrorCode ierr,(*f)(PC,PetscTruth);
516: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetShift_C",(void (**)(void))&f);
517: if (f) {
518: (*f)(pc,shift);
519: }
520: return(0);
521: }
525: /*@
526: PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
527: For dense matrices, this enables the solution of much larger problems.
528: For sparse matrices the factorization cannot be done truly in-place
529: so this does not save memory during the factorization, but after the matrix
530: is factored, the original unfactored matrix is freed, thus recovering that
531: space.
533: Collective on PC
535: Input Parameters:
536: . pc - the preconditioner context
538: Options Database Key:
539: . -pc_cholesky_in_place - Activates in-place factorization
541: Notes:
542: PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
543: a different matrix is provided for the multiply and the preconditioner in
544: a call to KSPSetOperators().
545: This is because the Krylov space methods require an application of the
546: matrix multiplication, which is not possible here because the matrix has
547: been factored in-place, replacing the original matrix.
549: Level: intermediate
551: .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
553: .seealso: PCICholeskySetUseInPlace()
554: @*/
555: PetscErrorCode PCCholeskySetUseInPlace(PC pc)
556: {
557: PetscErrorCode ierr,(*f)(PC);
561: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)(void))&f);
562: if (f) {
563: (*f)(pc);
564: }
565: return(0);
566: }
570: /*@
571: PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
572: be used it the Cholesky factorization.
574: Collective on PC
576: Input Parameters:
577: + pc - the preconditioner context
578: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
580: Options Database Key:
581: . -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
583: Level: intermediate
585: .seealso: PCICholeskySetMatOrdering()
586: @*/
587: PetscErrorCode PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
588: {
589: PetscErrorCode ierr,(*f)(PC,MatOrderingType);
592: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)(void))&f);
593: if (f) {
594: (*f)(pc,ordering);
595: }
596: return(0);
597: }
599: /*MC
600: PCCholesky - Uses a direct solver, based on Cholesky factorization, as a preconditioner
602: Options Database Keys:
603: + -pc_cholesky_reuse_ordering - Activate PCLUSetReuseOrdering()
604: . -pc_cholesky_reuse_fill - Activates PCLUSetReuseFill()
605: . -pc_cholesky_fill <fill> - Sets fill amount
606: . -pc_cholesky_damping <damping> - Sets damping amount
607: . -pc_cholesky_shift - Activates Manteuffel shift
608: . -pc_cholesky_in_place - Activates in-place factorization
609: - -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
611: Notes: Not all options work for all matrix formats
613: Level: beginner
615: Concepts: Cholesky factorization, direct solver
617: Notes: Usually this will compute an "exact" solution in one iteration and does
618: not need a Krylov method (i.e. you can use -ksp_type preonly, or
619: KSPSetType(ksp,KSPPREONLY) for the Krylov method
621: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
622: PCILU, PCLU, PCICC, PCCholeskySetReuseOrdering(), PCCholeskySetReuseFill(), PCGetFactoredMatrix(),
623: PCCholeskySetFill(), PCCholeskySetDamping(), PCCholeskySetShift(),
624: PCCholeskySetUseInPlace(), PCCholeskySetMatOrdering()
626: M*/
631: PetscErrorCode PCCreate_Cholesky(PC pc)
632: {
634: PC_Cholesky *dir;
637: PetscNew(PC_Cholesky,&dir);
638: PetscLogObjectMemory(pc,sizeof(PC_Cholesky));
640: dir->fact = 0;
641: dir->inplace = PETSC_FALSE;
642: MatFactorInfoInitialize(&dir->info);
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: }