xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 3ff96bc97e95df46fbdbdf211ddbac0ba3669e66)
1 #define PETSCKSP_DLL
2 
3 #include "../src/ksp/pc/impls/factor/icc/icc.h"   /*I "petscpc.h" I*/
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "PCSetup_ICC"
7 static PetscErrorCode PCSetup_ICC(PC pc)
8 {
9   PC_ICC         *icc = (PC_ICC*)pc->data;
10   IS             perm,cperm;
11   PetscErrorCode ierr;
12   MatInfo        info;
13 
14   PetscFunctionBegin;
15   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
16 
17   if (!pc->setupcalled) {
18     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr);
19     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
20   } else if (pc->flag != SAME_NONZERO_PATTERN) {
21     ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);
22     ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
23     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
24   }
25   ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
26   icc->actualfill = info.fill_ratio_needed;
27 
28   ierr = ISDestroy(cperm);CHKERRQ(ierr);
29   ierr = ISDestroy(perm);CHKERRQ(ierr);
30   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 #undef __FUNCT__
35 #define __FUNCT__ "PCDestroy_ICC"
36 static PetscErrorCode PCDestroy_ICC(PC pc)
37 {
38   PC_ICC         *icc = (PC_ICC*)pc->data;
39   PetscErrorCode ierr;
40 
41   PetscFunctionBegin;
42   if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);}
43   ierr = PetscStrfree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
44   ierr = PetscStrfree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
45   ierr = PetscFree(icc);CHKERRQ(ierr);
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "PCApply_ICC"
51 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
52 {
53   PC_ICC         *icc = (PC_ICC*)pc->data;
54   PetscErrorCode ierr;
55 
56   PetscFunctionBegin;
57   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "PCApplySymmetricLeft_ICC"
63 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
64 {
65   PetscErrorCode ierr;
66   PC_ICC         *icc = (PC_ICC*)pc->data;
67 
68   PetscFunctionBegin;
69   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "PCApplySymmetricRight_ICC"
75 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
76 {
77   PetscErrorCode ierr;
78   PC_ICC         *icc = (PC_ICC*)pc->data;
79 
80   PetscFunctionBegin;
81   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "PCSetFromOptions_ICC"
87 static PetscErrorCode PCSetFromOptions_ICC(PC pc)
88 {
89   PC_ICC         *icc = (PC_ICC*)pc->data;
90   PetscTruth     flg;
91   PetscReal      dt[3];
92   PetscErrorCode ierr;
93 
94   PetscFunctionBegin;
95   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
96     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
97 
98     ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
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     ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
105     if (flg) {
106       ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
107     }
108     */
109   ierr = PetscOptionsTail();CHKERRQ(ierr);
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNCT__
114 #define __FUNCT__ "PCView_ICC"
115 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
116 {
117   PetscErrorCode ierr;
118 
119   PetscFunctionBegin;
120   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 EXTERN_C_BEGIN
125 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
126 EXTERN_C_END
127 
128 /*MC
129      PCICC - Incomplete Cholesky factorization preconditioners.
130 
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
137 
138    Level: beginner
139 
140   Concepts: incomplete Cholesky factorization
141 
142    Notes: Only implemented for some matrix formats. Not implemented in parallel.
143 
144           For BAIJ matrices this implements a point block ICC.
145 
146           The Manteuffel shift is only implemented for matrices with block size 1
147 
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.
150 
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.
156 
157 
158 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
159            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
160            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
161            PCFactorSetLevels()
162 
163 M*/
164 
165 EXTERN_C_BEGIN
166 #undef __FUNCT__
167 #define __FUNCT__ "PCCreate_ICC"
168 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc)
169 {
170   PetscErrorCode ierr;
171   PC_ICC         *icc;
172 
173   PetscFunctionBegin;
174   ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr);
175 
176   ((PC_Factor*)icc)->fact	          = 0;
177   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
178   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
179   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
180   ((PC_Factor*)icc)->info.levels	  = 0.;
181   ((PC_Factor*)icc)->info.fill          = 1.0;
182   icc->implctx            = 0;
183 
184   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
185   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
186   ((PC_Factor*)icc)->info.shiftamount = 0.0;
187   ((PC_Factor*)icc)->info.zeropivot   = 1.e-12;
188 
189   pc->data	               = (void*)icc;
190   pc->ops->apply	       = PCApply_ICC;
191   pc->ops->setup               = PCSetup_ICC;
192   pc->ops->destroy	       = PCDestroy_ICC;
193   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
194   pc->ops->view                = PCView_ICC;
195   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
196   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
197   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
198 
199   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
200                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
201   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
202                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
203   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
204                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
205   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
206                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
207   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
208                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
209   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
210                     PCFactorSetFill_Factor);CHKERRQ(ierr);
211   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
212                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
213   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
214                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
215   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
216                     PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 EXTERN_C_END
220 
221 
222