Actual source code: cholesky.c
1: #define PETSCKSP_DLL
3: /*
4: Defines a direct factorization preconditioner for any Mat implementation
5: Note: this need not be consided a preconditioner since it supplies
6: a direct solver.
7: */
8: #include ../src/ksp/pc/impls/factor/factor.h
10: typedef struct {
11: PC_Factor hdr;
12: PetscReal actualfill; /* actual fill in factor */
13: PetscTruth inplace; /* flag indicating in-place factorization */
14: IS row,col; /* index sets used for reordering */
15: PetscTruth reuseordering; /* reuses previous reordering computed */
16: PetscTruth reusefill; /* reuse fill from previous Cholesky */
17: } PC_Cholesky;
22: PetscErrorCode PCFactorSetReuseOrdering_Cholesky(PC pc,PetscTruth flag)
23: {
24: PC_Cholesky *lu;
25:
27: lu = (PC_Cholesky*)pc->data;
28: lu->reuseordering = flag;
29: return(0);
30: }
36: PetscErrorCode PCFactorSetReuseFill_Cholesky(PC pc,PetscTruth flag)
37: {
38: PC_Cholesky *lu;
39:
41: lu = (PC_Cholesky*)pc->data;
42: lu->reusefill = flag;
43: return(0);
44: }
49: static PetscErrorCode PCSetFromOptions_Cholesky(PC pc)
50: {
51: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
53: PetscTruth flg;
54: char tname[256], solvertype[64];
55: PetscFList ordlist;
56:
58: if (!MatOrderingRegisterAllCalled) {MatOrderingRegisterAll(PETSC_NULL);}
59: PetscOptionsHead("Cholesky options");
60: PetscOptionsName("-pc_factor_in_place","Form Cholesky in the same memory as the matrix","PCFactorSetUseInPlace",&flg);
61: if (flg) {
62: PCFactorSetUseInPlace(pc);
63: }
64: PetscOptionsReal("-pc_factor_fill","Expected non-zeros in Cholesky/non-zeros in matrix","PCFactorSetFill",((PC_Factor*)lu)->info.fill,&((PC_Factor*)lu)->info.fill,0);
65:
66: PetscOptionsName("-pc_factor_reuse_fill","Use fill from previous factorization","PCFactorSetReuseFill",&flg);
67: if (flg) {
68: PCFactorSetReuseFill(pc,PETSC_TRUE);
69: }
70: PetscOptionsName("-pc_factor_reuse_ordering","Reuse ordering from previous factorization","PCFactorSetReuseOrdering",&flg);
71: if (flg) {
72: PCFactorSetReuseOrdering(pc,PETSC_TRUE);
73: }
74:
75: MatGetOrderingList(&ordlist);
76: PetscOptionsList("-pc_factor_mat_ordering_type","Reordering to reduce nonzeros in Cholesky","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)lu)->ordering,tname,256,&flg);
77: if (flg) {
78: PCFactorSetMatOrderingType(pc,tname);
79: }
81: /* maybe should have MatGetSolverTypes(Mat,&list) like the ordering list */
82: PetscOptionsString("-pc_factor_mat_solver_package","Specific Cholesky solver to use","MatGetFactor",((PC_Factor*)lu)->solvertype,solvertype,64,&flg);
83: if (flg) {
84: PCFactorSetMatSolverPackage(pc,solvertype);
85: }
87: PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);
88: if (flg) {
89: PCFactorSetShiftNonzero(pc,(PetscReal) PETSC_DECIDE);
90: }
91: PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)lu)->info.shiftnz,&((PC_Factor*)lu)->info.shiftnz,0);
92: PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);
93: if (flg) {
94: PCFactorSetShiftPd(pc,PETSC_TRUE);
95: }
96: PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)lu)->info.zeropivot,&((PC_Factor*)lu)->info.zeropivot,0);
98: PetscOptionsTail();
99: return(0);
100: }
104: static PetscErrorCode PCView_Cholesky(PC pc,PetscViewer viewer)
105: {
106: PC_Cholesky *lu = (PC_Cholesky*)pc->data;
108: PetscTruth iascii,isstring;
109:
111: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
112: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
113: if (iascii) {
114:
115: if (lu->inplace) {PetscViewerASCIIPrintf(viewer," Cholesky: in-place factorization\n");}
116: else {PetscViewerASCIIPrintf(viewer," Cholesky: out-of-place factorization\n");}
117: PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);
118: if (((PC_Factor*)lu)->fact) {
119: PetscViewerASCIIPrintf(viewer," Cholesky: factor fill ratio needed %G\n",lu->actualfill);
120: PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");
121: PetscViewerASCIIPushTab(viewer);
122: PetscViewerASCIIPushTab(viewer);
123: PetscViewerASCIIPushTab(viewer);
124: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
125: MatView(((PC_Factor*)lu)->fact,viewer);
126: PetscViewerPopFormat(viewer);
127: PetscViewerASCIIPopTab(viewer);
128: PetscViewerASCIIPopTab(viewer);
129: PetscViewerASCIIPopTab(viewer);
130: }
131: if (lu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");}
132: if (lu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");}
133: } else if (isstring) {
134: PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);
135: } else {
136: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCCHOLESKY",((PetscObject)viewer)->type_name);
137: }
138: return(0);
139: }
144: static PetscErrorCode PCSetUp_Cholesky(PC pc)
145: {
147: PetscTruth flg;
148: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
151: if (dir->reusefill && pc->setupcalled) ((PC_Factor*)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,((PC_Factor*)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,&((PC_Factor*)dir)->info);
169: ((PC_Factor*)dir)->fact = pc->pmat;
170: } else {
171: MatInfo info;
172: if (!pc->setupcalled) {
173: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
174: /* check if dir->row == dir->col */
175: ISEqual(dir->row,dir->col,&flg);
176: if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"row and column permutations must equal");
177: ISDestroy(dir->col); /* only pass one ordering into CholeskyFactor */
178: dir->col=0;
180: PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);
181: if (flg) {
182: PetscReal tol = 1.e-10;
183: PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
184: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
185: }
186: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
187: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
188: MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
189: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
190: dir->actualfill = info.fill_ratio_needed;
191: PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
192: } else if (pc->flag != SAME_NONZERO_PATTERN) {
193: if (!dir->reuseordering) {
194: if (dir->row && dir->col && (dir->row != dir->col)) {
195: ISDestroy(dir->row);
196: dir->row = 0;
197: }
198: if (dir->col) {
199: ISDestroy(dir->col);
200: dir->col =0;
201: }
202: MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);
203: if (dir->col && (dir->row != dir->col)) { /* only use row ordering for SBAIJ */
204: ISDestroy(dir->col);
205: dir->col=0;
206: }
207: PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&flg);
208: if (flg) {
209: PetscReal tol = 1.e-10;
210: PetscOptionsGetReal(((PetscObject)pc)->prefix,"-pc_factor_nonzeros_along_diagonal",&tol,PETSC_NULL);
211: MatReorderForNonzeroDiagonal(pc->pmat,tol,dir->row,dir->row);
212: }
213: if (dir->row) {PetscLogObjectParent(pc,dir->row);}
214: }
215: MatDestroy(((PC_Factor*)dir)->fact);
216: MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_CHOLESKY,&((PC_Factor*)dir)->fact);
217: MatCholeskyFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,&((PC_Factor*)dir)->info);
218: MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);
219: dir->actualfill = info.fill_ratio_needed;
220: PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);
221: }
222: MatCholeskyFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);
223: }
224: return(0);
225: }
229: static PetscErrorCode PCDestroy_Cholesky(PC pc)
230: {
231: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
235: if (!dir->inplace && ((PC_Factor*)dir)->fact) {MatDestroy(((PC_Factor*)dir)->fact);}
236: if (dir->row) {ISDestroy(dir->row);}
237: if (dir->col) {ISDestroy(dir->col);}
238: PetscStrfree(((PC_Factor*)dir)->ordering);
239: PetscStrfree(((PC_Factor*)dir)->solvertype);
240: PetscFree(dir);
241: return(0);
242: }
246: static PetscErrorCode PCApply_Cholesky(PC pc,Vec x,Vec y)
247: {
248: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
250:
252: if (dir->inplace) {MatSolve(pc->pmat,x,y);}
253: else {MatSolve(((PC_Factor*)dir)->fact,x,y);}
254: return(0);
255: }
259: static PetscErrorCode PCApplyTranspose_Cholesky(PC pc,Vec x,Vec y)
260: {
261: PC_Cholesky *dir = (PC_Cholesky*)pc->data;
265: if (dir->inplace) {MatSolveTranspose(pc->pmat,x,y);}
266: else {MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);}
267: return(0);
268: }
270: /* -----------------------------------------------------------------------------------*/
275: PetscErrorCode PCFactorSetUseInPlace_Cholesky(PC pc)
276: {
277: PC_Cholesky *dir;
280: dir = (PC_Cholesky*)pc->data;
281: dir->inplace = PETSC_TRUE;
282: return(0);
283: }
286: /* -----------------------------------------------------------------------------------*/
290: /*@
291: PCFactorSetReuseOrdering - When similar matrices are factored, this
292: causes the ordering computed in the first factor to be used for all
293: following factors.
295: Collective on PC
297: Input Parameters:
298: + pc - the preconditioner context
299: - flag - PETSC_TRUE to reuse else PETSC_FALSE
301: Options Database Key:
302: . -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
304: Level: intermediate
306: .keywords: PC, levels, reordering, factorization, incomplete, LU
308: .seealso: PCFactorSetReuseFill()
309: @*/
310: PetscErrorCode PCFactorSetReuseOrdering(PC pc,PetscTruth flag)
311: {
312: PetscErrorCode ierr,(*f)(PC,PetscTruth);
316: PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",(void (**)(void))&f);
317: if (f) {
318: (*f)(pc,flag);
319: }
320: return(0);
321: }
324: /*MC
325: PCCHOLESKY - Uses a direct solver, based on Cholesky factorization, as a preconditioner
327: Options Database Keys:
328: + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering()
329: . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles
330: . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill()
331: . -pc_factor_fill <fill> - Sets fill amount
332: . -pc_factor_in_place - Activates in-place factorization
333: . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine
334: . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
335: - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
336: is optional with PETSC_TRUE being the default
338: Notes: Not all options work for all matrix formats
340: Level: beginner
342: Concepts: Cholesky factorization, direct solver
344: Notes: Usually this will compute an "exact" solution in one iteration and does
345: not need a Krylov method (i.e. you can use -ksp_type preonly, or
346: KSPSetType(ksp,KSPPREONLY) for the Krylov method
348: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC,
349: PCILU, PCLU, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(),
350: PCFactorSetFill(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
351: PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd()
353: M*/
358: PetscErrorCode PCCreate_Cholesky(PC pc)
359: {
361: PC_Cholesky *dir;
364: PetscNewLog(pc,PC_Cholesky,&dir);
366: ((PC_Factor*)dir)->fact = 0;
367: dir->inplace = PETSC_FALSE;
368: MatFactorInfoInitialize(&((PC_Factor*)dir)->info);
369: ((PC_Factor*)dir)->info.fill = 5.0;
370: ((PC_Factor*)dir)->info.shiftnz = 0.0;
371: ((PC_Factor*)dir)->info.shiftpd = 0.0; /* false */
372: ((PC_Factor*)dir)->info.pivotinblocks = 1.0;
373: dir->col = 0;
374: dir->row = 0;
375: PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);
376: PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);
377: dir->reusefill = PETSC_FALSE;
378: dir->reuseordering = PETSC_FALSE;
379: pc->data = (void*)dir;
381: pc->ops->destroy = PCDestroy_Cholesky;
382: pc->ops->apply = PCApply_Cholesky;
383: pc->ops->applytranspose = PCApplyTranspose_Cholesky;
384: pc->ops->setup = PCSetUp_Cholesky;
385: pc->ops->setfromoptions = PCSetFromOptions_Cholesky;
386: pc->ops->view = PCView_Cholesky;
387: pc->ops->applyrichardson = 0;
388: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
390: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
391: PCFactorSetMatSolverPackage_Factor);
392: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
393: PCFactorGetMatSolverPackage_Factor);
394: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
395: PCFactorSetZeroPivot_Factor);
396: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
397: PCFactorSetShiftNonzero_Factor);
398: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
399: PCFactorSetShiftPd_Factor);
401: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
402: PCFactorSetFill_Factor);
403: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_Cholesky",
404: PCFactorSetUseInPlace_Cholesky);
405: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
406: PCFactorSetMatOrderingType_Factor);
407: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_Cholesky",
408: PCFactorSetReuseOrdering_Cholesky);
409: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_Cholesky",
410: PCFactorSetReuseFill_Cholesky);
411: return(0);
412: }