Actual source code: ilu.c
1: /*$Id: ilu.c,v 1.172 2001/08/07 21:30:27 bsmith Exp $*/
2: /*
3: Defines a ILU factorization preconditioner for any Mat implementation
4: */
5: #include src/sles/pc/pcimpl.h
6: #include src/sles/pc/impls/ilu/ilu.h
7: #include src/mat/matimpl.h
9: /* ------------------------------------------------------------------------------------------*/
10: EXTERN_C_BEGIN
11: int PCILUSetSinglePrecisionSolves_ILU(PC pc,PetscTruth flag)
12: {
13: PC_ILU *dir;
16: dir = (PC_ILU*)pc->data;
17: dir->single_precision_solve = flag;
18: return(0);
19: }
20: EXTERN_C_END
22: EXTERN_C_BEGIN
23: int PCILUSetDamping_ILU(PC pc,PetscReal damping)
24: {
25: PC_ILU *dir;
28: dir = (PC_ILU*)pc->data;
29: dir->info.damping = damping;
30: dir->info.damp = 1.0;
31: return(0);
32: }
33: EXTERN_C_END
35: EXTERN_C_BEGIN
36: int PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,int dtcount)
37: {
38: PC_ILU *ilu;
41: ilu = (PC_ILU*)pc->data;
42: ilu->usedt = PETSC_TRUE;
43: ilu->info.dt = dt;
44: ilu->info.dtcol = dtcol;
45: ilu->info.dtcount = dtcount;
46: ilu->info.fill = PETSC_DEFAULT;
47: return(0);
48: }
49: EXTERN_C_END
51: EXTERN_C_BEGIN
52: int PCILUSetFill_ILU(PC pc,PetscReal fill)
53: {
54: PC_ILU *dir;
57: dir = (PC_ILU*)pc->data;
58: dir->info.fill = fill;
59: return(0);
60: }
61: EXTERN_C_END
63: EXTERN_C_BEGIN
64: int PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering)
65: {
66: PC_ILU *dir = (PC_ILU*)pc->data;
67: int ierr;
68:
70: PetscStrfree(dir->ordering);
71: PetscStrallocpy(ordering,&dir->ordering);
72: return(0);
73: }
74: EXTERN_C_END
76: EXTERN_C_BEGIN
77: int PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag)
78: {
79: PC_ILU *ilu;
82: ilu = (PC_ILU*)pc->data;
83: ilu->reuseordering = flag;
84: return(0);
85: }
86: EXTERN_C_END
88: EXTERN_C_BEGIN
89: int PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag)
90: {
91: PC_ILU *ilu;
94: ilu = (PC_ILU*)pc->data;
95: ilu->reusefill = flag;
96: if (flag) SETERRQ(1,"Not yet supported");
97: return(0);
98: }
99: EXTERN_C_END
101: EXTERN_C_BEGIN
102: int PCILUSetLevels_ILU(PC pc,int levels)
103: {
104: PC_ILU *ilu;
107: ilu = (PC_ILU*)pc->data;
108: ilu->info.levels = levels;
109: return(0);
110: }
111: EXTERN_C_END
113: EXTERN_C_BEGIN
114: int PCILUSetUseInPlace_ILU(PC pc)
115: {
116: PC_ILU *dir;
119: dir = (PC_ILU*)pc->data;
120: dir->inplace = PETSC_TRUE;
121: return(0);
122: }
123: EXTERN_C_END
125: EXTERN_C_BEGIN
126: int PCILUSetAllowDiagonalFill_ILU(PC pc)
127: {
128: PC_ILU *dir;
131: dir = (PC_ILU*)pc->data;
132: dir->info.diagonal_fill = 1;
133: return(0);
134: }
135: EXTERN_C_END
137: /* ------------------------------------------------------------------------------------------*/
138: /*@
139: PCILUSetSinglePrecisionSolves - Sets precision of triangular solves to single precision
141: Collective on PC
142:
143: Input Parameters:
144: + pc - the preconditioner context
145: - flag - PETSC_TRUE to use single precision, default is PETSC_FALSE
147: Options Database Key:
148: . -pc_ilu_single_precision_solves - activates PCILUSetSinglePrecisionSolves
149: . -pc_ilu_single_precision_solves <flag> - flag of TRUE/FALSE will override this call in the code.
151: Level: intermediate
153: .keywords: PC, set, solve, precision
154: @*/
155: int PCILUSetSinglePrecisionSolves(PC pc,PetscTruth flag)
156: {
157: int ierr,(*f)(PC,PetscTruth);
161: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetSinglePrecisionSolves_C",(void (**)(void))&f);
162: if (f) {
163: (*f)(pc,flag);
164: }
165: return(0);
166: }
168: /*@
169: PCILUSetDamping - adds this quantity to the diagonal of the matrix during the
170: ILU numerical factorization
172: Collective on PC
173:
174: Input Parameters:
175: + pc - the preconditioner context
176: - damping - amount of damping
178: Options Database Key:
179: . -pc_ilu_damping <damping> - Sets damping amount
181: Level: intermediate
183: .keywords: PC, set, factorization, direct, fill
185: .seealso: PCILUSetFill(), PCLUSetDamp()
186: @*/
187: int PCILUSetDamping(PC pc,PetscReal damping)
188: {
189: int ierr,(*f)(PC,PetscReal);
193: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);
194: if (f) {
195: (*f)(pc,damping);
196: }
197: return(0);
198: }
200: /*@
201: PCILUSetUseDropTolerance - The preconditioner will use an ILU
202: based on a drop tolerance.
204: Collective on PC
206: Input Parameters:
207: + pc - the preconditioner context
208: . dt - the drop tolerance, try from 1.e-10 to .1
209: . dtcol - tolerance for column pivot, good values [0.1 to 0.01]
210: - maxrowcount - the max number of nonzeros allowed in a row, best value
211: depends on the number of nonzeros in row of original matrix
213: Options Database Key:
214: . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance
216: Level: intermediate
218: Notes:
219: This uses the iludt() code of Saad's SPARSKIT package
221: .keywords: PC, levels, reordering, factorization, incomplete, ILU
222: @*/
223: int PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,int maxrowcount)
224: {
225: int ierr,(*f)(PC,PetscReal,PetscReal,int);
229: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);
230: if (f) {
231: (*f)(pc,dt,dtcol,maxrowcount);
232: }
233: return(0);
234: }
236: /*@
237: PCILUSetFill - Indicate the amount of fill you expect in the factored matrix,
238: where fill = number nonzeros in factor/number nonzeros in original matrix.
240: Collective on PC
242: Input Parameters:
243: + pc - the preconditioner context
244: - fill - amount of expected fill
246: Options Database Key:
247: $ -pc_ilu_fill <fill>
249: Note:
250: For sparse matrix factorizations it is difficult to predict how much
251: fill to expect. By running with the option -log_info PETSc will print the
252: actual amount of fill used; allowing you to set the value accurately for
253: future runs. But default PETSc uses a value of 1.0
255: Level: intermediate
257: .keywords: PC, set, factorization, direct, fill
259: .seealso: PCLUSetFill()
260: @*/
261: int PCILUSetFill(PC pc,PetscReal fill)
262: {
263: int ierr,(*f)(PC,PetscReal);
267: if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
268: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);
269: if (f) {
270: (*f)(pc,fill);
271: }
272: return(0);
273: }
275: /*@C
276: PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to
277: be used in the ILU factorization.
279: Collective on PC
281: Input Parameters:
282: + pc - the preconditioner context
283: - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
285: Options Database Key:
286: . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine
288: Level: intermediate
290: Notes: natural ordering is used by default
292: .seealso: PCLUSetMatOrdering()
294: .keywords: PC, ILU, set, matrix, reordering
296: @*/
297: int PCILUSetMatOrdering(PC pc,MatOrderingType ordering)
298: {
299: int ierr,(*f)(PC,MatOrderingType);
303: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);
304: if (f) {
305: (*f)(pc,ordering);
306: }
307: return(0);
308: }
310: /*@
311: PCILUSetReuseOrdering - When similar matrices are factored, this
312: causes the ordering computed in the first factor to be used for all
313: following factors; applies to both fill and drop tolerance ILUs.
315: Collective on PC
317: Input Parameters:
318: + pc - the preconditioner context
319: - flag - PETSC_TRUE to reuse else PETSC_FALSE
321: Options Database Key:
322: . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering()
324: Level: intermediate
326: .keywords: PC, levels, reordering, factorization, incomplete, ILU
328: .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
329: @*/
330: int PCILUSetReuseOrdering(PC pc,PetscTruth flag)
331: {
332: int ierr,(*f)(PC,PetscTruth);
336: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);
337: if (f) {
338: (*f)(pc,flag);
339: }
340: return(0);
341: }
343: /*@
344: PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored,
345: this causes later ones to use the fill computed in the initial factorization.
347: Collective on PC
349: Input Parameters:
350: + pc - the preconditioner context
351: - flag - PETSC_TRUE to reuse else PETSC_FALSE
353: Options Database Key:
354: . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill()
356: Level: intermediate
358: .keywords: PC, levels, reordering, factorization, incomplete, ILU
360: .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill()
361: @*/
362: int PCILUDTSetReuseFill(PC pc,PetscTruth flag)
363: {
364: int ierr,(*f)(PC,PetscTruth);
368: PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);
369: if (f) {
370: (*f)(pc,flag);
371: }
372: return(0);
373: }
375: /*@
376: PCILUSetLevels - Sets the number of levels of fill to use.
378: Collective on PC
380: Input Parameters:
381: + pc - the preconditioner context
382: - levels - number of levels of fill
384: Options Database Key:
385: . -pc_ilu_levels <levels> - Sets fill level
387: Level: intermediate
389: .keywords: PC, levels, fill, factorization, incomplete, ILU
390: @*/
391: int PCILUSetLevels(PC pc,int levels)
392: {
393: int ierr,(*f)(PC,int);
397: if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
398: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);
399: if (f) {
400: (*f)(pc,levels);
401: }
402: return(0);
403: }
405: /*@
406: PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be
407: treated as level 0 fill even if there is no non-zero location.
409: Collective on PC
411: Input Parameters:
412: + pc - the preconditioner context
414: Options Database Key:
415: . -pc_ilu_diagonal_fill
417: Notes:
418: Does not apply with 0 fill.
420: Level: intermediate
422: .keywords: PC, levels, fill, factorization, incomplete, ILU
423: @*/
424: int PCILUSetAllowDiagonalFill(PC pc)
425: {
426: int ierr,(*f)(PC);
430: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);
431: if (f) {
432: (*f)(pc);
433: }
434: return(0);
435: }
437: /*@
438: PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization.
439: Collective on PC
441: Input Parameters:
442: . pc - the preconditioner context
444: Options Database Key:
445: . -pc_ilu_in_place - Activates in-place factorization
447: Notes:
448: PCILUSetUseInPlace() is intended for use with matrix-free variants of
449: Krylov methods, or when a different matrices are employed for the linear
450: system and preconditioner, or with ASM preconditioning. Do NOT use
451: this option if the linear system
452: matrix also serves as the preconditioning matrix, since the factored
453: matrix would then overwrite the original matrix.
455: Only works well with ILU(0).
457: Level: intermediate
459: .keywords: PC, set, factorization, inplace, in-place, ILU
461: .seealso: PCLUSetUseInPlace()
462: @*/
463: int PCILUSetUseInPlace(PC pc)
464: {
465: int ierr,(*f)(PC);
469: PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);
470: if (f) {
471: (*f)(pc);
472: }
473: return(0);
474: }
476: /* ------------------------------------------------------------------------------------------*/
478: static int PCSetFromOptions_ILU(PC pc)
479: {
480: int ierr,dtmax = 3,itmp;
481: PetscTruth flg,single_prec;
482: PetscReal dt[3];
483: char tname[256];
484: PC_ILU *ilu = (PC_ILU*)pc->data;
485: PetscFList ordlist;
486: PetscReal tol;
489: MatOrderingRegisterAll(PETSC_NULL);
490: PetscOptionsHead("ILU Options");
491: PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(int)ilu->info.levels,&itmp,&flg);
492: if (flg) ilu->info.levels = (double) itmp;
493: PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);
494: PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);
495: ilu->info.diagonal_fill = (double) flg;
496: PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);
497: PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);
498: PetscOptionsLogical("-pc_ilu_single_precision_solves","Precision for triangular solves","PCILUSetSinglePrecisionSolves",ilu->single_precision_solve,&single_prec,&flg);
499: if (flg) {
500: PCILUSetSinglePrecisionSolves(pc,single_prec);
501: }
502: PetscOptionsHasName(pc->prefix,"-pc_ilu_damping",&flg);
503: if (flg) {
504: PCILUSetDamping(pc,0.0);
505: }
506: PetscOptionsReal("-pc_ilu_damping","Damping added to diagonal","PCILUSetDamping",ilu->info.damping,&ilu->info.damping,0);
507: PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);
509: dt[0] = ilu->info.dt;
510: dt[1] = ilu->info.dtcol;
511: dt[2] = ilu->info.dtcount;
512: PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);
513: if (flg) {
514: PCILUSetUseDropTolerance(pc,dt[0],dt[1],(int)dt[2]);
515: }
516: PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);
517: PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","MatReorderForNonzeroDiagonal",0.0,&tol,0);
519: MatGetOrderingList(&ordlist);
520: PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);
521: if (flg) {
522: PCILUSetMatOrdering(pc,tname);
523: }
524: PetscOptionsTail();
525: return(0);
526: }
528: static int PCView_ILU(PC pc,PetscViewer viewer)
529: {
530: PC_ILU *ilu = (PC_ILU*)pc->data;
531: int ierr;
532: PetscTruth isstring,isascii;
535: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
536: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
537: if (isascii) {
538: if (ilu->usedt) {
539: PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %gn",ilu->info.dt);
540: PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %dn",(int)ilu->info.dtcount);
541: PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %gn",ilu->info.dtcol);
542: } else if (ilu->info.levels == 1) {
543: PetscViewerASCIIPrintf(viewer," ILU: %d level of filln",(int)ilu->info.levels);
544: } else {
545: PetscViewerASCIIPrintf(viewer," ILU: %d levels of filln",(int)ilu->info.levels);
546: }
547: PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %gn",ilu->info.fill);
548: PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %gn",ilu->info.zeropivot);
549: if (ilu->inplace) {PetscViewerASCIIPrintf(viewer," in-place factorizationn");}
550: else {PetscViewerASCIIPrintf(viewer," out-of-place factorizationn");}
551: PetscViewerASCIIPrintf(viewer," matrix ordering: %sn",ilu->ordering);
552: if (ilu->reusefill) {PetscViewerASCIIPrintf(viewer," Reusing fill from past factorizationn");}
553: if (ilu->reuseordering) {PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorizationn");}
554: if (ilu->fact) {
555: PetscViewerASCIIPrintf(viewer," Factored matrix followsn");
556: PetscViewerASCIIPushTab(viewer);
557: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
558: MatView(ilu->fact,viewer);
559: PetscViewerPopFormat(viewer);
560: PetscViewerASCIIPopTab(viewer);
561: }
562: } else if (isstring) {
563: PetscViewerStringSPrintf(viewer," lvls=%g,order=%s",ilu->info.levels,ilu->ordering);
564: } else {
565: SETERRQ1(1,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name);
566: }
567: return(0);
568: }
570: static int PCSetUp_ILU(PC pc)
571: {
572: int ierr;
573: PetscTruth flg;
574: PC_ILU *ilu = (PC_ILU*)pc->data;
577: if (ilu->inplace) {
578: if (!pc->setupcalled) {
580: /* In-place factorization only makes sense with the natural ordering,
581: so we only need to get the ordering once, even if nonzero structure changes */
582: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
583: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
584: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
585: }
587: /* In place ILU only makes sense with fill factor of 1.0 because
588: cannot have levels of fill */
589: ilu->info.fill = 1.0;
590: ilu->info.diagonal_fill = 0;
591: MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);
592: ilu->fact = pc->pmat;
593: } else if (ilu->usedt) {
594: if (!pc->setupcalled) {
595: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
596: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
597: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
598: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
599: PetscLogObjectParent(pc,ilu->fact);
600: } else if (pc->flag != SAME_NONZERO_PATTERN) {
601: MatDestroy(ilu->fact);
602: if (!ilu->reuseordering) {
603: if (ilu->row) {ISDestroy(ilu->row);}
604: if (ilu->col) {ISDestroy(ilu->col);}
605: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
606: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
607: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
608: }
609: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
610: PetscLogObjectParent(pc,ilu->fact);
611: } else if (!ilu->reusefill) {
612: MatDestroy(ilu->fact);
613: MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);
614: PetscLogObjectParent(pc,ilu->fact);
615: } else {
616: MatLUFactorNumeric(pc->pmat,&ilu->fact);
617: }
618: } else {
619: if (!pc->setupcalled) {
620: /* first time in so compute reordering and symbolic factorization */
621: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
622: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
623: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
624: /* Remove zeros along diagonal? */
625: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
626: if (flg) {
627: PetscReal ntol = 1.e-10;
628: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
629: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
630: }
632: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
633: PetscLogObjectParent(pc,ilu->fact);
634: } else if (pc->flag != SAME_NONZERO_PATTERN) {
635: if (!ilu->reuseordering) {
636: /* compute a new ordering for the ILU */
637: ISDestroy(ilu->row);
638: ISDestroy(ilu->col);
639: MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);
640: if (ilu->row) PetscLogObjectParent(pc,ilu->row);
641: if (ilu->col) PetscLogObjectParent(pc,ilu->col);
642: /* Remove zeros along diagonal? */
643: PetscOptionsHasName(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&flg);
644: if (flg) {
645: PetscReal ntol = 1.e-10;
646: PetscOptionsGetReal(pc->prefix,"-pc_ilu_nonzeros_along_diagonal",&ntol,PETSC_NULL);
647: MatReorderForNonzeroDiagonal(pc->pmat,ntol,ilu->row,ilu->col);
648: }
649: }
650: MatDestroy(ilu->fact);
651: MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);
652: PetscLogObjectParent(pc,ilu->fact);
653: }
654: MatLUFactorNumeric(pc->pmat,&ilu->fact);
655: }
656: if (ilu->single_precision_solve) {
657: MatSetOption(ilu->fact,MAT_USE_SINGLE_PRECISION_SOLVES);
658: }
659: return(0);
660: }
662: static int PCDestroy_ILU(PC pc)
663: {
664: PC_ILU *ilu = (PC_ILU*)pc->data;
665: int ierr;
668: if (!ilu->inplace && ilu->fact) {MatDestroy(ilu->fact);}
669: if (ilu->row && ilu->col && ilu->row != ilu->col) {ISDestroy(ilu->row);}
670: if (ilu->col) {ISDestroy(ilu->col);}
671: PetscStrfree(ilu->ordering);
672: PetscFree(ilu);
673: return(0);
674: }
676: static int PCApply_ILU(PC pc,Vec x,Vec y)
677: {
678: PC_ILU *ilu = (PC_ILU*)pc->data;
679: int ierr;
682: MatSolve(ilu->fact,x,y);
683: return(0);
684: }
686: static int PCApplyTranspose_ILU(PC pc,Vec x,Vec y)
687: {
688: PC_ILU *ilu = (PC_ILU*)pc->data;
689: int ierr;
692: MatSolveTranspose(ilu->fact,x,y);
693: return(0);
694: }
696: static int PCGetFactoredMatrix_ILU(PC pc,Mat *mat)
697: {
698: PC_ILU *ilu = (PC_ILU*)pc->data;
701: if (!ilu->fact) SETERRQ(1,"Matrix not yet factored; call after SLESSetUp() or PCSetUp()");
702: *mat = ilu->fact;
703: return(0);
704: }
706: EXTERN_C_BEGIN
707: int PCCreate_ILU(PC pc)
708: {
709: int ierr;
710: PC_ILU *ilu;
713: PetscNew(PC_ILU,&ilu);
714: PetscLogObjectMemory(pc,sizeof(PC_ILU));
716: ilu->fact = 0;
717: ilu->info.levels = 0;
718: ilu->info.fill = 1.0;
719: ilu->col = 0;
720: ilu->row = 0;
721: ilu->inplace = PETSC_FALSE;
722: PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);
723: ilu->reuseordering = PETSC_FALSE;
724: ilu->usedt = PETSC_FALSE;
725: ilu->info.dt = PETSC_DEFAULT;
726: ilu->info.dtcount = PETSC_DEFAULT;
727: ilu->info.dtcol = PETSC_DEFAULT;
728: ilu->info.damp = 0.0;
729: ilu->info.damping = 0.0;
730: ilu->info.zeropivot = 1.e-12;
731: ilu->reusefill = PETSC_FALSE;
732: ilu->info.diagonal_fill = 0;
733: ilu->single_precision_solve = PETSC_FALSE;
734: pc->data = (void*)ilu;
736: pc->ops->destroy = PCDestroy_ILU;
737: pc->ops->apply = PCApply_ILU;
738: pc->ops->applytranspose = PCApplyTranspose_ILU;
739: pc->ops->setup = PCSetUp_ILU;
740: pc->ops->setfromoptions = PCSetFromOptions_ILU;
741: pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU;
742: pc->ops->view = PCView_ILU;
743: pc->ops->applyrichardson = 0;
745: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetSinglePrecisionSolves_C",
746: "PCILUSetSinglePrecisionSolves_ILU",
747: PCILUSetSinglePrecisionSolves_ILU);
748: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU",
749: PCILUSetUseDropTolerance_ILU);
750: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU",
751: PCILUSetFill_ILU);
752: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU",
753: PCILUSetDamping_ILU);
754: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU",
755: PCILUSetMatOrdering_ILU);
756: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU",
757: PCILUSetReuseOrdering_ILU);
758: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT",
759: PCILUDTSetReuseFill_ILUDT);
760: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU",
761: PCILUSetLevels_ILU);
762: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU",
763: PCILUSetUseInPlace_ILU);
764: PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU",
765: PCILUSetAllowDiagonalFill_ILU);
767: return(0);
768: }
769: EXTERN_C_END