xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 
2 #include <../src/ksp/pc/impls/factor/icc/icc.h>   /*I "petscpc.h" I*/
3 
4 static PetscErrorCode PCSetUp_ICC(PC pc)
5 {
6   PC_ICC                 *icc = (PC_ICC*)pc->data;
7   IS                     perm = NULL,cperm = NULL;
8   MatInfo                info;
9   MatSolverType          stype;
10   MatFactorError         err;
11   PetscBool              canuseordering;
12 
13   PetscFunctionBegin;
14   pc->failedreason = PC_NOERROR;
15 
16   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
17   if (!pc->setupcalled) {
18     if (!((PC_Factor*)icc)->fact) {
19       PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact));
20     }
21     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering));
22     if (canuseordering) {
23       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
24       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm));
25     }
26     PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info));
27   } else if (pc->flag != SAME_NONZERO_PATTERN) {
28     PetscBool canuseordering;
29     PetscCall(MatDestroy(&((PC_Factor*)icc)->fact));
30     PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact));
31     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering));
32     if (canuseordering) {
33       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
34       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm));
35     }
36     PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info));
37   }
38   PetscCall(MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info));
39   icc->hdr.actualfill = info.fill_ratio_needed;
40 
41   PetscCall(ISDestroy(&cperm));
42   PetscCall(ISDestroy(&perm));
43 
44   PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err));
45   if (err) { /* FactorSymbolic() fails */
46     pc->failedreason = (PCFailedReason)err;
47     PetscFunctionReturn(0);
48   }
49 
50   PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info));
51   PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err));
52   if (err) { /* FactorNumeric() fails */
53     pc->failedreason = (PCFailedReason)err;
54   }
55 
56   PetscCall(PCFactorGetMatSolverType(pc,&stype));
57   if (!stype) {
58     MatSolverType solverpackage;
59     PetscCall(MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage));
60     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
61   }
62   PetscFunctionReturn(0);
63 }
64 
65 static PetscErrorCode PCReset_ICC(PC pc)
66 {
67   PC_ICC         *icc = (PC_ICC*)pc->data;
68 
69   PetscFunctionBegin;
70   PetscCall(MatDestroy(&((PC_Factor*)icc)->fact));
71   PetscFunctionReturn(0);
72 }
73 
74 static PetscErrorCode PCDestroy_ICC(PC pc)
75 {
76   PC_ICC         *icc = (PC_ICC*)pc->data;
77 
78   PetscFunctionBegin;
79   PetscCall(PCReset_ICC(pc));
80   PetscCall(PetscFree(((PC_Factor*)icc)->ordering));
81   PetscCall(PetscFree(((PC_Factor*)icc)->solvertype));
82   PetscCall(PetscFree(pc->data));
83   PetscFunctionReturn(0);
84 }
85 
86 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
87 {
88   PC_ICC         *icc = (PC_ICC*)pc->data;
89 
90   PetscFunctionBegin;
91   PetscCall(MatSolve(((PC_Factor*)icc)->fact,x,y));
92   PetscFunctionReturn(0);
93 }
94 
95 static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y)
96 {
97   PC_ICC         *icc = (PC_ICC*)pc->data;
98 
99   PetscFunctionBegin;
100   PetscCall(MatMatSolve(((PC_Factor*)icc)->fact,X,Y));
101   PetscFunctionReturn(0);
102 }
103 
104 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
105 {
106   PC_ICC         *icc = (PC_ICC*)pc->data;
107 
108   PetscFunctionBegin;
109   PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y));
110   PetscFunctionReturn(0);
111 }
112 
113 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
114 {
115   PC_ICC         *icc = (PC_ICC*)pc->data;
116 
117   PetscFunctionBegin;
118   PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y));
119   PetscFunctionReturn(0);
120 }
121 
122 static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
123 {
124   PC_ICC         *icc = (PC_ICC*)pc->data;
125   PetscBool      flg;
126   /* PetscReal      dt[3];*/
127 
128   PetscFunctionBegin;
129   PetscCall(PetscOptionsHead(PetscOptionsObject,"ICC Options"));
130   PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc));
131 
132   PetscCall(PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg));
133   /*dt[0] = ((PC_Factor*)icc)->info.dt;
134   dt[1] = ((PC_Factor*)icc)->info.dtcol;
135   dt[2] = ((PC_Factor*)icc)->info.dtcount;
136   PetscInt       dtmax = 3;
137   PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg));
138   if (flg) {
139     PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]));
140   }
141   */
142   PetscCall(PetscOptionsTail());
143   PetscFunctionReturn(0);
144 }
145 
146 extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
147 
148 /*MC
149      PCICC - Incomplete Cholesky factorization preconditioners.
150 
151    Options Database Keys:
152 +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
153 .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
154                       its factorization (overwrites original matrix)
155 .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
156 -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
157 
158    Level: beginner
159 
160    Notes:
161     Only implemented for some matrix formats. Not implemented in parallel.
162 
163           For BAIJ matrices this implements a point block ICC.
164 
165           The Manteuffel shift is only implemented for matrices with block size 1
166 
167           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
168           to turn off the shift.
169 
170    References:
171 .  * - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
172       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
173       Science and Engineering, Kluwer.
174 
175 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
176            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
177            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
178            PCFactorSetLevels()
179 
180 M*/
181 
182 PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
183 {
184   PC_ICC         *icc;
185 
186   PetscFunctionBegin;
187   PetscCall(PetscNewLog(pc,&icc));
188   pc->data = (void*)icc;
189   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
190 
191   ((PC_Factor*)icc)->info.fill      = 1.0;
192   ((PC_Factor*)icc)->info.dtcol     = PETSC_DEFAULT;
193   ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
194 
195   pc->ops->apply               = PCApply_ICC;
196   pc->ops->matapply            = PCMatApply_ICC;
197   pc->ops->applytranspose      = PCApply_ICC;
198   pc->ops->setup               = PCSetUp_ICC;
199   pc->ops->reset               = PCReset_ICC;
200   pc->ops->destroy             = PCDestroy_ICC;
201   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
202   pc->ops->view                = PCView_Factor;
203   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
204   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
205   PetscFunctionReturn(0);
206 }
207