xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 52e6d16ba989af9362d0fcebb12e73dd7c9666ed)
19b54502bSHong Zhang /*
29b54502bSHong Zhang    Defines a Cholesky factorization preconditioner for any Mat implementation.
39b54502bSHong Zhang   Presently only provided for MPIRowbs format (i.e. BlockSolve).
49b54502bSHong Zhang */
59b54502bSHong Zhang 
685b398ccSKris Buschelman #include "src/ksp/pc/impls/factor/icc/icc.h"   /*I "petscpc.h" I*/
79b54502bSHong Zhang 
89b54502bSHong Zhang EXTERN_C_BEGIN
99b54502bSHong Zhang #undef __FUNCT__
10afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ICC"
11afaefe49SHong Zhang PetscErrorCode PCFactorSetZeroPivot_ICC(PC pc,PetscReal z)
12afaefe49SHong Zhang {
13afaefe49SHong Zhang   PC_ICC *icc;
14afaefe49SHong Zhang 
15afaefe49SHong Zhang   PetscFunctionBegin;
16afaefe49SHong Zhang   icc                 = (PC_ICC*)pc->data;
17afaefe49SHong Zhang   icc->info.zeropivot = z;
18afaefe49SHong Zhang   PetscFunctionReturn(0);
19afaefe49SHong Zhang }
20afaefe49SHong Zhang EXTERN_C_END
21afaefe49SHong Zhang 
22afaefe49SHong Zhang EXTERN_C_BEGIN
23afaefe49SHong Zhang #undef __FUNCT__
24afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ICC"
25afaefe49SHong Zhang PetscErrorCode PCFactorSetShiftNonzero_ICC(PC pc,PetscReal shift)
26afaefe49SHong Zhang {
27afaefe49SHong Zhang   PC_ICC *dir;
28afaefe49SHong Zhang 
29afaefe49SHong Zhang   PetscFunctionBegin;
30afaefe49SHong Zhang   dir = (PC_ICC*)pc->data;
31afaefe49SHong Zhang   if (shift == (PetscReal) PETSC_DECIDE) {
32afaefe49SHong Zhang     dir->info.shiftnz = 1.e-12;
33afaefe49SHong Zhang   } else {
34afaefe49SHong Zhang     dir->info.shiftnz = shift;
35afaefe49SHong Zhang   }
36afaefe49SHong Zhang   PetscFunctionReturn(0);
37afaefe49SHong Zhang }
38afaefe49SHong Zhang EXTERN_C_END
39afaefe49SHong Zhang 
40afaefe49SHong Zhang EXTERN_C_BEGIN
41afaefe49SHong Zhang #undef __FUNCT__
42afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ICC"
43afaefe49SHong Zhang PetscErrorCode PCFactorSetShiftPd_ICC(PC pc,PetscTruth shift)
44afaefe49SHong Zhang {
45afaefe49SHong Zhang   PC_ICC *dir;
46afaefe49SHong Zhang 
47afaefe49SHong Zhang   PetscFunctionBegin;
48afaefe49SHong Zhang   dir = (PC_ICC*)pc->data;
49afaefe49SHong Zhang   dir->info.shiftpd = shift;
50afaefe49SHong Zhang   if (shift) dir->info.shift_fraction = 0.0;
51afaefe49SHong Zhang   PetscFunctionReturn(0);
52afaefe49SHong Zhang }
53afaefe49SHong Zhang EXTERN_C_END
54afaefe49SHong Zhang 
55afaefe49SHong Zhang EXTERN_C_BEGIN
56afaefe49SHong Zhang #undef __FUNCT__
579b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering_ICC"
589b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering)
599b54502bSHong Zhang {
609b54502bSHong Zhang   PC_ICC         *dir = (PC_ICC*)pc->data;
619b54502bSHong Zhang   PetscErrorCode ierr;
629b54502bSHong Zhang 
639b54502bSHong Zhang   PetscFunctionBegin;
649b54502bSHong Zhang   ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr);
659b54502bSHong Zhang   ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr);
669b54502bSHong Zhang   PetscFunctionReturn(0);
679b54502bSHong Zhang }
689b54502bSHong Zhang EXTERN_C_END
699b54502bSHong Zhang 
709b54502bSHong Zhang EXTERN_C_BEGIN
719b54502bSHong Zhang #undef __FUNCT__
729b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill_ICC"
739b54502bSHong Zhang PetscErrorCode PCICCSetFill_ICC(PC pc,PetscReal fill)
749b54502bSHong Zhang {
759b54502bSHong Zhang   PC_ICC *dir;
769b54502bSHong Zhang 
779b54502bSHong Zhang   PetscFunctionBegin;
789b54502bSHong Zhang   dir            = (PC_ICC*)pc->data;
799b54502bSHong Zhang   dir->info.fill = fill;
809b54502bSHong Zhang   PetscFunctionReturn(0);
819b54502bSHong Zhang }
829b54502bSHong Zhang EXTERN_C_END
839b54502bSHong Zhang 
849b54502bSHong Zhang EXTERN_C_BEGIN
859b54502bSHong Zhang #undef __FUNCT__
869b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels_ICC"
879b54502bSHong Zhang PetscErrorCode PCICCSetLevels_ICC(PC pc,PetscInt levels)
889b54502bSHong Zhang {
899b54502bSHong Zhang   PC_ICC *icc;
909b54502bSHong Zhang 
919b54502bSHong Zhang   PetscFunctionBegin;
929b54502bSHong Zhang   icc = (PC_ICC*)pc->data;
939b54502bSHong Zhang   icc->info.levels = levels;
949b54502bSHong Zhang   PetscFunctionReturn(0);
959b54502bSHong Zhang }
969b54502bSHong Zhang EXTERN_C_END
979b54502bSHong Zhang 
989b54502bSHong Zhang #undef __FUNCT__
999b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering"
1009b54502bSHong Zhang /*@
1019b54502bSHong Zhang     PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to
1029b54502bSHong Zhang     be used it the ICC factorization.
1039b54502bSHong Zhang 
1049b54502bSHong Zhang     Collective on PC
1059b54502bSHong Zhang 
1069b54502bSHong Zhang     Input Parameters:
1079b54502bSHong Zhang +   pc - the preconditioner context
1089b54502bSHong Zhang -   ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM
1099b54502bSHong Zhang 
1109b54502bSHong Zhang     Options Database Key:
1119b54502bSHong Zhang .   -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine
1129b54502bSHong Zhang 
1139b54502bSHong Zhang     Level: intermediate
1149b54502bSHong Zhang 
1159b54502bSHong Zhang .seealso: PCLUSetMatOrdering()
1169b54502bSHong Zhang 
1179b54502bSHong Zhang .keywords: PC, ICC, set, matrix, reordering
1189b54502bSHong Zhang 
1199b54502bSHong Zhang @*/
1209b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering(PC pc,MatOrderingType ordering)
1219b54502bSHong Zhang {
1229b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,MatOrderingType);
1239b54502bSHong Zhang 
1249b54502bSHong Zhang   PetscFunctionBegin;
1259b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1269b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr);
1279b54502bSHong Zhang   if (f) {
1289b54502bSHong Zhang     ierr = (*f)(pc,ordering);CHKERRQ(ierr);
1299b54502bSHong Zhang   }
1309b54502bSHong Zhang   PetscFunctionReturn(0);
1319b54502bSHong Zhang }
1329b54502bSHong Zhang 
1339b54502bSHong Zhang #undef __FUNCT__
1349b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels"
1359b54502bSHong Zhang /*@
1369b54502bSHong Zhang    PCICCSetLevels - Sets the number of levels of fill to use.
1379b54502bSHong Zhang 
1389b54502bSHong Zhang    Collective on PC
1399b54502bSHong Zhang 
1409b54502bSHong Zhang    Input Parameters:
1419b54502bSHong Zhang +  pc - the preconditioner context
1429b54502bSHong Zhang -  levels - number of levels of fill
1439b54502bSHong Zhang 
1449b54502bSHong Zhang    Options Database Key:
1459b54502bSHong Zhang .  -pc_icc_levels <levels> - Sets fill level
1469b54502bSHong Zhang 
1479b54502bSHong Zhang    Level: intermediate
1489b54502bSHong Zhang 
1499b54502bSHong Zhang    Concepts: ICC^setting levels of fill
1509b54502bSHong Zhang 
1519b54502bSHong Zhang @*/
1529b54502bSHong Zhang PetscErrorCode PCICCSetLevels(PC pc,PetscInt levels)
1539b54502bSHong Zhang {
1549b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscInt);
1559b54502bSHong Zhang 
1569b54502bSHong Zhang   PetscFunctionBegin;
1579b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1589b54502bSHong Zhang   if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels");
1599b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr);
1609b54502bSHong Zhang   if (f) {
1619b54502bSHong Zhang     ierr = (*f)(pc,levels);CHKERRQ(ierr);
1629b54502bSHong Zhang   }
1639b54502bSHong Zhang   PetscFunctionReturn(0);
1649b54502bSHong Zhang }
1659b54502bSHong Zhang 
1669b54502bSHong Zhang #undef __FUNCT__
1679b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill"
1689b54502bSHong Zhang /*@
1699b54502bSHong Zhang    PCICCSetFill - Indicate the amount of fill you expect in the factored matrix,
1709b54502bSHong Zhang    where fill = number nonzeros in factor/number nonzeros in original matrix.
1719b54502bSHong Zhang 
1729b54502bSHong Zhang    Collective on PC
1739b54502bSHong Zhang 
1749b54502bSHong Zhang    Input Parameters:
1759b54502bSHong Zhang +  pc - the preconditioner context
1769b54502bSHong Zhang -  fill - amount of expected fill
1779b54502bSHong Zhang 
1789b54502bSHong Zhang    Options Database Key:
1799b54502bSHong Zhang $  -pc_icc_fill <fill>
1809b54502bSHong Zhang 
1819b54502bSHong Zhang    Note:
1829b54502bSHong Zhang    For sparse matrix factorizations it is difficult to predict how much
1839b54502bSHong Zhang    fill to expect. By running with the option -log_info PETSc will print the
1849b54502bSHong Zhang    actual amount of fill used; allowing you to set the value accurately for
1859b54502bSHong Zhang    future runs. But default PETSc uses a value of 1.0
1869b54502bSHong Zhang 
1879b54502bSHong Zhang    Level: intermediate
1889b54502bSHong Zhang 
1899b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill
1909b54502bSHong Zhang 
1919b54502bSHong Zhang .seealso: PCLUSetFill()
1929b54502bSHong Zhang @*/
1939b54502bSHong Zhang PetscErrorCode PCICCSetFill(PC pc,PetscReal fill)
1949b54502bSHong Zhang {
1959b54502bSHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscReal);
1969b54502bSHong Zhang 
1979b54502bSHong Zhang   PetscFunctionBegin;
1989b54502bSHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
1999b54502bSHong Zhang   if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0");
2009b54502bSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr);
2019b54502bSHong Zhang   if (f) {
2029b54502bSHong Zhang     ierr = (*f)(pc,fill);CHKERRQ(ierr);
2039b54502bSHong Zhang   }
2049b54502bSHong Zhang   PetscFunctionReturn(0);
2059b54502bSHong Zhang }
2069b54502bSHong Zhang 
2079b54502bSHong Zhang #undef __FUNCT__
2089b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC"
2099b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc)
2109b54502bSHong Zhang {
2119b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2129b54502bSHong Zhang   IS             perm,cperm;
2139b54502bSHong Zhang   PetscErrorCode ierr;
2149b54502bSHong Zhang 
2159b54502bSHong Zhang   PetscFunctionBegin;
2169b54502bSHong Zhang   ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr);
2179b54502bSHong Zhang 
2189b54502bSHong Zhang   if (!pc->setupcalled) {
2199b54502bSHong Zhang     ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr);
2209b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
2219b54502bSHong Zhang     ierr = MatDestroy(icc->fact);CHKERRQ(ierr);
2229b54502bSHong Zhang     ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr);
2239b54502bSHong Zhang   }
2249b54502bSHong Zhang   ierr = ISDestroy(cperm);CHKERRQ(ierr);
2259b54502bSHong Zhang   ierr = ISDestroy(perm);CHKERRQ(ierr);
2269b54502bSHong Zhang   ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr);
2279b54502bSHong Zhang   PetscFunctionReturn(0);
2289b54502bSHong Zhang }
2299b54502bSHong Zhang 
2309b54502bSHong Zhang #undef __FUNCT__
2319b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
2329b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
2339b54502bSHong Zhang {
2349b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2359b54502bSHong Zhang   PetscErrorCode ierr;
2369b54502bSHong Zhang 
2379b54502bSHong Zhang   PetscFunctionBegin;
2389b54502bSHong Zhang   if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);}
2399b54502bSHong Zhang   ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr);
2409b54502bSHong Zhang   ierr = PetscFree(icc);CHKERRQ(ierr);
2419b54502bSHong Zhang   PetscFunctionReturn(0);
2429b54502bSHong Zhang }
2439b54502bSHong Zhang 
2449b54502bSHong Zhang #undef __FUNCT__
2459b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
2469b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
2479b54502bSHong Zhang {
2489b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2499b54502bSHong Zhang   PetscErrorCode ierr;
2509b54502bSHong Zhang 
2519b54502bSHong Zhang   PetscFunctionBegin;
2529b54502bSHong Zhang   ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr);
2539b54502bSHong Zhang   PetscFunctionReturn(0);
2549b54502bSHong Zhang }
2559b54502bSHong Zhang 
2569b54502bSHong Zhang #undef __FUNCT__
2579b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
2589b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
2599b54502bSHong Zhang {
2609b54502bSHong Zhang   PetscErrorCode ierr;
2619b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2629b54502bSHong Zhang 
2639b54502bSHong Zhang   PetscFunctionBegin;
2649b54502bSHong Zhang   ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr);
2659b54502bSHong Zhang   PetscFunctionReturn(0);
2669b54502bSHong Zhang }
2679b54502bSHong Zhang 
2689b54502bSHong Zhang #undef __FUNCT__
2699b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
2709b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
2719b54502bSHong Zhang {
2729b54502bSHong Zhang   PetscErrorCode ierr;
2739b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2749b54502bSHong Zhang 
2759b54502bSHong Zhang   PetscFunctionBegin;
2769b54502bSHong Zhang   ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr);
2779b54502bSHong Zhang   PetscFunctionReturn(0);
2789b54502bSHong Zhang }
2799b54502bSHong Zhang 
2809b54502bSHong Zhang #undef __FUNCT__
2819b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ICC"
2829b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat)
2839b54502bSHong Zhang {
2849b54502bSHong Zhang   PC_ICC *icc = (PC_ICC*)pc->data;
2859b54502bSHong Zhang 
2869b54502bSHong Zhang   PetscFunctionBegin;
2879b54502bSHong Zhang   *mat = icc->fact;
2889b54502bSHong Zhang   PetscFunctionReturn(0);
2899b54502bSHong Zhang }
2909b54502bSHong Zhang 
2919b54502bSHong Zhang #undef __FUNCT__
2929b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
2939b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc)
2949b54502bSHong Zhang {
2959b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
2969b54502bSHong Zhang   char           tname[256];
2979b54502bSHong Zhang   PetscTruth     flg;
2989b54502bSHong Zhang   PetscErrorCode ierr;
2999b54502bSHong Zhang   PetscFList     ordlist;
3009b54502bSHong Zhang 
3019b54502bSHong Zhang   PetscFunctionBegin;
3029b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
3039b54502bSHong Zhang   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
3049b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr);
3059b54502bSHong Zhang     ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr);
3069b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
3079b54502bSHong Zhang     ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr);
3089b54502bSHong Zhang     if (flg) {
3099b54502bSHong Zhang       ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr);
3109b54502bSHong Zhang     }
3110a29876aSHong Zhang     ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
3129b54502bSHong Zhang     if (flg) {
313afaefe49SHong Zhang       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
3149b54502bSHong Zhang     }
3150a29876aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",icc->info.shiftnz,&icc->info.shiftnz,0);CHKERRQ(ierr);
3160a29876aSHong Zhang     ierr = PetscOptionsName("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr);
3179b54502bSHong Zhang     if (flg) {
318afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr);
3199b54502bSHong Zhang     } else {
320afaefe49SHong Zhang       ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr);
3219b54502bSHong Zhang     }
322ee45ca4aSHong Zhang     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr);
3239b54502bSHong Zhang 
3249b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3259b54502bSHong Zhang   PetscFunctionReturn(0);
3269b54502bSHong Zhang }
3279b54502bSHong Zhang 
3289b54502bSHong Zhang #undef __FUNCT__
3299b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
3309b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
3319b54502bSHong Zhang {
3329b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
3339b54502bSHong Zhang   PetscErrorCode ierr;
3349b54502bSHong Zhang   PetscTruth     isstring,iascii;
3359b54502bSHong Zhang 
3369b54502bSHong Zhang   PetscFunctionBegin;
3379b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3389b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
3399b54502bSHong Zhang   if (iascii) {
3409b54502bSHong Zhang     if (icc->info.levels == 1) {
3419b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr);
3429b54502bSHong Zhang     } else {
3439b54502bSHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr);
3449b54502bSHong Zhang     }
3459b54502bSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr);
3460a29876aSHong Zhang     if (icc->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using Manteuffel shift\n");CHKERRQ(ierr);}
3479b54502bSHong Zhang   } else if (isstring) {
3489b54502bSHong Zhang     ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr);
3499b54502bSHong Zhang   } else {
3509b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name);
3519b54502bSHong Zhang   }
3529b54502bSHong Zhang   PetscFunctionReturn(0);
3539b54502bSHong Zhang }
3549b54502bSHong Zhang 
3559b54502bSHong Zhang /*MC
3569b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
3579b54502bSHong Zhang 
3589b54502bSHong Zhang    Options Database Keys:
3599b54502bSHong Zhang +  -pc_icc_levels <k> - number of levels of fill for ICC(k)
3609b54502bSHong Zhang .  -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
3619b54502bSHong Zhang                       its factorization (overwrites original matrix)
3629b54502bSHong Zhang .  -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
3639b54502bSHong Zhang -  -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
3649b54502bSHong Zhang 
3659b54502bSHong Zhang    Level: beginner
3669b54502bSHong Zhang 
3679b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
3689b54502bSHong Zhang 
3699b54502bSHong Zhang    Notes: Only implemented for some matrix formats. Not implemented in parallel
3709b54502bSHong Zhang 
3719b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
3729b54502bSHong Zhang 
3739b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
3749b54502bSHong Zhang 
3759b54502bSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE);
3769b54502bSHong Zhang           to turn off the shift.
3779b54502bSHong Zhang 
3789b54502bSHong Zhang 
3799b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
380ee45ca4aSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
3819b54502bSHong Zhang            PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(),
3829b54502bSHong Zhang            PCICCSetLevels()
3839b54502bSHong Zhang 
3849b54502bSHong Zhang M*/
3859b54502bSHong Zhang 
3869b54502bSHong Zhang EXTERN_C_BEGIN
3879b54502bSHong Zhang #undef __FUNCT__
3889b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
3899b54502bSHong Zhang PetscErrorCode PCCreate_ICC(PC pc)
3909b54502bSHong Zhang {
3919b54502bSHong Zhang   PetscErrorCode ierr;
3929b54502bSHong Zhang   PC_ICC         *icc;
3939b54502bSHong Zhang 
3949b54502bSHong Zhang   PetscFunctionBegin;
3959b54502bSHong Zhang   ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr);
396*52e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_ICC));CHKERRQ(ierr);
3979b54502bSHong Zhang 
3989b54502bSHong Zhang   icc->fact	          = 0;
3999b54502bSHong Zhang   ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr);
4009b54502bSHong Zhang   ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr);
4019b54502bSHong Zhang   icc->info.levels	  = 0;
4029b54502bSHong Zhang   icc->info.fill          = 1.0;
4039b54502bSHong Zhang   icc->implctx            = 0;
4049b54502bSHong Zhang 
4059b54502bSHong Zhang   icc->info.dtcol              = PETSC_DEFAULT;
4060a29876aSHong Zhang   icc->info.shiftnz            = 0.0;
4070a29876aSHong Zhang   icc->info.shiftpd            = PETSC_TRUE;
4089b54502bSHong Zhang   icc->info.shift_fraction     = 0.0;
4099b54502bSHong Zhang   icc->info.zeropivot          = 1.e-12;
4109b54502bSHong Zhang   pc->data	               = (void*)icc;
4119b54502bSHong Zhang 
4129b54502bSHong Zhang   pc->ops->apply	       = PCApply_ICC;
4139b54502bSHong Zhang   pc->ops->setup               = PCSetup_ICC;
4149b54502bSHong Zhang   pc->ops->destroy	       = PCDestroy_ICC;
4159b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
4169b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
4179b54502bSHong Zhang   pc->ops->getfactoredmatrix   = PCGetFactoredMatrix_ICC;
4189b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
4199b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
4209b54502bSHong Zhang 
421afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ICC",
422afaefe49SHong Zhang                     PCFactorSetZeroPivot_ICC);CHKERRQ(ierr);
423afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ICC",
424afaefe49SHong Zhang                     PCFactorSetShiftNonzero_ICC);CHKERRQ(ierr);
425afaefe49SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ICC",
426afaefe49SHong Zhang                     PCFactorSetShiftPd_ICC);CHKERRQ(ierr);
427afaefe49SHong Zhang 
4289b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC",
4299b54502bSHong Zhang                     PCICCSetLevels_ICC);CHKERRQ(ierr);
4309b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC",
4319b54502bSHong Zhang                     PCICCSetFill_ICC);CHKERRQ(ierr);
4329b54502bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC",
4339b54502bSHong Zhang                     PCICCSetMatOrdering_ICC);CHKERRQ(ierr);
4349b54502bSHong Zhang   PetscFunctionReturn(0);
4359b54502bSHong Zhang }
4369b54502bSHong Zhang EXTERN_C_END
4379b54502bSHong Zhang 
4389b54502bSHong Zhang 
439