xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 307d0106c7ffa1e3caba7e700e4bcbb71ca7fbbe)
1 
2 #include <../src/ksp/pc/impls/factor/icc/icc.h>   /*I "petscpc.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PCSetUp_ICC"
6 static PetscErrorCode PCSetUp_ICC(PC pc)
7 {
8   PC_ICC         *icc = (PC_ICC*)pc->data;
9   IS             perm,cperm;
10   PetscErrorCode ierr;
11   MatInfo        info;
12   Mat            F;
13   const MatSolverPackage stype;
14 
15   PetscFunctionBegin;
16   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
17 
18   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
19   if (!pc->setupcalled) {
20     if (!((PC_Factor*)icc)->fact) {
21       ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
22     }
23     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
24   } else if (pc->flag != SAME_NONZERO_PATTERN) {
25     ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
26     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
27     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
28   }
29   ierr            = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
30   icc->actualfill = info.fill_ratio_needed;
31 
32   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
33   ierr = ISDestroy(&perm);CHKERRQ(ierr);
34 
35   F = ((PC_Factor*)icc)->fact;
36   if (F->errortype) { /* FactorSymbolic() fails */
37     pc->failedreason = (PCFailedReason)F->errortype;
38     PetscFunctionReturn(0);
39   }
40 
41   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
42   if (F->errortype) { /* FactorNumeric() fails */
43     pc->failedreason = (PCFailedReason)F->errortype;
44   }
45 
46   ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
47   if (!stype) {
48     ierr = PCFactorSetMatSolverPackage(pc,((PC_Factor*)icc)->fact->solvertype);CHKERRQ(ierr);
49   }
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "PCReset_ICC"
55 static PetscErrorCode PCReset_ICC(PC pc)
56 {
57   PC_ICC         *icc = (PC_ICC*)pc->data;
58   PetscErrorCode ierr;
59 
60   PetscFunctionBegin;
61   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "PCDestroy_ICC"
67 static PetscErrorCode PCDestroy_ICC(PC pc)
68 {
69   PC_ICC         *icc = (PC_ICC*)pc->data;
70   PetscErrorCode ierr;
71 
72   PetscFunctionBegin;
73   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
74   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
75   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
76   ierr = PetscFree(pc->data);CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "PCApply_ICC"
82 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
83 {
84   PC_ICC         *icc = (PC_ICC*)pc->data;
85   PetscErrorCode ierr;
86 
87   PetscFunctionBegin;
88   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "PCApplySymmetricLeft_ICC"
94 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
95 {
96   PetscErrorCode ierr;
97   PC_ICC         *icc = (PC_ICC*)pc->data;
98 
99   PetscFunctionBegin;
100   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "PCApplySymmetricRight_ICC"
106 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
107 {
108   PetscErrorCode ierr;
109   PC_ICC         *icc = (PC_ICC*)pc->data;
110 
111   PetscFunctionBegin;
112   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 #undef __FUNCT__
117 #define __FUNCT__ "PCSetFromOptions_ICC"
118 static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
119 {
120   PC_ICC         *icc = (PC_ICC*)pc->data;
121   PetscBool      flg;
122   PetscErrorCode ierr;
123   /* PetscReal      dt[3];*/
124 
125   PetscFunctionBegin;
126   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
127   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
128 
129   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
130   /*dt[0] = ((PC_Factor*)icc)->info.dt;
131   dt[1] = ((PC_Factor*)icc)->info.dtcol;
132   dt[2] = ((PC_Factor*)icc)->info.dtcount;
133   PetscInt       dtmax = 3;
134   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
135   if (flg) {
136     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
137   }
138   */
139   ierr = PetscOptionsTail();CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "PCView_ICC"
145 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
146 {
147   PetscErrorCode ierr;
148 
149   PetscFunctionBegin;
150   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
151   PetscFunctionReturn(0);
152 }
153 
154 extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "PCFactorGetUseInPlace_ICC"
158 PetscErrorCode  PCFactorGetUseInPlace_ICC(PC pc,PetscBool *flg)
159 {
160   PetscFunctionBegin;
161   *flg = PETSC_FALSE;
162   PetscFunctionReturn(0);
163 }
164 
165 /*MC
166      PCICC - Incomplete Cholesky factorization preconditioners.
167 
168    Options Database Keys:
169 +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
170 .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
171                       its factorization (overwrites original matrix)
172 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
173 -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
174 
175    Level: beginner
176 
177   Concepts: incomplete Cholesky factorization
178 
179    Notes: Only implemented for some matrix formats. Not implemented in parallel.
180 
181           For BAIJ matrices this implements a point block ICC.
182 
183           The Manteuffel shift is only implemented for matrices with block size 1
184 
185           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
186           to turn off the shift.
187 
188    References:
189 .  1. - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
190       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
191       Science and Engineering, Kluwer.
192 
193 
194 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
195            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
196            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
197            PCFactorSetLevels()
198 
199 M*/
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PCCreate_ICC"
203 PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
204 {
205   PetscErrorCode ierr;
206   PC_ICC         *icc;
207 
208   PetscFunctionBegin;
209   ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr);
210 
211   ((PC_Factor*)icc)->fact = 0;
212 
213   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
214   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
215 
216   ((PC_Factor*)icc)->factortype  = MAT_FACTOR_ICC;
217   ((PC_Factor*)icc)->info.levels = 0.;
218   ((PC_Factor*)icc)->info.fill   = 1.0;
219   icc->implctx                   = 0;
220 
221   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
222   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
223   ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
224   ((PC_Factor*)icc)->info.zeropivot   = 100.0*PETSC_MACHINE_EPSILON;
225 
226   pc->data                     = (void*)icc;
227   pc->ops->apply               = PCApply_ICC;
228   pc->ops->applytranspose      = PCApply_ICC;
229   pc->ops->setup               = PCSetUp_ICC;
230   pc->ops->reset               = PCReset_ICC;
231   pc->ops->destroy             = PCDestroy_ICC;
232   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
233   pc->ops->view                = PCView_ICC;
234   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
235   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
236   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
237 
238   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
239   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor);CHKERRQ(ierr);
240   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
241   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor);CHKERRQ(ierr);
242   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
243   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor);CHKERRQ(ierr);
244   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
245   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
246   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
247   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
248   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
249   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
250   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
251   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
252   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ICC);CHKERRQ(ierr);
253   PetscFunctionReturn(0);
254 }
255 
256 
257