Actual source code: cholesky.c
1: /*$Id: cholesky.c,v 1.10 2001/04/10 19:36:21 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 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: PetscOptionsDouble("-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: PetscOptionsDouble("-pc_cholesky_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,0,0);
76:
77: PetscOptionsTail();
78: return(0);
79: }
81: static int PCView_Cholesky(PC pc,PetscViewer viewer)
82: {
83: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
84: int ierr;
85: PetscTruth isascii,isstring;
86:
88: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
89: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
90: if (isascii) {
91: MatInfo info;
92:
93: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorizationn");}
94: else {PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorizationn");}
95: PetscViewerASCIIPrintf(viewer," matrix ordering: %sn",lu->ordering);
96: if (lu->fact) {
97: MatGetInfo(lu->fact,MAT_LOCAL,&info);
98: PetscViewerASCIIPrintf(viewer," Cholesky nonzeros %gn",info.nz_used);
99: }
100: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorizationn");}
101: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorizationn");}
102: } else if (isstring) {
103: PetscViewerStringSPrintf(viewer," order=%s",lu->ordering);
104: } else {
105: SETERRQ1(1,"Viewer type %s not supported for PCCholesky",((PetscObject)viewer)->type_name);
106: }
107: return(0);
108: }
110: static int PCGetFactoredMatrix_Cholesky(PC pc,Mat *mat)
111: {
112: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
113:
115: if (!dir->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
116: *mat = dir->fact;
117: return(0);
118: }
120: static int PCSetUp_Cholesky(PC pc)
121: {
122: int ierr;
123: PetscTruth flg;
124: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
127: if (dir->reusefill && pc->setupcalled) dir->info.fill = dir->actualfill;
128:
129: if (dir->inplace) {
130: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
131: if (dir->col) {ISDestroy(dir->col);}
132: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
133: if (dir->col && dir->row != dir->col)
134: {ISDestroy(dir->col); dir->col=0;} /* only use row ordering for SBAIJ */
135: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
136: MatCholeskyFactor(pc->pmat,dir->row,dir->info.fill);
137: dir->fact = pc->pmat;
138: } else {
139: MatInfo info;
140: if (!pc->setupcalled) {
141: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
142: if (dir->col && dir->row != dir->col)
143: {ISDestroy(dir->col); dir->col=0;} /* only use row ordering for SBAIJ */
144: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
145: if (flg) {
146: PetscReal tol = 1.e-10;
147: PetscOptionsGetDouble(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
148: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
149: }
150: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
151: MatCholeskyFactorSymbolic(pc->pmat,dir->row,dir->info.fill,&dir->fact);
152: MatGetInfo(dir->fact,MAT_LOCAL,&info);
153: dir->actualfill = info.fill_ratio_needed;
154: PetscLogObjectParent(pc,dir->fact);
155: } else if (pc->flag != SAME_NONZERO_PATTERN) {
156: if (!dir->reuseordering) {
157: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
158: if (dir->col) {ISDestroy(dir->col);}
159: MatGetOrdering(pc->pmat,dir->ordering,&dir->row,&dir->col);
160: if (dir->col && dir->row != dir->col)
161: {ISDestroy(dir->col); dir->col=0;} /* only use row ordering for SBAIJ */
162: PetscOptionsHasName(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&flg);
163: if (flg) {
164: PetscReal tol = 1.e-10;
165: PetscOptionsGetDouble(pc->prefix,"-pc_cholesky_nonzeros_along_diagonal",&tol,PETSC_NULL);
166: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
167: }
168: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
169: }
170: MatDestroy(dir->fact);
171: MatCholeskyFactorSymbolic(pc->pmat,dir->row,dir->info.fill,&dir->fact);
172: MatGetInfo(dir->fact,MAT_LOCAL,&info);
173: dir->actualfill = info.fill_ratio_needed;
174: PetscLogObjectParent(pc,dir->fact);
175: }
176: MatCholeskyFactorNumeric(pc->pmat,&dir->fact);
177: }
178: return(0);
179: }
181: static int PCDestroy_Cholesky(PC pc)
182: {
183: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
184: int ierr;
187: if (!dir->inplace && dir->fact) {MatDestroy(dir->fact);}
188: if (dir->row && dir->col && dir->row != dir->col) {ISDestroy(dir->row);}
189: if (dir->col) {ISDestroy(dir->col);}
190: PetscStrfree(dir->ordering);
191: PetscFree(dir);
192: return(0);
193: }
195: static int PCApply_Cholesky(PC pc,Vec x,Vec y)
196: {
197: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
198: int ierr;
199:
201: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
202: else {MatSolve(dir->fact,x,y);}
203: return(0);
204: }
206: static int PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
207: {
208: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
209: int ierr;
212: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
213: else {MatSolveTranspose(dir->fact,x,y);}
214: return(0);
215: }
217: /* -----------------------------------------------------------------------------------*/
219: EXTERN_C_BEGIN
220: int PCCholeskySetFill_Cholesky(PC pc,PetscReal fill)
221: {
222: PC_Cholesky *dir;
223:
225: dir = (PC_Cholesky*)pc->data;
226: dir->info.fill = fill;
227: return(0);
228: }
229: EXTERN_C_END
231: EXTERN_C_BEGIN
232: int PCCholeskySetDamping_Cholesky(PC pc,PetscReal damping)
233: {
234: PC_Cholesky *dir;
235:
237: dir = (PC_Cholesky*)pc->data;
238: dir->info.damping = damping;
239: dir->info.damp = 1.0;
240: return(0);
241: }
242: EXTERN_C_END
244: EXTERN_C_BEGIN
245: int PCCholeskySetUseInPlace_Cholesky(PC pc)
246: {
247: PC_Cholesky *dir;
250: dir = (PC_Cholesky*)pc->data;
251: dir->inplace = 1;
252: return(0);
253: }
254: EXTERN_C_END
256: EXTERN_C_BEGIN
257: int PCCholeskySetMatOrdering_Cholesky(PC pc,MatOrderingType ordering)
258: {
259: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
260: int ierr;
263: PetscStrfree(dir->ordering);
264: PetscStrallocpy(ordering,&dir->ordering);
265: return(0);
266: }
267: EXTERN_C_END
269: /* -----------------------------------------------------------------------------------*/
271: /*@
272: PCCholeskySetReuseOrdering - When similar matrices are factored, this
273: causes the ordering computed in the first factor to be used for all
274: following factors.
276: Collective on PC
278: Input Parameters:
279: + pc - the preconditioner context
280: - flag - PETSC_TRUE to reuse else PETSC_FALSE
282: Options Database Key:
283: . -pc_cholesky_reuse_ordering - Activate PCCholeskySetReuseOrdering()
285: Level: intermediate
287: .keywords: PC, levels, reordering, factorization, incomplete, LU
289: .seealso: PCCholeskySetReuseFill(), PCICholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
290: @*/
291: int PCCholeskySetReuseOrdering(PC pc,PetscTruth flag)
292: {
293: int ierr,(*f)(PC,PetscTruth);
297: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseOrdering_C",(void (**)())&f);
298: if (f) {
299: (*f)(pc,flag);
300: }
301: return(0);
302: }
304: /*@
305: PCCholeskySetReuseFill - When matrices with same nonzero structure are Cholesky factored,
306: this causes later ones to use the fill computed in the initial factorization.
308: Collective on PC
310: Input Parameters:
311: + pc - the preconditioner context
312: - flag - PETSC_TRUE to reuse else PETSC_FALSE
314: Options Database Key:
315: . -pc_cholesky_reuse_fill - Activates PCCholeskySetReuseFill()
317: Level: intermediate
319: .keywords: PC, levels, reordering, factorization, incomplete, Cholesky
321: .seealso: PCICholeskySetReuseOrdering(), PCCholeskySetReuseOrdering(), PCICholeskyDTSetReuseFill()
322: @*/
323: int PCCholeskySetReuseFill(PC pc,PetscTruth flag)
324: {
325: int ierr,(*f)(PC,PetscTruth);
329: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetReuseFill_C",(void (**)())&f);
330: if (f) {
331: (*f)(pc,flag);
332: }
333: return(0);
334: }
336: /*@
337: PCCholeskySetFill - Indicates the amount of fill you expect in the factored matrix,
338: fill = number nonzeros in factor/number nonzeros in original matrix.
340: Collective on PC
341:
342: Input Parameters:
343: + pc - the preconditioner context
344: - fill - amount of expected fill
346: Options Database Key:
347: . -pc_cholesky_fill <fill> - Sets fill amount
349: Level: intermediate
351: Note:
352: For sparse matrix factorizations it is difficult to predict how much
353: fill to expect. By running with the option -log_info PETSc will print the
354: actual amount of fill used; allowing you to set the value accurately for
355: future runs. Default PETSc uses a value of 5.0
357: .keywords: PC, set, factorization, direct, fill
359: .seealso: PCILUSetFill()
360: @*/
361: int PCCholeskySetFill(PC pc,PetscReal fill)
362: {
363: int ierr,(*f)(PC,PetscReal);
367: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0");
368: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetFill_C",(void (**)())&f);
369: if (f) {
370: (*f)(pc,fill);
371: }
372: return(0);
373: }
375: /*@
376: PCCholeskySetDamping - Adds this quantity to the diagonal of the matrix during the
377: Cholesky numerical factorization.
379: Collective on PC
380:
381: Input Parameters:
382: + pc - the preconditioner context
383: - damping - amount of damping
385: Options Database Key:
386: . -pc_cholesky_damping <damping> - Sets damping amount
388: Level: intermediate
390: .keywords: PC, set, factorization, direct, fill
392: .seealso: PCICholeskySetFill(), PCILUSetDamp()
393: @*/
394: int PCCholeskySetDamping(PC pc,PetscReal damping)
395: {
396: int ierr,(*f)(PC,PetscReal);
400: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetDamping_C",(void (**)())&f);
401: if (f) {
402: (*f)(pc,damping);
403: }
404: return(0);
405: }
407: /*@
408: PCCholeskySetUseInPlace - Tells the system to do an in-place factorization.
409: For dense matrices, this enables the solution of much larger problems.
410: For sparse matrices the factorization cannot be done truly in-place
411: so this does not save memory during the factorization, but after the matrix
412: is factored, the original unfactored matrix is freed, thus recovering that
413: space.
415: Collective on PC
417: Input Parameters:
418: . pc - the preconditioner context
420: Options Database Key:
421: . -pc_cholesky_in_place - Activates in-place factorization
423: Notes:
424: PCCholeskySetUseInplace() can only be used with the KSP method KSPPREONLY or when
425: a different matrix is provided for the multiply and the preconditioner in
426: a call to SLESSetOperators().
427: This is because the Krylov space methods require an application of the
428: matrix multiplication, which is not possible here because the matrix has
429: been factored in-place, replacing the original matrix.
431: Level: intermediate
433: .keywords: PC, set, factorization, direct, inplace, in-place, Cholesky
435: .seealso: PCICholeskySetUseInPlace()
436: @*/
437: int PCCholeskySetUseInPlace(PC pc)
438: {
439: int ierr,(*f)(PC);
443: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetUseInPlace_C",(void (**)())&f);
444: if (f) {
445: (*f)(pc);
446: }
447: return(0);
448: }
450: /*@
451: PCCholeskySetMatOrdering - Sets the ordering routine (to reduce fill) to
452: be used it the Cholesky factorization.
454: Collective on PC
456: Input Parameters:
457: + pc - the preconditioner context
458: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
460: Options Database Key:
461: . -pc_cholesky_mat_ordering_type <nd,rcm,...> - Sets ordering routine
463: Level: intermediate
465: .seealso: PCICholeskySetMatOrdering()
466: @*/
467: int PCCholeskySetMatOrdering(PC pc,MatOrderingType ordering)
468: {
469: int ierr,(*f)(PC,MatOrderingType);
472: PetscObjectQueryFunction((PetscObject)pc,"PCCholeskySetMatOrdering_C",(void (**)())&f);
473: if (f) {
474: (*f)(pc,ordering);
475: }
476: return(0);
477: }
479: /* ------------------------------------------------------------------------ */
481: EXTERN_C_BEGIN
482: int PCCreate_Cholesky(PC pc)
483: {
484: int ierr;
485: PC_Cholesky *dir;
488: PetscNew(PC_Cholesky,&dir);
489: PetscLogObjectMemory(pc,sizeof(PC_Cholesky));
491: dir->fact = 0;
492: dir->inplace = 0;
493: dir->info.fill = 5.0;
494: dir->info.damping = 0.0;
495: dir->info.damp = 0.0;
496: dir->col = 0;
497: dir->row = 0;
498: PetscStrallocpy(MATORDERING_NATURAL,&dir->ordering);
499: dir->reusefill = PETSC_FALSE;
500: dir->reuseordering = PETSC_FALSE;
501: pc->data = (void*)dir;
503: pc->ops->destroy = PCDestroy_Cholesky;
504: pc->ops->apply = PCApply_Cholesky;
505: pc->ops->applytranspose = PCApplyTranspose_Cholesky;
506: pc->ops->setup = PCSetUp_Cholesky;
507: pc->ops->setfromoptions = PCSetFromOptions_Cholesky;
508: pc->ops->view = PCView_Cholesky;
509: pc->ops->applyrichardson = 0;
510: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_Cholesky;
512: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetFill_C","PCCholeskySetFill_Cholesky",
513: PCCholeskySetFill_Cholesky);
514: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetDamping_C","PCCholeskySetDamping_Cholesky",
515: PCCholeskySetDamping_Cholesky);
516: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetUseInPlace_C","PCCholeskySetUseInPlace_Cholesky",
517: PCCholeskySetUseInPlace_Cholesky);
518: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetMatOrdering_C","PCCholeskySetMatOrdering_Cholesky",
519: PCCholeskySetMatOrdering_Cholesky);
520: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseOrdering_C","PCCholeskySetReuseOrdering_Cholesky",
521: PCCholeskySetReuseOrdering_Cholesky);
522: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCholeskySetReuseFill_C","PCCholeskySetReuseFill_Cholesky",
523: PCCholeskySetReuseFill_Cholesky);
524: return(0);
525: }
526: EXTERN_C_END