xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 657ed8aec61a800d97f0df7c487e5eaa2941ceea)
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   PetscInt       dtmax = 3;
92   PetscReal      dt[3];
93   PetscErrorCode ierr;
94 
95   PetscFunctionBegin;
96   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
97     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
98 
99     ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
100     dt[0] = ((PC_Factor*)icc)->info.dt;
101     dt[1] = ((PC_Factor*)icc)->info.dtcol;
102     dt[2] = ((PC_Factor*)icc)->info.dtcount;
103     ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
104     if (flg) {
105       ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
106     }
107 
108   ierr = PetscOptionsTail();CHKERRQ(ierr);
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "PCView_ICC"
114 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
115 {
116   PC_ICC         *icc = (PC_ICC*)pc->data;
117   PetscErrorCode ierr;
118   PetscTruth     isstring,iascii;
119 
120   PetscFunctionBegin;
121   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
122   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
123   if (iascii) {
124     if (((PC_Factor*)icc)->info.levels == 1) {
125         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D level of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
126     } else {
127         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D levels of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
128     }
129     ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio allocated %G, ordering used %s\n",((PC_Factor*)icc)->info.fill,((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
130     if (((PC_Factor*)icc)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using Manteuffel shift\n");CHKERRQ(ierr);}
131     if (((PC_Factor*)icc)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
132     if (((PC_Factor*)icc)->fact) {
133       ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr);
134       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
135       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
136       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
137       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
138       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
139       ierr = MatView(((PC_Factor*)icc)->fact,viewer);CHKERRQ(ierr);
140       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
141       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
142       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
143       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
144     }
145   } else if (isstring) {
146     ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);CHKERRQ(ierr);
147   } else {
148     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name);
149   }
150   PetscFunctionReturn(0);
151 }
152 
153 extern "C" PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
154 
155 /*MC
156      PCICC - Incomplete Cholesky factorization preconditioners.
157 
158    Options Database Keys:
159 +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
160 .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
161                       its factorization (overwrites original matrix)
162 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
163 .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
164 .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
165 -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
166    is optional with PETSC_TRUE being the default
167 
168    Level: beginner
169 
170   Concepts: incomplete Cholesky factorization
171 
172    Notes: Only implemented for some matrix formats. Not implemented in parallel.
173 
174           For BAIJ matrices this implements a point block ICC.
175 
176           The Manteuffel shift is only implemented for matrices with block size 1
177 
178           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE);
179           to turn off the shift.
180 
181    References:
182    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
183       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
184       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
185       Science and Engineering, Kluwer, pp. 167--202.
186 
187 
188 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
189            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),PCFactorSetShiftInBlocks(),
190            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
191            PCFactorSetLevels()
192 
193 M*/
194 
195 EXTERN_C_BEGIN
196 #undef __FUNCT__
197 #define __FUNCT__ "PCCreate_ICC"
198 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc)
199 {
200   PetscErrorCode ierr;
201   PC_ICC         *icc;
202 
203   PetscFunctionBegin;
204   ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr);
205 
206   ((PC_Factor*)icc)->fact	          = 0;
207   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
208   ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
209   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
210   ((PC_Factor*)icc)->info.levels	  = 0.;
211   ((PC_Factor*)icc)->info.fill          = 1.0;
212   icc->implctx            = 0;
213 
214   ((PC_Factor*)icc)->info.dtcol              = PETSC_DEFAULT;
215   ((PC_Factor*)icc)->info.shiftnz            = 0.0;
216   ((PC_Factor*)icc)->info.shiftpd            = 1.0; /* true */
217   ((PC_Factor*)icc)->info.zeropivot          = 1.e-12;
218   pc->data	               = (void*)icc;
219 
220   pc->ops->apply	       = PCApply_ICC;
221   pc->ops->setup               = PCSetup_ICC;
222   pc->ops->destroy	       = PCDestroy_ICC;
223   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
224   pc->ops->view                = PCView_ICC;
225   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
226   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
227   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
228 
229   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
230                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
231   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
232                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
233   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
234                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
235   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
236                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
237   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor",
238                     PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr);
239 
240   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
241                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
242   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
243                     PCFactorSetFill_Factor);CHKERRQ(ierr);
244   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
245                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
246   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
247                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
248   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU",
249                     PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 EXTERN_C_END
253 
254 
255