xref: /petsc/src/ksp/pc/impls/factor/ilu/ilu.c (revision 9b54502bff450e1352e581e897e40ddb4a064c7d)
1*9b54502bSHong Zhang /*
2*9b54502bSHong Zhang    Defines a ILU factorization preconditioner for any Mat implementation
3*9b54502bSHong Zhang */
4*9b54502bSHong Zhang #include "src/ksp/pc/pcimpl.h"                 /*I "petscpc.h"  I*/
5*9b54502bSHong Zhang #include "src/ksp/pc/impls/ilu/ilu.h"
6*9b54502bSHong Zhang #include "src/mat/matimpl.h"
7*9b54502bSHong Zhang 
8*9b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/
9*9b54502bSHong Zhang 
10*9b54502bSHong Zhang EXTERN_C_BEGIN
11*9b54502bSHong Zhang #undef __FUNCT__
12*9b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal_ILU"
13*9b54502bSHong Zhang PetscErrorCode PCILUReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z)
14*9b54502bSHong Zhang {
15*9b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU*)pc->data;
16*9b54502bSHong Zhang 
17*9b54502bSHong Zhang   PetscFunctionBegin;
18*9b54502bSHong Zhang   ilu->nonzerosalongdiagonal = PETSC_TRUE;
19*9b54502bSHong Zhang   if (z == PETSC_DECIDE) {
20*9b54502bSHong Zhang     ilu->nonzerosalongdiagonaltol = 1.e-10;
21*9b54502bSHong Zhang   } else {
22*9b54502bSHong Zhang     ilu->nonzerosalongdiagonaltol = z;
23*9b54502bSHong Zhang   }
24*9b54502bSHong Zhang   PetscFunctionReturn(0);
25*9b54502bSHong Zhang }
26*9b54502bSHong Zhang EXTERN_C_END
27*9b54502bSHong Zhang 
28*9b54502bSHong Zhang EXTERN_C_BEGIN
29*9b54502bSHong Zhang #undef __FUNCT__
30*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetSetZeroPivot_ILU"
31*9b54502bSHong Zhang PetscErrorCode PCILUSetZeroPivot_ILU(PC pc,PetscReal z)
32*9b54502bSHong Zhang {
33*9b54502bSHong Zhang   PC_ILU *lu;
34*9b54502bSHong Zhang 
35*9b54502bSHong Zhang   PetscFunctionBegin;
36*9b54502bSHong Zhang   lu                 = (PC_ILU*)pc->data;
37*9b54502bSHong Zhang   lu->info.zeropivot = z;
38*9b54502bSHong Zhang   PetscFunctionReturn(0);
39*9b54502bSHong Zhang }
40*9b54502bSHong Zhang EXTERN_C_END
41*9b54502bSHong Zhang 
42*9b54502bSHong Zhang #undef __FUNCT__
43*9b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal"
44*9b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc)
45*9b54502bSHong Zhang {
46*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
47*9b54502bSHong Zhang   PetscErrorCode ierr;
48*9b54502bSHong Zhang 
49*9b54502bSHong Zhang   PetscFunctionBegin;
50*9b54502bSHong Zhang   if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);}
51*9b54502bSHong Zhang   if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
52*9b54502bSHong Zhang   if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
53*9b54502bSHong Zhang   PetscFunctionReturn(0);
54*9b54502bSHong Zhang }
55*9b54502bSHong Zhang 
56*9b54502bSHong Zhang EXTERN_C_BEGIN
57*9b54502bSHong Zhang #undef __FUNCT__
58*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetDamping_ILU"
59*9b54502bSHong Zhang PetscErrorCode PCILUSetDamping_ILU(PC pc,PetscReal damping)
60*9b54502bSHong Zhang {
61*9b54502bSHong Zhang   PC_ILU *dir;
62*9b54502bSHong Zhang 
63*9b54502bSHong Zhang   PetscFunctionBegin;
64*9b54502bSHong Zhang   dir = (PC_ILU*)pc->data;
65*9b54502bSHong Zhang   if (damping == (PetscReal) PETSC_DECIDE) {
66*9b54502bSHong Zhang     dir->info.damping = 1.e-12;
67*9b54502bSHong Zhang   } else {
68*9b54502bSHong Zhang     dir->info.damping = damping;
69*9b54502bSHong Zhang   }
70*9b54502bSHong Zhang   PetscFunctionReturn(0);
71*9b54502bSHong Zhang }
72*9b54502bSHong Zhang EXTERN_C_END
73*9b54502bSHong Zhang 
74*9b54502bSHong Zhang EXTERN_C_BEGIN
75*9b54502bSHong Zhang #undef __FUNCT__
76*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetShift_ILU"
77*9b54502bSHong Zhang PetscErrorCode PCILUSetShift_ILU(PC pc,PetscTruth shift)
78*9b54502bSHong Zhang {
79*9b54502bSHong Zhang   PC_ILU *dir;
80*9b54502bSHong Zhang 
81*9b54502bSHong Zhang   PetscFunctionBegin;
82*9b54502bSHong Zhang   dir = (PC_ILU*)pc->data;
83*9b54502bSHong Zhang   dir->info.shift = shift;
84*9b54502bSHong Zhang   if (shift) dir->info.shift_fraction = 0.0;
85*9b54502bSHong Zhang   PetscFunctionReturn(0);
86*9b54502bSHong Zhang }
87*9b54502bSHong Zhang EXTERN_C_END
88*9b54502bSHong Zhang 
89*9b54502bSHong Zhang EXTERN_C_BEGIN
90*9b54502bSHong Zhang #undef __FUNCT__
91*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance_ILU"
92*9b54502bSHong Zhang PetscErrorCode PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount)
93*9b54502bSHong Zhang {
94*9b54502bSHong Zhang   PC_ILU         *ilu;
95*9b54502bSHong Zhang   PetscErrorCode ierr;
96*9b54502bSHong Zhang 
97*9b54502bSHong Zhang   PetscFunctionBegin;
98*9b54502bSHong Zhang   ilu = (PC_ILU*)pc->data;
99*9b54502bSHong Zhang   if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) {
100*9b54502bSHong Zhang     pc->setupcalled   = 0;
101*9b54502bSHong Zhang     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
102*9b54502bSHong Zhang   }
103*9b54502bSHong Zhang   ilu->usedt        = PETSC_TRUE;
104*9b54502bSHong Zhang   ilu->info.dt      = dt;
105*9b54502bSHong Zhang   ilu->info.dtcol   = dtcol;
106*9b54502bSHong Zhang   ilu->info.dtcount = dtcount;
107*9b54502bSHong Zhang   ilu->info.fill    = PETSC_DEFAULT;
108*9b54502bSHong Zhang   PetscFunctionReturn(0);
109*9b54502bSHong Zhang }
110*9b54502bSHong Zhang EXTERN_C_END
111*9b54502bSHong Zhang 
112*9b54502bSHong Zhang EXTERN_C_BEGIN
113*9b54502bSHong Zhang #undef __FUNCT__
114*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetFill_ILU"
115*9b54502bSHong Zhang PetscErrorCode PCILUSetFill_ILU(PC pc,PetscReal fill)
116*9b54502bSHong Zhang {
117*9b54502bSHong Zhang   PC_ILU *dir;
118*9b54502bSHong Zhang 
119*9b54502bSHong Zhang   PetscFunctionBegin;
120*9b54502bSHong Zhang   dir            = (PC_ILU*)pc->data;
121*9b54502bSHong Zhang   dir->info.fill = fill;
122*9b54502bSHong Zhang   PetscFunctionReturn(0);
123*9b54502bSHong Zhang }
124*9b54502bSHong Zhang EXTERN_C_END
125*9b54502bSHong Zhang 
126*9b54502bSHong Zhang EXTERN_C_BEGIN
127*9b54502bSHong Zhang #undef __FUNCT__
128*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering_ILU"
129*9b54502bSHong Zhang PetscErrorCode PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
130*9b54502bSHong Zhang {
131*9b54502bSHong Zhang   PC_ILU         *dir = (PC_ILU*)pc->data;
132*9b54502bSHong Zhang   PetscErrorCode ierr;
133*9b54502bSHong Zhang   PetscTruth     flg;
134*9b54502bSHong Zhang 
135*9b54502bSHong Zhang   PetscFunctionBegin;
136*9b54502bSHong Zhang   if (!pc->setupcalled) {
137*9b54502bSHong Zhang      ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
138*9b54502bSHong Zhang      ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
139*9b54502bSHong Zhang   } else {
140*9b54502bSHong Zhang     ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr);
141*9b54502bSHong Zhang     if (!flg) {
142*9b54502bSHong Zhang       pc->setupcalled = 0;
143*9b54502bSHong Zhang       ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
144*9b54502bSHong Zhang       ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
145*9b54502bSHong Zhang       /* free the data structures, then create them again */
146*9b54502bSHong Zhang       ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
147*9b54502bSHong Zhang     }
148*9b54502bSHong Zhang   }
149*9b54502bSHong Zhang   PetscFunctionReturn(0);
150*9b54502bSHong Zhang }
151*9b54502bSHong Zhang EXTERN_C_END
152*9b54502bSHong Zhang 
153*9b54502bSHong Zhang EXTERN_C_BEGIN
154*9b54502bSHong Zhang #undef __FUNCT__
155*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering_ILU"
156*9b54502bSHong Zhang PetscErrorCode PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
157*9b54502bSHong Zhang {
158*9b54502bSHong Zhang   PC_ILU *ilu;
159*9b54502bSHong Zhang 
160*9b54502bSHong Zhang   PetscFunctionBegin;
161*9b54502bSHong Zhang   ilu                = (PC_ILU*)pc->data;
162*9b54502bSHong Zhang   ilu->reuseordering = flag;
163*9b54502bSHong Zhang   PetscFunctionReturn(0);
164*9b54502bSHong Zhang }
165*9b54502bSHong Zhang EXTERN_C_END
166*9b54502bSHong Zhang 
167*9b54502bSHong Zhang EXTERN_C_BEGIN
168*9b54502bSHong Zhang #undef __FUNCT__
169*9b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill_ILUDT"
170*9b54502bSHong Zhang PetscErrorCode PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
171*9b54502bSHong Zhang {
172*9b54502bSHong Zhang   PC_ILU *ilu;
173*9b54502bSHong Zhang 
174*9b54502bSHong Zhang   PetscFunctionBegin;
175*9b54502bSHong Zhang   ilu = (PC_ILU*)pc->data;
176*9b54502bSHong Zhang   ilu->reusefill = flag;
177*9b54502bSHong Zhang   if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported");
178*9b54502bSHong Zhang   PetscFunctionReturn(0);
179*9b54502bSHong Zhang }
180*9b54502bSHong Zhang EXTERN_C_END
181*9b54502bSHong Zhang 
182*9b54502bSHong Zhang EXTERN_C_BEGIN
183*9b54502bSHong Zhang #undef __FUNCT__
184*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels_ILU"
185*9b54502bSHong Zhang PetscErrorCode PCILUSetLevels_ILU(PC pc,PetscInt levels)
186*9b54502bSHong Zhang {
187*9b54502bSHong Zhang   PC_ILU         *ilu;
188*9b54502bSHong Zhang   PetscErrorCode ierr;
189*9b54502bSHong Zhang 
190*9b54502bSHong Zhang   PetscFunctionBegin;
191*9b54502bSHong Zhang   ilu = (PC_ILU*)pc->data;
192*9b54502bSHong Zhang 
193*9b54502bSHong Zhang   if (!pc->setupcalled) {
194*9b54502bSHong Zhang     ilu->info.levels = levels;
195*9b54502bSHong Zhang   } else if (ilu->usedt || ilu->info.levels != levels) {
196*9b54502bSHong Zhang     ilu->info.levels = levels;
197*9b54502bSHong Zhang     pc->setupcalled  = 0;
198*9b54502bSHong Zhang     ilu->usedt       = PETSC_FALSE;
199*9b54502bSHong Zhang     ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
200*9b54502bSHong Zhang   }
201*9b54502bSHong Zhang   PetscFunctionReturn(0);
202*9b54502bSHong Zhang }
203*9b54502bSHong Zhang EXTERN_C_END
204*9b54502bSHong Zhang 
205*9b54502bSHong Zhang EXTERN_C_BEGIN
206*9b54502bSHong Zhang #undef __FUNCT__
207*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace_ILU"
208*9b54502bSHong Zhang PetscErrorCode PCILUSetUseInPlace_ILU(PC pc)
209*9b54502bSHong Zhang {
210*9b54502bSHong Zhang   PC_ILU *dir;
211*9b54502bSHong Zhang 
212*9b54502bSHong Zhang   PetscFunctionBegin;
213*9b54502bSHong Zhang   dir          = (PC_ILU*)pc->data;
214*9b54502bSHong Zhang   dir->inplace = PETSC_TRUE;
215*9b54502bSHong Zhang   PetscFunctionReturn(0);
216*9b54502bSHong Zhang }
217*9b54502bSHong Zhang EXTERN_C_END
218*9b54502bSHong Zhang 
219*9b54502bSHong Zhang EXTERN_C_BEGIN
220*9b54502bSHong Zhang #undef __FUNCT__
221*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill"
222*9b54502bSHong Zhang PetscErrorCode PCILUSetAllowDiagonalFill_ILU(PC pc)
223*9b54502bSHong Zhang {
224*9b54502bSHong Zhang   PC_ILU *dir;
225*9b54502bSHong Zhang 
226*9b54502bSHong Zhang   PetscFunctionBegin;
227*9b54502bSHong Zhang   dir                 = (PC_ILU*)pc->data;
228*9b54502bSHong Zhang   dir->info.diagonal_fill = 1;
229*9b54502bSHong Zhang   PetscFunctionReturn(0);
230*9b54502bSHong Zhang }
231*9b54502bSHong Zhang EXTERN_C_END
232*9b54502bSHong Zhang 
233*9b54502bSHong Zhang #undef __FUNCT__
234*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetZeroPivot"
235*9b54502bSHong Zhang /*@
236*9b54502bSHong Zhang    PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero
237*9b54502bSHong Zhang 
238*9b54502bSHong Zhang    Collective on PC
239*9b54502bSHong Zhang 
240*9b54502bSHong Zhang    Input Parameters:
241*9b54502bSHong Zhang +  pc - the preconditioner context
242*9b54502bSHong Zhang -  zero - all pivots smaller than this will be considered zero
243*9b54502bSHong Zhang 
244*9b54502bSHong Zhang    Options Database Key:
245*9b54502bSHong Zhang .  -pc_ilu_zeropivot <zero> - Sets the zero pivot size
246*9b54502bSHong Zhang 
247*9b54502bSHong Zhang    Level: intermediate
248*9b54502bSHong Zhang 
249*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
250*9b54502bSHong Zhang 
251*9b54502bSHong Zhang .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot()
252*9b54502bSHong Zhang @*/
253*9b54502bSHong Zhang PetscErrorCode PCILUSetZeroPivot(PC pc,PetscReal zero)
254*9b54502bSHong Zhang {
255*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
256*9b54502bSHong Zhang 
257*9b54502bSHong Zhang   PetscFunctionBegin;
258*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
259*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr);
260*9b54502bSHong Zhang   if (f) {
261*9b54502bSHong Zhang     ierr = (*f)(pc,zero);CHKERRQ(ierr);
262*9b54502bSHong Zhang   }
263*9b54502bSHong Zhang   PetscFunctionReturn(0);
264*9b54502bSHong Zhang }
265*9b54502bSHong Zhang 
266*9b54502bSHong Zhang #undef __FUNCT__
267*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetDamping"
268*9b54502bSHong Zhang /*@
269*9b54502bSHong Zhang    PCILUSetDamping - adds this quantity to the diagonal of the matrix during the
270*9b54502bSHong Zhang      ILU numerical factorization
271*9b54502bSHong Zhang 
272*9b54502bSHong Zhang    Collective on PC
273*9b54502bSHong Zhang 
274*9b54502bSHong Zhang    Input Parameters:
275*9b54502bSHong Zhang +  pc - the preconditioner context
276*9b54502bSHong Zhang -  damping - amount of damping
277*9b54502bSHong Zhang 
278*9b54502bSHong Zhang    Options Database Key:
279*9b54502bSHong Zhang .  -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default
280*9b54502bSHong Zhang 
281*9b54502bSHong Zhang    Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero
282*9b54502bSHong Zhang          pivot, then the damping is doubled until this is alleviated.
283*9b54502bSHong Zhang 
284*9b54502bSHong Zhang    Level: intermediate
285*9b54502bSHong Zhang 
286*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
287*9b54502bSHong Zhang 
288*9b54502bSHong Zhang .seealso: PCILUSetFill(), PCLUSetDamping(), PCILUSetShift()
289*9b54502bSHong Zhang @*/
290*9b54502bSHong Zhang PetscErrorCode PCILUSetDamping(PC pc,PetscReal damping)
291*9b54502bSHong Zhang {
292*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
293*9b54502bSHong Zhang 
294*9b54502bSHong Zhang   PetscFunctionBegin;
295*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
296*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
297*9b54502bSHong Zhang   if (f) {
298*9b54502bSHong Zhang     ierr = (*f)(pc,damping);CHKERRQ(ierr);
299*9b54502bSHong Zhang   }
300*9b54502bSHong Zhang   PetscFunctionReturn(0);
301*9b54502bSHong Zhang }
302*9b54502bSHong Zhang 
303*9b54502bSHong Zhang #undef __FUNCT__
304*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetShift"
305*9b54502bSHong Zhang /*@
306*9b54502bSHong Zhang    PCILUSetShift - specify whether to use Manteuffel shifting of ILU.
307*9b54502bSHong Zhang    If an ILU factorisation breaks down because of nonpositive pivots,
308*9b54502bSHong Zhang    adding sufficient identity to the diagonal will remedy this.
309*9b54502bSHong Zhang    Setting this causes a bisection method to find the minimum shift that
310*9b54502bSHong Zhang    will lead to a well-defined ILU.
311*9b54502bSHong Zhang 
312*9b54502bSHong Zhang    Input parameters:
313*9b54502bSHong Zhang +  pc - the preconditioner context
314*9b54502bSHong Zhang -  shifting - PETSC_TRUE to set shift else PETSC_FALSE
315*9b54502bSHong Zhang 
316*9b54502bSHong Zhang    Options Database Key:
317*9b54502bSHong Zhang .  -pc_ilu_shift [1/0] - Activate/Deactivate PCILUSetShift(); the value
318*9b54502bSHong Zhang    is optional with 1 being the default
319*9b54502bSHong Zhang 
320*9b54502bSHong Zhang    Level: intermediate
321*9b54502bSHong Zhang 
322*9b54502bSHong Zhang .keywords: PC, indefinite, factorization, incomplete, ILU
323*9b54502bSHong Zhang 
324*9b54502bSHong Zhang .seealso: PCILUSetDamping()
325*9b54502bSHong Zhang @*/
326*9b54502bSHong Zhang PetscErrorCode PCILUSetShift(PC pc,PetscTruth shifting)
327*9b54502bSHong Zhang {
328*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
329*9b54502bSHong Zhang 
330*9b54502bSHong Zhang   PetscFunctionBegin;
331*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
332*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr);
333*9b54502bSHong Zhang   if (f) {
334*9b54502bSHong Zhang     ierr = (*f)(pc,shifting);CHKERRQ(ierr);
335*9b54502bSHong Zhang   }
336*9b54502bSHong Zhang   PetscFunctionReturn(0);
337*9b54502bSHong Zhang }
338*9b54502bSHong Zhang 
339*9b54502bSHong Zhang 
340*9b54502bSHong Zhang #undef __FUNCT__
341*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance"
342*9b54502bSHong Zhang /*@
343*9b54502bSHong Zhang    PCILUSetUseDropTolerance - The preconditioner will use an ILU
344*9b54502bSHong Zhang    based on a drop tolerance.
345*9b54502bSHong Zhang 
346*9b54502bSHong Zhang    Collective on PC
347*9b54502bSHong Zhang 
348*9b54502bSHong Zhang    Input Parameters:
349*9b54502bSHong Zhang +  pc - the preconditioner context
350*9b54502bSHong Zhang .  dt - the drop tolerance, try from 1.e-10 to .1
351*9b54502bSHong Zhang .  dtcol - tolerance for column pivot, good values [0.1 to 0.01]
352*9b54502bSHong Zhang -  maxrowcount - the max number of nonzeros allowed in a row, best value
353*9b54502bSHong Zhang                  depends on the number of nonzeros in row of original matrix
354*9b54502bSHong Zhang 
355*9b54502bSHong Zhang    Options Database Key:
356*9b54502bSHong Zhang .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
357*9b54502bSHong Zhang 
358*9b54502bSHong Zhang    Level: intermediate
359*9b54502bSHong Zhang 
360*9b54502bSHong Zhang     Notes:
361*9b54502bSHong Zhang       This uses the iludt() code of Saad's SPARSKIT package
362*9b54502bSHong Zhang 
363*9b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU
364*9b54502bSHong Zhang @*/
365*9b54502bSHong Zhang PetscErrorCode PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount)
366*9b54502bSHong Zhang {
367*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt);
368*9b54502bSHong Zhang 
369*9b54502bSHong Zhang   PetscFunctionBegin;
370*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
371*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
372*9b54502bSHong Zhang   if (f) {
373*9b54502bSHong Zhang     ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr);
374*9b54502bSHong Zhang   }
375*9b54502bSHong Zhang   PetscFunctionReturn(0);
376*9b54502bSHong Zhang }
377*9b54502bSHong Zhang 
378*9b54502bSHong Zhang #undef __FUNCT__
379*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetFill"
380*9b54502bSHong Zhang /*@
381*9b54502bSHong Zhang    PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
382*9b54502bSHong Zhang    where fill = number nonzeros in factor/number nonzeros in original matrix.
383*9b54502bSHong Zhang 
384*9b54502bSHong Zhang    Collective on PC
385*9b54502bSHong Zhang 
386*9b54502bSHong Zhang    Input Parameters:
387*9b54502bSHong Zhang +  pc - the preconditioner context
388*9b54502bSHong Zhang -  fill - amount of expected fill
389*9b54502bSHong Zhang 
390*9b54502bSHong Zhang    Options Database Key:
391*9b54502bSHong Zhang $  -pc_ilu_fill <fill>
392*9b54502bSHong Zhang 
393*9b54502bSHong Zhang    Note:
394*9b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
395*9b54502bSHong Zhang    fill to expect. By running with the option -log_info PETSc will print the
396*9b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
397*9b54502bSHong Zhang    future runs. But default PETSc uses a value of 1.0
398*9b54502bSHong Zhang 
399*9b54502bSHong Zhang    Level: intermediate
400*9b54502bSHong Zhang 
401*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
402*9b54502bSHong Zhang 
403*9b54502bSHong Zhang .seealso: PCLUSetFill()
404*9b54502bSHong Zhang @*/
405*9b54502bSHong Zhang PetscErrorCode PCILUSetFill(PC pc,PetscReal fill)
406*9b54502bSHong Zhang {
407*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
408*9b54502bSHong Zhang 
409*9b54502bSHong Zhang   PetscFunctionBegin;
410*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
411*9b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
412*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
413*9b54502bSHong Zhang   if (f) {
414*9b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
415*9b54502bSHong Zhang   }
416*9b54502bSHong Zhang   PetscFunctionReturn(0);
417*9b54502bSHong Zhang }
418*9b54502bSHong Zhang 
419*9b54502bSHong Zhang #undef __FUNCT__
420*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering"
421*9b54502bSHong Zhang /*@C
422*9b54502bSHong Zhang     PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to
423*9b54502bSHong Zhang     be used in the ILU factorization.
424*9b54502bSHong Zhang 
425*9b54502bSHong Zhang     Collective on PC
426*9b54502bSHong Zhang 
427*9b54502bSHong Zhang     Input Parameters:
428*9b54502bSHong Zhang +   pc - the preconditioner context
429*9b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
430*9b54502bSHong Zhang 
431*9b54502bSHong Zhang     Options Database Key:
432*9b54502bSHong Zhang .   -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
433*9b54502bSHong Zhang 
434*9b54502bSHong Zhang     Level: intermediate
435*9b54502bSHong Zhang 
436*9b54502bSHong Zhang     Notes: natural ordering is used by default
437*9b54502bSHong Zhang 
438*9b54502bSHong Zhang .seealso: PCLUSetMatOrdering()
439*9b54502bSHong Zhang 
440*9b54502bSHong Zhang .keywords: PC, ILU, set, matrix, reordering
441*9b54502bSHong Zhang 
442*9b54502bSHong Zhang @*/
443*9b54502bSHong Zhang PetscErrorCode PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
444*9b54502bSHong Zhang {
445*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
446*9b54502bSHong Zhang 
447*9b54502bSHong Zhang   PetscFunctionBegin;
448*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
449*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
450*9b54502bSHong Zhang   if (f) {
451*9b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
452*9b54502bSHong Zhang   }
453*9b54502bSHong Zhang   PetscFunctionReturn(0);
454*9b54502bSHong Zhang }
455*9b54502bSHong Zhang 
456*9b54502bSHong Zhang #undef __FUNCT__
457*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering"
458*9b54502bSHong Zhang /*@
459*9b54502bSHong Zhang    PCILUSetReuseOrdering - When similar matrices are factored, this
460*9b54502bSHong Zhang    causes the ordering computed in the first factor to be used for all
461*9b54502bSHong Zhang    following factors; applies to both fill and drop tolerance ILUs.
462*9b54502bSHong Zhang 
463*9b54502bSHong Zhang    Collective on PC
464*9b54502bSHong Zhang 
465*9b54502bSHong Zhang    Input Parameters:
466*9b54502bSHong Zhang +  pc - the preconditioner context
467*9b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
468*9b54502bSHong Zhang 
469*9b54502bSHong Zhang    Options Database Key:
470*9b54502bSHong Zhang .  -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()
471*9b54502bSHong Zhang 
472*9b54502bSHong Zhang    Level: intermediate
473*9b54502bSHong Zhang 
474*9b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU
475*9b54502bSHong Zhang 
476*9b54502bSHong Zhang .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
477*9b54502bSHong Zhang @*/
478*9b54502bSHong Zhang PetscErrorCode PCILUSetReuseOrdering(PC pc,PetscTruth flag)
479*9b54502bSHong Zhang {
480*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
481*9b54502bSHong Zhang 
482*9b54502bSHong Zhang   PetscFunctionBegin;
483*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
484*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
485*9b54502bSHong Zhang   if (f) {
486*9b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
487*9b54502bSHong Zhang   }
488*9b54502bSHong Zhang   PetscFunctionReturn(0);
489*9b54502bSHong Zhang }
490*9b54502bSHong Zhang 
491*9b54502bSHong Zhang #undef __FUNCT__
492*9b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill"
493*9b54502bSHong Zhang /*@
494*9b54502bSHong Zhang    PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
495*9b54502bSHong Zhang    this causes later ones to use the fill computed in the initial factorization.
496*9b54502bSHong Zhang 
497*9b54502bSHong Zhang    Collective on PC
498*9b54502bSHong Zhang 
499*9b54502bSHong Zhang    Input Parameters:
500*9b54502bSHong Zhang +  pc - the preconditioner context
501*9b54502bSHong Zhang -  flag - PETSC_TRUE to reuse else PETSC_FALSE
502*9b54502bSHong Zhang 
503*9b54502bSHong Zhang    Options Database Key:
504*9b54502bSHong Zhang .  -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()
505*9b54502bSHong Zhang 
506*9b54502bSHong Zhang    Level: intermediate
507*9b54502bSHong Zhang 
508*9b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU
509*9b54502bSHong Zhang 
510*9b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
511*9b54502bSHong Zhang @*/
512*9b54502bSHong Zhang PetscErrorCode PCILUDTSetReuseFill(PC pc,PetscTruth flag)
513*9b54502bSHong Zhang {
514*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
515*9b54502bSHong Zhang 
516*9b54502bSHong Zhang   PetscFunctionBegin;
517*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
518*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr);
519*9b54502bSHong Zhang   if (f) {
520*9b54502bSHong Zhang     ierr = (*f)(pc,flag);CHKERRQ(ierr);
521*9b54502bSHong Zhang   }
522*9b54502bSHong Zhang   PetscFunctionReturn(0);
523*9b54502bSHong Zhang }
524*9b54502bSHong Zhang 
525*9b54502bSHong Zhang #undef __FUNCT__
526*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels"
527*9b54502bSHong Zhang /*@
528*9b54502bSHong Zhang    PCILUSetLevels - Sets the number of levels of fill to use.
529*9b54502bSHong Zhang 
530*9b54502bSHong Zhang    Collective on PC
531*9b54502bSHong Zhang 
532*9b54502bSHong Zhang    Input Parameters:
533*9b54502bSHong Zhang +  pc - the preconditioner context
534*9b54502bSHong Zhang -  levels - number of levels of fill
535*9b54502bSHong Zhang 
536*9b54502bSHong Zhang    Options Database Key:
537*9b54502bSHong Zhang .  -pc_ilu_levels <levels> - Sets fill level
538*9b54502bSHong Zhang 
539*9b54502bSHong Zhang    Level: intermediate
540*9b54502bSHong Zhang 
541*9b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU
542*9b54502bSHong Zhang @*/
543*9b54502bSHong Zhang PetscErrorCode PCILUSetLevels(PC pc,PetscInt levels)
544*9b54502bSHong Zhang {
545*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscInt);
546*9b54502bSHong Zhang 
547*9b54502bSHong Zhang   PetscFunctionBegin;
548*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
549*9b54502bSHong Zhang   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
550*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
551*9b54502bSHong Zhang   if (f) {
552*9b54502bSHong Zhang     ierr = (*f)(pc,levels);CHKERRQ(ierr);
553*9b54502bSHong Zhang   }
554*9b54502bSHong Zhang   PetscFunctionReturn(0);
555*9b54502bSHong Zhang }
556*9b54502bSHong Zhang 
557*9b54502bSHong Zhang #undef __FUNCT__
558*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill"
559*9b54502bSHong Zhang /*@
560*9b54502bSHong Zhang    PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be
561*9b54502bSHong Zhang    treated as level 0 fill even if there is no non-zero location.
562*9b54502bSHong Zhang 
563*9b54502bSHong Zhang    Collective on PC
564*9b54502bSHong Zhang 
565*9b54502bSHong Zhang    Input Parameters:
566*9b54502bSHong Zhang +  pc - the preconditioner context
567*9b54502bSHong Zhang 
568*9b54502bSHong Zhang    Options Database Key:
569*9b54502bSHong Zhang .  -pc_ilu_diagonal_fill
570*9b54502bSHong Zhang 
571*9b54502bSHong Zhang    Notes:
572*9b54502bSHong Zhang    Does not apply with 0 fill.
573*9b54502bSHong Zhang 
574*9b54502bSHong Zhang    Level: intermediate
575*9b54502bSHong Zhang 
576*9b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU
577*9b54502bSHong Zhang @*/
578*9b54502bSHong Zhang PetscErrorCode PCILUSetAllowDiagonalFill(PC pc)
579*9b54502bSHong Zhang {
580*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
581*9b54502bSHong Zhang 
582*9b54502bSHong Zhang   PetscFunctionBegin;
583*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
584*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr);
585*9b54502bSHong Zhang   if (f) {
586*9b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
587*9b54502bSHong Zhang   }
588*9b54502bSHong Zhang   PetscFunctionReturn(0);
589*9b54502bSHong Zhang }
590*9b54502bSHong Zhang 
591*9b54502bSHong Zhang #undef __FUNCT__
592*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace"
593*9b54502bSHong Zhang /*@
594*9b54502bSHong Zhang    PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
595*9b54502bSHong Zhang    Collective on PC
596*9b54502bSHong Zhang 
597*9b54502bSHong Zhang    Input Parameters:
598*9b54502bSHong Zhang .  pc - the preconditioner context
599*9b54502bSHong Zhang 
600*9b54502bSHong Zhang    Options Database Key:
601*9b54502bSHong Zhang .  -pc_ilu_in_place - Activates in-place factorization
602*9b54502bSHong Zhang 
603*9b54502bSHong Zhang    Notes:
604*9b54502bSHong Zhang    PCILUSetUseInPlace() is intended for use with matrix-free variants of
605*9b54502bSHong Zhang    Krylov methods, or when a different matrices are employed for the linear
606*9b54502bSHong Zhang    system and preconditioner, or with ASM preconditioning.  Do NOT use
607*9b54502bSHong Zhang    this option if the linear system
608*9b54502bSHong Zhang    matrix also serves as the preconditioning matrix, since the factored
609*9b54502bSHong Zhang    matrix would then overwrite the original matrix.
610*9b54502bSHong Zhang 
611*9b54502bSHong Zhang    Only works well with ILU(0).
612*9b54502bSHong Zhang 
613*9b54502bSHong Zhang    Level: intermediate
614*9b54502bSHong Zhang 
615*9b54502bSHong Zhang .keywords: PC, set, factorization, inplace, in-place, ILU
616*9b54502bSHong Zhang 
617*9b54502bSHong Zhang .seealso:  PCLUSetUseInPlace()
618*9b54502bSHong Zhang @*/
619*9b54502bSHong Zhang PetscErrorCode PCILUSetUseInPlace(PC pc)
620*9b54502bSHong Zhang {
621*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC);
622*9b54502bSHong Zhang 
623*9b54502bSHong Zhang   PetscFunctionBegin;
624*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
625*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr);
626*9b54502bSHong Zhang   if (f) {
627*9b54502bSHong Zhang     ierr = (*f)(pc);CHKERRQ(ierr);
628*9b54502bSHong Zhang   }
629*9b54502bSHong Zhang   PetscFunctionReturn(0);
630*9b54502bSHong Zhang }
631*9b54502bSHong Zhang 
632*9b54502bSHong Zhang #undef __FUNCT__
633*9b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal"
634*9b54502bSHong Zhang /*@
635*9b54502bSHong Zhang    PCILUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal
636*9b54502bSHong Zhang 
637*9b54502bSHong Zhang    Collective on PC
638*9b54502bSHong Zhang 
639*9b54502bSHong Zhang    Input Parameters:
640*9b54502bSHong Zhang +  pc - the preconditioner context
641*9b54502bSHong Zhang -  tol - diagonal entries smaller than this in absolute value are considered zero
642*9b54502bSHong Zhang 
643*9b54502bSHong Zhang    Options Database Key:
644*9b54502bSHong Zhang .  -pc_lu_nonzeros_along_diagonal
645*9b54502bSHong Zhang 
646*9b54502bSHong Zhang    Level: intermediate
647*9b54502bSHong Zhang 
648*9b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
649*9b54502bSHong Zhang 
650*9b54502bSHong Zhang .seealso: PCILUSetFill(), PCILUSetDamp(), PCIILUSetZeroPivot(), MatReorderForNonzeroDiagonal()
651*9b54502bSHong Zhang @*/
652*9b54502bSHong Zhang PetscErrorCode PCILUReorderForNonzeroDiagonal(PC pc,PetscReal rtol)
653*9b54502bSHong Zhang {
654*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
655*9b54502bSHong Zhang 
656*9b54502bSHong Zhang   PetscFunctionBegin;
657*9b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
658*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
659*9b54502bSHong Zhang   if (f) {
660*9b54502bSHong Zhang     ierr = (*f)(pc,rtol);CHKERRQ(ierr);
661*9b54502bSHong Zhang   }
662*9b54502bSHong Zhang   PetscFunctionReturn(0);
663*9b54502bSHong Zhang }
664*9b54502bSHong Zhang 
665*9b54502bSHong Zhang #undef __FUNCT__
666*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks"
667*9b54502bSHong Zhang /*@
668*9b54502bSHong Zhang     PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block
669*9b54502bSHong Zhang       with BAIJ or SBAIJ matrices
670*9b54502bSHong Zhang 
671*9b54502bSHong Zhang     Collective on PC
672*9b54502bSHong Zhang 
673*9b54502bSHong Zhang     Input Parameters:
674*9b54502bSHong Zhang +   pc - the preconditioner context
675*9b54502bSHong Zhang -   pivot - PETSC_TRUE or PETSC_FALSE
676*9b54502bSHong Zhang 
677*9b54502bSHong Zhang     Options Database Key:
678*9b54502bSHong Zhang .   -pc_ilu_pivot_in_blocks <true,false>
679*9b54502bSHong Zhang 
680*9b54502bSHong Zhang     Level: intermediate
681*9b54502bSHong Zhang 
682*9b54502bSHong Zhang .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting()
683*9b54502bSHong Zhang @*/
684*9b54502bSHong Zhang PetscErrorCode PCILUSetPivotInBlocks(PC pc,PetscTruth pivot)
685*9b54502bSHong Zhang {
686*9b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscTruth);
687*9b54502bSHong Zhang 
688*9b54502bSHong Zhang   PetscFunctionBegin;
689*9b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr);
690*9b54502bSHong Zhang   if (f) {
691*9b54502bSHong Zhang     ierr = (*f)(pc,pivot);CHKERRQ(ierr);
692*9b54502bSHong Zhang   }
693*9b54502bSHong Zhang   PetscFunctionReturn(0);
694*9b54502bSHong Zhang }
695*9b54502bSHong Zhang 
696*9b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/
697*9b54502bSHong Zhang 
698*9b54502bSHong Zhang EXTERN_C_BEGIN
699*9b54502bSHong Zhang #undef __FUNCT__
700*9b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks_ILU"
701*9b54502bSHong Zhang PetscErrorCode PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot)
702*9b54502bSHong Zhang {
703*9b54502bSHong Zhang   PC_ILU *dir = (PC_ILU*)pc->data;
704*9b54502bSHong Zhang 
705*9b54502bSHong Zhang   PetscFunctionBegin;
706*9b54502bSHong Zhang   dir->info.pivotinblocks = pivot ? 1.0 : 0.0;
707*9b54502bSHong Zhang   PetscFunctionReturn(0);
708*9b54502bSHong Zhang }
709*9b54502bSHong Zhang EXTERN_C_END
710*9b54502bSHong Zhang 
711*9b54502bSHong Zhang #undef __FUNCT__
712*9b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU"
713*9b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc)
714*9b54502bSHong Zhang {
715*9b54502bSHong Zhang   PetscErrorCode ierr;
716*9b54502bSHong Zhang   PetscInt       dtmax = 3,itmp;
717*9b54502bSHong Zhang   PetscTruth     flg,set;
718*9b54502bSHong Zhang   PetscReal      dt[3];
719*9b54502bSHong Zhang   char           tname[256];
720*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
721*9b54502bSHong Zhang   PetscFList     ordlist;
722*9b54502bSHong Zhang   PetscReal      tol;
723*9b54502bSHong Zhang 
724*9b54502bSHong Zhang   PetscFunctionBegin;
725*9b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
726*9b54502bSHong Zhang   ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr);
727*9b54502bSHong Zhang     ierr = PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr);
728*9b54502bSHong Zhang     if (flg) ilu->info.levels = itmp;
729*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);CHKERRQ(ierr);
730*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);CHKERRQ(ierr);
731*9b54502bSHong Zhang     ilu->info.diagonal_fill = (double) flg;
732*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);CHKERRQ(ierr);
733*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr);
734*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",&flg);CHKERRQ(ierr);
735*9b54502bSHong Zhang     if (flg) {
736*9b54502bSHong Zhang       ierr = PCILUSetDamping(pc,(PetscReal) PETSC_DECIDE);CHKERRQ(ierr);
737*9b54502bSHong Zhang     }
738*9b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);CHKERRQ(ierr);
739*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);CHKERRQ(ierr);
740*9b54502bSHong Zhang     if (flg) {
741*9b54502bSHong Zhang       ierr = PetscOptionsInt("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",(PetscInt)ilu->info.shift,&itmp,&flg); CHKERRQ(ierr);
742*9b54502bSHong Zhang       if (flg && !itmp) {
743*9b54502bSHong Zhang 	ierr = PCILUSetShift(pc,PETSC_FALSE);CHKERRQ(ierr);
744*9b54502bSHong Zhang       } else {
745*9b54502bSHong Zhang 	ierr = PCILUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr);
746*9b54502bSHong Zhang       }
747*9b54502bSHong Zhang     }
748*9b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr);
749*9b54502bSHong Zhang 
750*9b54502bSHong Zhang     dt[0] = ilu->info.dt;
751*9b54502bSHong Zhang     dt[1] = ilu->info.dtcol;
752*9b54502bSHong Zhang     dt[2] = ilu->info.dtcount;
753*9b54502bSHong Zhang     ierr = PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
754*9b54502bSHong Zhang     if (flg) {
755*9b54502bSHong Zhang       ierr = PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
756*9b54502bSHong Zhang     }
757*9b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr);
758*9b54502bSHong Zhang     ierr = PetscOptionsName("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr);
759*9b54502bSHong Zhang     if (flg) {
760*9b54502bSHong Zhang       tol = PETSC_DECIDE;
761*9b54502bSHong Zhang       ierr = PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr);
762*9b54502bSHong Zhang       ierr = PCILUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr);
763*9b54502bSHong Zhang     }
764*9b54502bSHong Zhang 
765*9b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
766*9b54502bSHong Zhang     ierr = PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr);
767*9b54502bSHong Zhang     if (flg) {
768*9b54502bSHong Zhang       ierr = PCILUSetMatOrdering(pc,tname);CHKERRQ(ierr);
769*9b54502bSHong Zhang     }
770*9b54502bSHong Zhang     flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE;
771*9b54502bSHong Zhang     ierr = PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr);
772*9b54502bSHong Zhang     if (set) {
773*9b54502bSHong Zhang       ierr = PCILUSetPivotInBlocks(pc,flg);CHKERRQ(ierr);
774*9b54502bSHong Zhang     }
775*9b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
776*9b54502bSHong Zhang   PetscFunctionReturn(0);
777*9b54502bSHong Zhang }
778*9b54502bSHong Zhang 
779*9b54502bSHong Zhang #undef __FUNCT__
780*9b54502bSHong Zhang #define __FUNCT__ "PCView_ILU"
781*9b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer)
782*9b54502bSHong Zhang {
783*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
784*9b54502bSHong Zhang   PetscErrorCode ierr;
785*9b54502bSHong Zhang   PetscTruth     isstring,iascii;
786*9b54502bSHong Zhang 
787*9b54502bSHong Zhang   PetscFunctionBegin;
788*9b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
789*9b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
790*9b54502bSHong Zhang   if (iascii) {
791*9b54502bSHong Zhang     if (ilu->usedt) {
792*9b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: drop tolerance %g\n",ilu->info.dt);CHKERRQ(ierr);
793*9b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr);
794*9b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: column permutation tolerance %g\n",ilu->info.dtcol);CHKERRQ(ierr);
795*9b54502bSHong Zhang     } else if (ilu->info.levels == 1) {
796*9b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
797*9b54502bSHong Zhang     } else {
798*9b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr);
799*9b54502bSHong Zhang     }
800*9b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: max fill ratio allocated %g\n",ilu->info.fill);CHKERRQ(ierr);
801*9b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);CHKERRQ(ierr);
802*9b54502bSHong Zhang     if (ilu->info.shift) {ierr = PetscViewerASCIIPrintf(viewer,"  ILU: using Manteuffel shift\n");CHKERRQ(ierr);}
803*9b54502bSHong Zhang     if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer,"       in-place factorization\n");CHKERRQ(ierr);}
804*9b54502bSHong Zhang     else              {ierr = PetscViewerASCIIPrintf(viewer,"       out-of-place factorization\n");CHKERRQ(ierr);}
805*9b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"       matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr);
806*9b54502bSHong Zhang     if (ilu->reusefill)     {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing fill from past factorization\n");CHKERRQ(ierr);}
807*9b54502bSHong Zhang     if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer,"       Reusing reordering from past factorization\n");CHKERRQ(ierr);}
808*9b54502bSHong Zhang     if (ilu->fact) {
809*9b54502bSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
810*9b54502bSHong Zhang       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
811*9b54502bSHong Zhang       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
812*9b54502bSHong Zhang       ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr);
813*9b54502bSHong Zhang       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
814*9b54502bSHong Zhang       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
815*9b54502bSHong Zhang     }
816*9b54502bSHong Zhang   } else if (isstring) {
817*9b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr);
818*9b54502bSHong Zhang   } else {
819*9b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
820*9b54502bSHong Zhang   }
821*9b54502bSHong Zhang   PetscFunctionReturn(0);
822*9b54502bSHong Zhang }
823*9b54502bSHong Zhang 
824*9b54502bSHong Zhang #undef __FUNCT__
825*9b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU"
826*9b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc)
827*9b54502bSHong Zhang {
828*9b54502bSHong Zhang   PetscErrorCode ierr;
829*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
830*9b54502bSHong Zhang 
831*9b54502bSHong Zhang   PetscFunctionBegin;
832*9b54502bSHong Zhang   if (ilu->inplace) {
833*9b54502bSHong Zhang     if (!pc->setupcalled) {
834*9b54502bSHong Zhang 
835*9b54502bSHong Zhang       /* In-place factorization only makes sense with the natural ordering,
836*9b54502bSHong Zhang          so we only need to get the ordering once, even if nonzero structure changes */
837*9b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
838*9b54502bSHong Zhang       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
839*9b54502bSHong Zhang       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
840*9b54502bSHong Zhang     }
841*9b54502bSHong Zhang 
842*9b54502bSHong Zhang     /* In place ILU only makes sense with fill factor of 1.0 because
843*9b54502bSHong Zhang        cannot have levels of fill */
844*9b54502bSHong Zhang     ilu->info.fill          = 1.0;
845*9b54502bSHong Zhang     ilu->info.diagonal_fill = 0;
846*9b54502bSHong Zhang     ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr);
847*9b54502bSHong Zhang     ilu->fact = pc->pmat;
848*9b54502bSHong Zhang   } else if (ilu->usedt) {
849*9b54502bSHong Zhang     if (!pc->setupcalled) {
850*9b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
851*9b54502bSHong Zhang       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
852*9b54502bSHong Zhang       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
853*9b54502bSHong Zhang       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
854*9b54502bSHong Zhang       PetscLogObjectParent(pc,ilu->fact);
855*9b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
856*9b54502bSHong Zhang       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
857*9b54502bSHong Zhang       if (!ilu->reuseordering) {
858*9b54502bSHong Zhang         if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);}
859*9b54502bSHong Zhang         if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);}
860*9b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
861*9b54502bSHong Zhang         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
862*9b54502bSHong Zhang         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
863*9b54502bSHong Zhang       }
864*9b54502bSHong Zhang       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
865*9b54502bSHong Zhang       PetscLogObjectParent(pc,ilu->fact);
866*9b54502bSHong Zhang     } else if (!ilu->reusefill) {
867*9b54502bSHong Zhang       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
868*9b54502bSHong Zhang       ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr);
869*9b54502bSHong Zhang       PetscLogObjectParent(pc,ilu->fact);
870*9b54502bSHong Zhang     } else {
871*9b54502bSHong Zhang       ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
872*9b54502bSHong Zhang     }
873*9b54502bSHong Zhang    } else {
874*9b54502bSHong Zhang     if (!pc->setupcalled) {
875*9b54502bSHong Zhang       /* first time in so compute reordering and symbolic factorization */
876*9b54502bSHong Zhang       ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
877*9b54502bSHong Zhang       if (ilu->row) PetscLogObjectParent(pc,ilu->row);
878*9b54502bSHong Zhang       if (ilu->col) PetscLogObjectParent(pc,ilu->col);
879*9b54502bSHong Zhang       /*  Remove zeros along diagonal?     */
880*9b54502bSHong Zhang       if (ilu->nonzerosalongdiagonal) {
881*9b54502bSHong Zhang         ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
882*9b54502bSHong Zhang       }
883*9b54502bSHong Zhang       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
884*9b54502bSHong Zhang       PetscLogObjectParent(pc,ilu->fact);
885*9b54502bSHong Zhang     } else if (pc->flag != SAME_NONZERO_PATTERN) {
886*9b54502bSHong Zhang       if (!ilu->reuseordering) {
887*9b54502bSHong Zhang         /* compute a new ordering for the ILU */
888*9b54502bSHong Zhang         ierr = ISDestroy(ilu->row);CHKERRQ(ierr);
889*9b54502bSHong Zhang         ierr = ISDestroy(ilu->col);CHKERRQ(ierr);
890*9b54502bSHong Zhang         ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr);
891*9b54502bSHong Zhang         if (ilu->row) PetscLogObjectParent(pc,ilu->row);
892*9b54502bSHong Zhang         if (ilu->col) PetscLogObjectParent(pc,ilu->col);
893*9b54502bSHong Zhang         /*  Remove zeros along diagonal?     */
894*9b54502bSHong Zhang         if (ilu->nonzerosalongdiagonal) {
895*9b54502bSHong Zhang           ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr);
896*9b54502bSHong Zhang         }
897*9b54502bSHong Zhang       }
898*9b54502bSHong Zhang       ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);
899*9b54502bSHong Zhang       ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr);
900*9b54502bSHong Zhang       PetscLogObjectParent(pc,ilu->fact);
901*9b54502bSHong Zhang     }
902*9b54502bSHong Zhang     ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr);
903*9b54502bSHong Zhang   }
904*9b54502bSHong Zhang   PetscFunctionReturn(0);
905*9b54502bSHong Zhang }
906*9b54502bSHong Zhang 
907*9b54502bSHong Zhang #undef __FUNCT__
908*9b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU"
909*9b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc)
910*9b54502bSHong Zhang {
911*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
912*9b54502bSHong Zhang   PetscErrorCode ierr;
913*9b54502bSHong Zhang 
914*9b54502bSHong Zhang   PetscFunctionBegin;
915*9b54502bSHong Zhang   ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr);
916*9b54502bSHong Zhang   ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr);
917*9b54502bSHong Zhang   ierr = PetscFree(ilu);CHKERRQ(ierr);
918*9b54502bSHong Zhang   PetscFunctionReturn(0);
919*9b54502bSHong Zhang }
920*9b54502bSHong Zhang 
921*9b54502bSHong Zhang #undef __FUNCT__
922*9b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU"
923*9b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y)
924*9b54502bSHong Zhang {
925*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
926*9b54502bSHong Zhang   PetscErrorCode ierr;
927*9b54502bSHong Zhang 
928*9b54502bSHong Zhang   PetscFunctionBegin;
929*9b54502bSHong Zhang   ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr);
930*9b54502bSHong Zhang   PetscFunctionReturn(0);
931*9b54502bSHong Zhang }
932*9b54502bSHong Zhang 
933*9b54502bSHong Zhang #undef __FUNCT__
934*9b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU"
935*9b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
936*9b54502bSHong Zhang {
937*9b54502bSHong Zhang   PC_ILU         *ilu = (PC_ILU*)pc->data;
938*9b54502bSHong Zhang   PetscErrorCode ierr;
939*9b54502bSHong Zhang 
940*9b54502bSHong Zhang   PetscFunctionBegin;
941*9b54502bSHong Zhang   ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr);
942*9b54502bSHong Zhang   PetscFunctionReturn(0);
943*9b54502bSHong Zhang }
944*9b54502bSHong Zhang 
945*9b54502bSHong Zhang #undef __FUNCT__
946*9b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ILU"
947*9b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
948*9b54502bSHong Zhang {
949*9b54502bSHong Zhang   PC_ILU *ilu = (PC_ILU*)pc->data;
950*9b54502bSHong Zhang 
951*9b54502bSHong Zhang   PetscFunctionBegin;
952*9b54502bSHong Zhang   if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()");
953*9b54502bSHong Zhang   *mat = ilu->fact;
954*9b54502bSHong Zhang   PetscFunctionReturn(0);
955*9b54502bSHong Zhang }
956*9b54502bSHong Zhang 
957*9b54502bSHong Zhang /*MC
958*9b54502bSHong Zhang      PCILU - Incomplete factorization preconditioners.
959*9b54502bSHong Zhang 
960*9b54502bSHong Zhang    Options Database Keys:
961*9b54502bSHong Zhang +  -pc_ilu_levels <k> - number of levels of fill for ILU(k)
962*9b54502bSHong Zhang .  -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for
963*9b54502bSHong Zhang                       its factorization (overwrites original matrix)
964*9b54502bSHong Zhang .  -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill
965*9b54502bSHong Zhang .  -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization
966*9b54502bSHong Zhang .  -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots
967*9b54502bSHong Zhang .  -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner
968*9b54502bSHong Zhang .  -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot
969*9b54502bSHong Zhang .  -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt
970*9b54502bSHong Zhang .  -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
971*9b54502bSHong Zhang .  -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal,
972*9b54502bSHong Zhang                                    this decreases the chance of getting a zero pivot
973*9b54502bSHong Zhang .  -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
974*9b54502bSHong Zhang -  -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger
975*9b54502bSHong Zhang                              than 1 the diagonal blocks are factored with partial pivoting (this increases the
976*9b54502bSHong Zhang                              stability of the ILU factorization
977*9b54502bSHong Zhang 
978*9b54502bSHong Zhang    Level: beginner
979*9b54502bSHong Zhang 
980*9b54502bSHong Zhang   Concepts: incomplete factorization
981*9b54502bSHong Zhang 
982*9b54502bSHong Zhang    Notes: Only implemented for some matrix formats. Not implemented in parallel
983*9b54502bSHong Zhang 
984*9b54502bSHong Zhang           For BAIJ matrices this implements a point block ILU
985*9b54502bSHong Zhang 
986*9b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
987*9b54502bSHong Zhang            PCILUSetSetZeroPivot(), PCILUSetDamping(), PCILUSetShift(), PCILUSetUseDropTolerance(),
988*9b54502bSHong Zhang            PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(),
989*9b54502bSHong Zhang            PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(),
990*9b54502bSHong Zhang 
991*9b54502bSHong Zhang M*/
992*9b54502bSHong Zhang 
993*9b54502bSHong Zhang EXTERN_C_BEGIN
994*9b54502bSHong Zhang #undef __FUNCT__
995*9b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU"
996*9b54502bSHong Zhang PetscErrorCode PCCreate_ILU(PC pc)
997*9b54502bSHong Zhang {
998*9b54502bSHong Zhang   PetscErrorCode ierr;
999*9b54502bSHong Zhang   PC_ILU         *ilu;
1000*9b54502bSHong Zhang 
1001*9b54502bSHong Zhang   PetscFunctionBegin;
1002*9b54502bSHong Zhang   ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr);
1003*9b54502bSHong Zhang   PetscLogObjectMemory(pc,sizeof(PC_ILU));
1004*9b54502bSHong Zhang 
1005*9b54502bSHong Zhang   ilu->fact                    = 0;
1006*9b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr);
1007*9b54502bSHong Zhang   ilu->info.levels             = 0;
1008*9b54502bSHong Zhang   ilu->info.fill               = 1.0;
1009*9b54502bSHong Zhang   ilu->col                     = 0;
1010*9b54502bSHong Zhang   ilu->row                     = 0;
1011*9b54502bSHong Zhang   ilu->inplace                 = PETSC_FALSE;
1012*9b54502bSHong Zhang   ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr);
1013*9b54502bSHong Zhang   ilu->reuseordering           = PETSC_FALSE;
1014*9b54502bSHong Zhang   ilu->usedt                   = PETSC_FALSE;
1015*9b54502bSHong Zhang   ilu->info.dt                 = PETSC_DEFAULT;
1016*9b54502bSHong Zhang   ilu->info.dtcount            = PETSC_DEFAULT;
1017*9b54502bSHong Zhang   ilu->info.dtcol              = PETSC_DEFAULT;
1018*9b54502bSHong Zhang   ilu->info.damping            = 0.0;
1019*9b54502bSHong Zhang   ilu->info.shift              = PETSC_FALSE;
1020*9b54502bSHong Zhang   ilu->info.shift_fraction     = 0.0;
1021*9b54502bSHong Zhang   ilu->info.zeropivot          = 1.e-12;
1022*9b54502bSHong Zhang   ilu->info.pivotinblocks      = 1.0;
1023*9b54502bSHong Zhang   ilu->reusefill               = PETSC_FALSE;
1024*9b54502bSHong Zhang   ilu->info.diagonal_fill      = 0;
1025*9b54502bSHong Zhang   pc->data                     = (void*)ilu;
1026*9b54502bSHong Zhang 
1027*9b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ILU;
1028*9b54502bSHong Zhang   pc->ops->apply               = PCApply_ILU;
1029*9b54502bSHong Zhang   pc->ops->applytranspose      = PCApplyTranspose_ILU;
1030*9b54502bSHong Zhang   pc->ops->setup               = PCSetUp_ILU;
1031*9b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ILU;
1032*9b54502bSHong Zhang   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ILU;
1033*9b54502bSHong Zhang   pc->ops->view                = PCView_ILU;
1034*9b54502bSHong Zhang   pc->ops->applyrichardson     = 0;
1035*9b54502bSHong Zhang 
1036*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
1037*9b54502bSHong Zhang                     PCILUSetUseDropTolerance_ILU);CHKERRQ(ierr);
1038*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
1039*9b54502bSHong Zhang                     PCILUSetFill_ILU);CHKERRQ(ierr);
1040*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
1041*9b54502bSHong Zhang                     PCILUSetDamping_ILU);CHKERRQ(ierr);
1042*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU",
1043*9b54502bSHong Zhang  		    PCILUSetShift_ILU);CHKERRQ(ierr);
1044*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
1045*9b54502bSHong Zhang                     PCILUSetMatOrdering_ILU);CHKERRQ(ierr);
1046*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
1047*9b54502bSHong Zhang                     PCILUSetReuseOrdering_ILU);CHKERRQ(ierr);
1048*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
1049*9b54502bSHong Zhang                     PCILUDTSetReuseFill_ILUDT);CHKERRQ(ierr);
1050*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
1051*9b54502bSHong Zhang                     PCILUSetLevels_ILU);CHKERRQ(ierr);
1052*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
1053*9b54502bSHong Zhang                     PCILUSetUseInPlace_ILU);CHKERRQ(ierr);
1054*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
1055*9b54502bSHong Zhang                     PCILUSetAllowDiagonalFill_ILU);CHKERRQ(ierr);
1056*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU",
1057*9b54502bSHong Zhang                     PCILUSetPivotInBlocks_ILU);CHKERRQ(ierr);
1058*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU",
1059*9b54502bSHong Zhang                     PCILUSetZeroPivot_ILU);CHKERRQ(ierr);
1060*9b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C","PCILUReorderForNonzeroDiagonal_ILU",
1061*9b54502bSHong Zhang                     PCILUReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr);
1062*9b54502bSHong Zhang   PetscFunctionReturn(0);
1063*9b54502bSHong Zhang }
1064*9b54502bSHong Zhang EXTERN_C_END
1065