xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 075768bc809001225c61399c80901bb8d8a18b20)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
385b398ccSKris Buschelman #include "src/ksp/pc/impls/factor/icc/icc.h"   /*I "petscpc.h" I*/
49b54502bSHong Zhang 
59b54502bSHong Zhang EXTERN_C_BEGIN
69b54502bSHong Zhang #undef __FUNCT__
7afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ICC"
8dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ICC(PC pc,PetscReal z)
9afaefe49SHong Zhang {
10*075768bcSBarry Smith   PC_ICC *icc  = (PC_ICC*)pc->data;
11afaefe49SHong Zhang 
12afaefe49SHong Zhang   PetscFunctionBegin;
13*075768bcSBarry Smith   ((PC_Factor*)icc)->info.zeropivot = z;
14afaefe49SHong Zhang   PetscFunctionReturn(0);
15afaefe49SHong Zhang }
16afaefe49SHong Zhang EXTERN_C_END
17afaefe49SHong Zhang 
18afaefe49SHong Zhang EXTERN_C_BEGIN
19afaefe49SHong Zhang #undef __FUNCT__
20afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ICC"
21dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ICC(PC pc,PetscReal shift)
22afaefe49SHong Zhang {
23*075768bcSBarry Smith   PC_ICC *dir  = (PC_ICC*)pc->data;
24afaefe49SHong Zhang 
25afaefe49SHong Zhang   PetscFunctionBegin;
26afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
27*075768bcSBarry Smith      ((PC_Factor*)dir)->info.shiftnz = 1.e-12;
28afaefe49SHong Zhang   } else {
29*075768bcSBarry Smith      ((PC_Factor*)dir)->info.shiftnz = shift;
30afaefe49SHong Zhang   }
31afaefe49SHong Zhang   PetscFunctionReturn(0);
32afaefe49SHong Zhang }
33afaefe49SHong Zhang EXTERN_C_END
34afaefe49SHong Zhang 
35afaefe49SHong Zhang EXTERN_C_BEGIN
36afaefe49SHong Zhang #undef __FUNCT__
37afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ICC"
38dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ICC(PC pc,PetscTruth shift)
39afaefe49SHong Zhang {
40*075768bcSBarry Smith   PC_ICC *dir = (PC_ICC*)pc->data;
41afaefe49SHong Zhang 
42afaefe49SHong Zhang   PetscFunctionBegin;
43fbf22428SSatish Balay   if (shift) {
44*075768bcSBarry Smith      ((PC_Factor*)dir)->info.shiftpd = 1.0;
45fbf22428SSatish Balay   } else {
46*075768bcSBarry Smith      ((PC_Factor*)dir)->info.shiftpd = 0.0;
47fbf22428SSatish Balay   }
48afaefe49SHong Zhang   PetscFunctionReturn(0);
49afaefe49SHong Zhang }
50afaefe49SHong Zhang EXTERN_C_END
51afaefe49SHong Zhang 
52afaefe49SHong Zhang EXTERN_C_BEGIN
53afaefe49SHong Zhang #undef __FUNCT__
54e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ICC"
55e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ICC(PC pc,const MatOrderingType ordering)
569b54502bSHong Zhang {
579b54502bSHong Zhang   PC_ICC         *dir = (PC_ICC*)pc->data;
589b54502bSHong Zhang   PetscErrorCode ierr;
599b54502bSHong Zhang 
609b54502bSHong Zhang   PetscFunctionBegin;
61*075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
62*075768bcSBarry Smith   ierr = PetscStrallocpy(ordering,& ((PC_Factor*)dir)->ordering);CHKERRQ(ierr);
639b54502bSHong Zhang   PetscFunctionReturn(0);
649b54502bSHong Zhang }
659b54502bSHong Zhang EXTERN_C_END
669b54502bSHong Zhang 
679b54502bSHong Zhang EXTERN_C_BEGIN
689b54502bSHong Zhang #undef __FUNCT__
6955ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ICC"
7055ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ICC(PC pc,PetscReal fill)
719b54502bSHong Zhang {
72*075768bcSBarry Smith   PC_ICC *dir = (PC_ICC*)pc->data;
739b54502bSHong Zhang 
749b54502bSHong Zhang   PetscFunctionBegin;
75*075768bcSBarry Smith   ((PC_Factor*)dir)->info.fill = fill;
769b54502bSHong Zhang   PetscFunctionReturn(0);
779b54502bSHong Zhang }
789b54502bSHong Zhang EXTERN_C_END
799b54502bSHong Zhang 
809b54502bSHong Zhang EXTERN_C_BEGIN
819b54502bSHong Zhang #undef __FUNCT__
822401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels_ICC"
832401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ICC(PC pc,PetscInt levels)
849b54502bSHong Zhang {
85*075768bcSBarry Smith   PC_ICC *icc = (PC_ICC*)pc->data;
869b54502bSHong Zhang 
879b54502bSHong Zhang   PetscFunctionBegin;
88*075768bcSBarry Smith    ((PC_Factor*)icc)->info.levels = levels;
899b54502bSHong Zhang   PetscFunctionReturn(0);
909b54502bSHong Zhang }
919b54502bSHong Zhang EXTERN_C_END
929b54502bSHong Zhang 
939b54502bSHong Zhang #undef __FUNCT__
949b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC"
959b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc)
969b54502bSHong Zhang {
979b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
989b54502bSHong Zhang   IS             perm,cperm;
999b54502bSHong Zhang   PetscErrorCode ierr;
100f3a39becSBarry Smith   MatInfo        info;
1019b54502bSHong Zhang 
1029b54502bSHong Zhang   PetscFunctionBegin;
103*075768bcSBarry Smith   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
1049b54502bSHong Zhang 
1059b54502bSHong Zhang   if (!pc->setupcalled) {
106*075768bcSBarry Smith     ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr);
107*075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
1089b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
109*075768bcSBarry Smith     ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);
110*075768bcSBarry Smith     ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
111*075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
1129b54502bSHong Zhang   }
113*075768bcSBarry Smith   ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
114f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
115f3a39becSBarry Smith 
1169b54502bSHong Zhang   ierr = ISDestroy(cperm);CHKERRQ(ierr);
1179b54502bSHong Zhang   ierr = ISDestroy(perm);CHKERRQ(ierr);
118*075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
1199b54502bSHong Zhang   PetscFunctionReturn(0);
1209b54502bSHong Zhang }
1219b54502bSHong Zhang 
1229b54502bSHong Zhang #undef __FUNCT__
1239b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
1249b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
1259b54502bSHong Zhang {
1269b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1279b54502bSHong Zhang   PetscErrorCode ierr;
1289b54502bSHong Zhang 
1299b54502bSHong Zhang   PetscFunctionBegin;
130*075768bcSBarry Smith   if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);}
131*075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
1329b54502bSHong Zhang   ierr = PetscFree(icc);CHKERRQ(ierr);
1339b54502bSHong Zhang   PetscFunctionReturn(0);
1349b54502bSHong Zhang }
1359b54502bSHong Zhang 
1369b54502bSHong Zhang #undef __FUNCT__
1379b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
1389b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
1399b54502bSHong Zhang {
1409b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1419b54502bSHong Zhang   PetscErrorCode ierr;
1429b54502bSHong Zhang 
1439b54502bSHong Zhang   PetscFunctionBegin;
144*075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1459b54502bSHong Zhang   PetscFunctionReturn(0);
1469b54502bSHong Zhang }
1479b54502bSHong Zhang 
1489b54502bSHong Zhang #undef __FUNCT__
1499b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
1509b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
1519b54502bSHong Zhang {
1529b54502bSHong Zhang   PetscErrorCode ierr;
1539b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1549b54502bSHong Zhang 
1559b54502bSHong Zhang   PetscFunctionBegin;
156*075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1579b54502bSHong Zhang   PetscFunctionReturn(0);
1589b54502bSHong Zhang }
1599b54502bSHong Zhang 
1609b54502bSHong Zhang #undef __FUNCT__
1619b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
1629b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1639b54502bSHong Zhang {
1649b54502bSHong Zhang   PetscErrorCode ierr;
1659b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1669b54502bSHong Zhang 
1679b54502bSHong Zhang   PetscFunctionBegin;
168*075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1699b54502bSHong Zhang   PetscFunctionReturn(0);
1709b54502bSHong Zhang }
1719b54502bSHong Zhang 
1729b54502bSHong Zhang #undef __FUNCT__
173a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_ICC"
174a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_ICC(PC pc,Mat *mat)
1759b54502bSHong Zhang {
1769b54502bSHong Zhang   PC_ICC *icc = (PC_ICC*)pc->data;
1779b54502bSHong Zhang 
1789b54502bSHong Zhang   PetscFunctionBegin;
179*075768bcSBarry Smith   *mat = ((PC_Factor*)icc)->fact;
1809b54502bSHong Zhang   PetscFunctionReturn(0);
1819b54502bSHong Zhang }
1829b54502bSHong Zhang 
1839b54502bSHong Zhang #undef __FUNCT__
1849b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1859b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc)
1869b54502bSHong Zhang {
1879b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1889b54502bSHong Zhang   char           tname[256];
1899b54502bSHong Zhang   PetscTruth     flg;
1909b54502bSHong Zhang   PetscErrorCode ierr;
1919b54502bSHong Zhang   PetscFList     ordlist;
1929b54502bSHong Zhang 
1939b54502bSHong Zhang   PetscFunctionBegin;
1949b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
1959b54502bSHong Zhang   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
196*075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
197*075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)icc)->info.fill,&((PC_Factor*)icc)->info.fill,&flg);CHKERRQ(ierr);
1989b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
199*075768bcSBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)icc)->ordering,tname,256,&flg);CHKERRQ(ierr);
2009b54502bSHong Zhang     if (flg) {
201e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
2029b54502bSHong Zhang     }
2039f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
2049b54502bSHong Zhang     if (flg) {
205afaefe49SHong Zhang       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
2069b54502bSHong Zhang     }
207*075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)icc)->info.shiftnz,&((PC_Factor*)icc)->info.shiftnz,0);CHKERRQ(ierr);
208*075768bcSBarry Smith     flg = (((PC_Factor*)icc)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE;
20962bba022SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
21062bba022SBarry Smith     ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr);
211*075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)icc)->info.zeropivot,&((PC_Factor*)icc)->info.zeropivot,0);CHKERRQ(ierr);
2129b54502bSHong Zhang 
2139b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
2149b54502bSHong Zhang   PetscFunctionReturn(0);
2159b54502bSHong Zhang }
2169b54502bSHong Zhang 
2179b54502bSHong Zhang #undef __FUNCT__
2189b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
2199b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
2209b54502bSHong Zhang {
2219b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2229b54502bSHong Zhang   PetscErrorCode ierr;
2239b54502bSHong Zhang   PetscTruth     isstring,iascii;
2249b54502bSHong Zhang 
2259b54502bSHong Zhang   PetscFunctionBegin;
2269b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
2279b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
2289b54502bSHong Zhang   if (iascii) {
229*075768bcSBarry Smith     if (((PC_Factor*)icc)->info.levels == 1) {
230*075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D level of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
2319b54502bSHong Zhang     } else {
232*075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D levels of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
2339b54502bSHong Zhang     }
234*075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio allocated %G, ordering used %s\n",((PC_Factor*)icc)->info.fill,((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
235*075768bcSBarry Smith     if (((PC_Factor*)icc)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using Manteuffel shift\n");CHKERRQ(ierr);}
236*075768bcSBarry Smith     if (((PC_Factor*)icc)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
237*075768bcSBarry Smith     if (((PC_Factor*)icc)->fact) {
238f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr);
239f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
240f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
241f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
242f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
243f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
244*075768bcSBarry Smith       ierr = MatView(((PC_Factor*)icc)->fact,viewer);CHKERRQ(ierr);
245f3a39becSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
246f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
247f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
248f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
249f3a39becSBarry Smith     }
2509b54502bSHong Zhang   } else if (isstring) {
251*075768bcSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);CHKERRQ(ierr);
2529b54502bSHong Zhang   } else {
2539b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name);
2549b54502bSHong Zhang   }
2559b54502bSHong Zhang   PetscFunctionReturn(0);
2569b54502bSHong Zhang }
2579b54502bSHong Zhang 
2589b54502bSHong Zhang /*MC
2599b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
2609b54502bSHong Zhang 
2619b54502bSHong Zhang    Options Database Keys:
2622401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
2632401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
2649b54502bSHong Zhang                       its factorization (overwrites original matrix)
26555ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
2662401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
267f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
268f251bdbdSHong Zhang -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
269f251bdbdSHong Zhang    is optional with PETSC_TRUE being the default
2709b54502bSHong Zhang 
2719b54502bSHong Zhang    Level: beginner
2729b54502bSHong Zhang 
2739b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
2749b54502bSHong Zhang 
275c582cd25SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel (for parallel use you
276c582cd25SBarry Smith              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ICC(0) and this is not recommended
277c582cd25SBarry Smith              unless you really want a parallel ICC).
2789b54502bSHong Zhang 
2799b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
2809b54502bSHong Zhang 
2819b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
2829b54502bSHong Zhang 
2832401956bSBarry Smith           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE);
2849b54502bSHong Zhang           to turn off the shift.
2859b54502bSHong Zhang 
286c582cd25SBarry Smith    References:
287c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
288c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
289c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
290c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
291c582cd25SBarry Smith 
2929b54502bSHong Zhang 
2939b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
294ee45ca4aSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
295e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
2962401956bSBarry Smith            PCFactorSetLevels(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(),
2979b54502bSHong Zhang 
2989b54502bSHong Zhang M*/
2999b54502bSHong Zhang 
3009b54502bSHong Zhang EXTERN_C_BEGIN
3019b54502bSHong Zhang #undef __FUNCT__
3029b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
303dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc)
3049b54502bSHong Zhang {
3059b54502bSHong Zhang   PetscErrorCode ierr;
3069b54502bSHong Zhang   PC_ICC         *icc;
3079b54502bSHong Zhang 
3089b54502bSHong Zhang   PetscFunctionBegin;
30938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr);
3109b54502bSHong Zhang 
311*075768bcSBarry Smith   ((PC_Factor*)icc)->fact	          = 0;
312*075768bcSBarry Smith   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
313*075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
314*075768bcSBarry Smith   ((PC_Factor*)icc)->info.levels	  = 0;
315*075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill          = 1.0;
3169b54502bSHong Zhang   icc->implctx            = 0;
3179b54502bSHong Zhang 
318*075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol              = PETSC_DEFAULT;
319*075768bcSBarry Smith   ((PC_Factor*)icc)->info.shiftnz            = 0.0;
320*075768bcSBarry Smith   ((PC_Factor*)icc)->info.shiftpd            = 1.0; /* true */
321*075768bcSBarry Smith   ((PC_Factor*)icc)->info.zeropivot          = 1.e-12;
3229b54502bSHong Zhang   pc->data	               = (void*)icc;
3239b54502bSHong Zhang 
3249b54502bSHong Zhang   pc->ops->apply	       = PCApply_ICC;
3259b54502bSHong Zhang   pc->ops->setup               = PCSetup_ICC;
3269b54502bSHong Zhang   pc->ops->destroy	       = PCDestroy_ICC;
3279b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
3289b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
329a4fd02acSBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_ICC;
3309b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
3319b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
3329b54502bSHong Zhang 
333afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ICC",
334afaefe49SHong Zhang                     PCFactorSetZeroPivot_ICC);CHKERRQ(ierr);
335afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ICC",
336afaefe49SHong Zhang                     PCFactorSetShiftNonzero_ICC);CHKERRQ(ierr);
337afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ICC",
338afaefe49SHong Zhang                     PCFactorSetShiftPd_ICC);CHKERRQ(ierr);
339afaefe49SHong Zhang 
3402401956bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ICC",
3412401956bSBarry Smith                     PCFactorSetLevels_ICC);CHKERRQ(ierr);
34255ba2a51SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ICC",
34355ba2a51SBarry Smith                     PCFactorSetFill_ICC);CHKERRQ(ierr);
344e5a9bf91SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ICC",
345e5a9bf91SBarry Smith                     PCFactorSetMatOrderingType_ICC);CHKERRQ(ierr);
3469b54502bSHong Zhang   PetscFunctionReturn(0);
3479b54502bSHong Zhang }
3489b54502bSHong Zhang EXTERN_C_END
3499b54502bSHong Zhang 
3509b54502bSHong Zhang 
351