Actual source code: icc.c
1: #define PETSCKSP_DLL
3: #include ../src/ksp/pc/impls/factor/icc/icc.h
7: static PetscErrorCode PCSetup_ICC(PC pc)
8: {
9: PC_ICC *icc = (PC_ICC*)pc->data;
10: IS perm,cperm;
12: MatInfo info;
15: MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);
17: if (!pc->setupcalled) {
18: MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);
19: MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);
20: } else if (pc->flag != SAME_NONZERO_PATTERN) {
21: MatDestroy(((PC_Factor*)icc)->fact);
22: MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);
23: MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);
24: }
25: MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);
26: icc->actualfill = info.fill_ratio_needed;
28: ISDestroy(cperm);
29: ISDestroy(perm);
30: MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);
31: return(0);
32: }
36: static PetscErrorCode PCDestroy_ICC(PC pc)
37: {
38: PC_ICC *icc = (PC_ICC*)pc->data;
42: if (((PC_Factor*)icc)->fact) {MatDestroy(((PC_Factor*)icc)->fact);}
43: PetscStrfree(((PC_Factor*)icc)->ordering);
44: PetscStrfree(((PC_Factor*)icc)->solvertype);
45: PetscFree(icc);
46: return(0);
47: }
51: static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
52: {
53: PC_ICC *icc = (PC_ICC*)pc->data;
57: MatSolve(((PC_Factor*)icc)->fact,x,y);
58: return(0);
59: }
63: static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
64: {
66: PC_ICC *icc = (PC_ICC*)pc->data;
69: MatForwardSolve(((PC_Factor*)icc)->fact,x,y);
70: return(0);
71: }
75: static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
76: {
78: PC_ICC *icc = (PC_ICC*)pc->data;
81: MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);
82: return(0);
83: }
87: static PetscErrorCode PCSetFromOptions_ICC(PC pc)
88: {
89: PC_ICC *icc = (PC_ICC*)pc->data;
90: PetscTruth flg;
91: PetscReal dt[3];
95: PetscOptionsHead("ICC Options");
96: PCSetFromOptions_Factor(pc);
98: PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);
99: dt[0] = ((PC_Factor*)icc)->info.dt;
100: dt[1] = ((PC_Factor*)icc)->info.dtcol;
101: dt[2] = ((PC_Factor*)icc)->info.dtcount;
102: /*
103: PetscInt dtmax = 3;
104: PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);
105: if (flg) {
106: PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);
107: }
108: */
109: PetscOptionsTail();
110: return(0);
111: }
115: static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
116: {
118:
120: PCView_Factor(pc,viewer);
121: return(0);
122: }
128: /*MC
129: PCICC - Incomplete Cholesky factorization preconditioners.
131: Options Database Keys:
132: + -pc_factor_levels <k> - number of levels of fill for ICC(k)
133: . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
134: its factorization (overwrites original matrix)
135: . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
136: - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
138: Level: beginner
140: Concepts: incomplete Cholesky factorization
142: Notes: Only implemented for some matrix formats. Not implemented in parallel.
144: For BAIJ matrices this implements a point block ICC.
146: The Manteuffel shift is only implemented for matrices with block size 1
148: By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
149: to turn off the shift.
151: References:
152: Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
153: http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
154: Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
155: Science and Engineering, Kluwer, pp. 167--202.
158: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
159: PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
160: PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
161: PCFactorSetLevels()
163: M*/
168: PetscErrorCode PCCreate_ICC(PC pc)
169: {
171: PC_ICC *icc;
174: PetscNewLog(pc,PC_ICC,&icc);
176: ((PC_Factor*)icc)->fact = 0;
177: PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);
178: PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)icc)->solvertype);
179: MatFactorInfoInitialize(&((PC_Factor*)icc)->info);
180: ((PC_Factor*)icc)->factortype = MAT_FACTOR_ICC;
181: ((PC_Factor*)icc)->info.levels = 0.;
182: ((PC_Factor*)icc)->info.fill = 1.0;
183: icc->implctx = 0;
185: ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT;
186: ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
187: ((PC_Factor*)icc)->info.shiftamount = 0.0;
188: ((PC_Factor*)icc)->info.zeropivot = 1.e-12;
190: pc->data = (void*)icc;
191: pc->ops->apply = PCApply_ICC;
192: pc->ops->setup = PCSetup_ICC;
193: pc->ops->destroy = PCDestroy_ICC;
194: pc->ops->setfromoptions = PCSetFromOptions_ICC;
195: pc->ops->view = PCView_ICC;
196: pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor;
197: pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC;
198: pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
200: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
201: PCFactorGetMatSolverPackage_Factor);
202: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
203: PCFactorSetZeroPivot_Factor);
204: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
205: PCFactorSetShiftType_Factor);
206: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
207: PCFactorSetShiftAmount_Factor);
208: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
209: PCFactorSetLevels_Factor);
210: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
211: PCFactorSetFill_Factor);
212: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
213: PCFactorSetMatOrderingType_Factor);
214: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
215: PCFactorSetMatSolverPackage_Factor);
216: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
217: PCFactorSetDropTolerance_ILU);
218: return(0);
219: }