1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 39b54502bSHong Zhang /* 49b54502bSHong Zhang Defines a Cholesky factorization preconditioner for any Mat implementation. 59b54502bSHong Zhang Presently only provided for MPIRowbs format (i.e. BlockSolve). 69b54502bSHong Zhang */ 79b54502bSHong Zhang 885b398ccSKris Buschelman #include "src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 99b54502bSHong Zhang 109b54502bSHong Zhang EXTERN_C_BEGIN 119b54502bSHong Zhang #undef __FUNCT__ 12afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ICC" 13dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ICC(PC pc,PetscReal z) 14afaefe49SHong Zhang { 15afaefe49SHong Zhang PC_ICC *icc; 16afaefe49SHong Zhang 17afaefe49SHong Zhang PetscFunctionBegin; 18afaefe49SHong Zhang icc = (PC_ICC*)pc->data; 19afaefe49SHong Zhang icc->info.zeropivot = z; 20afaefe49SHong Zhang PetscFunctionReturn(0); 21afaefe49SHong Zhang } 22afaefe49SHong Zhang EXTERN_C_END 23afaefe49SHong Zhang 24afaefe49SHong Zhang EXTERN_C_BEGIN 25afaefe49SHong Zhang #undef __FUNCT__ 26afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ICC" 27dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ICC(PC pc,PetscReal shift) 28afaefe49SHong Zhang { 29afaefe49SHong Zhang PC_ICC *dir; 30afaefe49SHong Zhang 31afaefe49SHong Zhang PetscFunctionBegin; 32afaefe49SHong Zhang dir = (PC_ICC*)pc->data; 33afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 34afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 35afaefe49SHong Zhang } else { 36afaefe49SHong Zhang dir->info.shiftnz = shift; 37afaefe49SHong Zhang } 38afaefe49SHong Zhang PetscFunctionReturn(0); 39afaefe49SHong Zhang } 40afaefe49SHong Zhang EXTERN_C_END 41afaefe49SHong Zhang 42afaefe49SHong Zhang EXTERN_C_BEGIN 43afaefe49SHong Zhang #undef __FUNCT__ 44afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ICC" 45dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ICC(PC pc,PetscTruth shift) 46afaefe49SHong Zhang { 47afaefe49SHong Zhang PC_ICC *dir; 48afaefe49SHong Zhang 49afaefe49SHong Zhang PetscFunctionBegin; 50afaefe49SHong Zhang dir = (PC_ICC*)pc->data; 51fbf22428SSatish Balay if (shift) { 52fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 53fbf22428SSatish Balay dir->info.shiftpd = 1.0; 54fbf22428SSatish Balay } else { 55fbf22428SSatish Balay dir->info.shiftpd = 0.0; 56fbf22428SSatish Balay } 57afaefe49SHong Zhang PetscFunctionReturn(0); 58afaefe49SHong Zhang } 59afaefe49SHong Zhang EXTERN_C_END 60afaefe49SHong Zhang 61afaefe49SHong Zhang EXTERN_C_BEGIN 62afaefe49SHong Zhang #undef __FUNCT__ 63e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ICC" 64e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ICC(PC pc,MatOrderingType ordering) 659b54502bSHong Zhang { 669b54502bSHong Zhang PC_ICC *dir = (PC_ICC*)pc->data; 679b54502bSHong Zhang PetscErrorCode ierr; 689b54502bSHong Zhang 699b54502bSHong Zhang PetscFunctionBegin; 709b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 719b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 729b54502bSHong Zhang PetscFunctionReturn(0); 739b54502bSHong Zhang } 749b54502bSHong Zhang EXTERN_C_END 759b54502bSHong Zhang 769b54502bSHong Zhang EXTERN_C_BEGIN 779b54502bSHong Zhang #undef __FUNCT__ 7855ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ICC" 7955ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ICC(PC pc,PetscReal fill) 809b54502bSHong Zhang { 819b54502bSHong Zhang PC_ICC *dir; 829b54502bSHong Zhang 839b54502bSHong Zhang PetscFunctionBegin; 849b54502bSHong Zhang dir = (PC_ICC*)pc->data; 859b54502bSHong Zhang dir->info.fill = fill; 869b54502bSHong Zhang PetscFunctionReturn(0); 879b54502bSHong Zhang } 889b54502bSHong Zhang EXTERN_C_END 899b54502bSHong Zhang 909b54502bSHong Zhang EXTERN_C_BEGIN 919b54502bSHong Zhang #undef __FUNCT__ 922401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels_ICC" 932401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ICC(PC pc,PetscInt levels) 949b54502bSHong Zhang { 959b54502bSHong Zhang PC_ICC *icc; 969b54502bSHong Zhang 979b54502bSHong Zhang PetscFunctionBegin; 989b54502bSHong Zhang icc = (PC_ICC*)pc->data; 999b54502bSHong Zhang icc->info.levels = levels; 1009b54502bSHong Zhang PetscFunctionReturn(0); 1019b54502bSHong Zhang } 1029b54502bSHong Zhang EXTERN_C_END 1039b54502bSHong Zhang 1049b54502bSHong Zhang #undef __FUNCT__ 1059b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 1069b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 1079b54502bSHong Zhang { 1089b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1099b54502bSHong Zhang IS perm,cperm; 1109b54502bSHong Zhang PetscErrorCode ierr; 111f3a39becSBarry Smith MatInfo info; 1129b54502bSHong Zhang 1139b54502bSHong Zhang PetscFunctionBegin; 1149b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 1159b54502bSHong Zhang 1169b54502bSHong Zhang if (!pc->setupcalled) { 1179b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1189b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1199b54502bSHong Zhang ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 1209b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1219b54502bSHong Zhang } 122f3a39becSBarry Smith ierr = MatGetInfo(icc->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 123f3a39becSBarry Smith icc->actualfill = info.fill_ratio_needed; 124f3a39becSBarry Smith 1259b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 1269b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 1279b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 1289b54502bSHong Zhang PetscFunctionReturn(0); 1299b54502bSHong Zhang } 1309b54502bSHong Zhang 1319b54502bSHong Zhang #undef __FUNCT__ 1329b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 1339b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 1349b54502bSHong Zhang { 1359b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1369b54502bSHong Zhang PetscErrorCode ierr; 1379b54502bSHong Zhang 1389b54502bSHong Zhang PetscFunctionBegin; 1399b54502bSHong Zhang if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 1409b54502bSHong Zhang ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 1419b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 1429b54502bSHong Zhang PetscFunctionReturn(0); 1439b54502bSHong Zhang } 1449b54502bSHong Zhang 1459b54502bSHong Zhang #undef __FUNCT__ 1469b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 1479b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 1489b54502bSHong Zhang { 1499b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1509b54502bSHong Zhang PetscErrorCode ierr; 1519b54502bSHong Zhang 1529b54502bSHong Zhang PetscFunctionBegin; 1539b54502bSHong Zhang ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 1549b54502bSHong Zhang PetscFunctionReturn(0); 1559b54502bSHong Zhang } 1569b54502bSHong Zhang 1579b54502bSHong Zhang #undef __FUNCT__ 1589b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 1599b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 1609b54502bSHong Zhang { 1619b54502bSHong Zhang PetscErrorCode ierr; 1629b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1639b54502bSHong Zhang 1649b54502bSHong Zhang PetscFunctionBegin; 1659b54502bSHong Zhang ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 1669b54502bSHong Zhang PetscFunctionReturn(0); 1679b54502bSHong Zhang } 1689b54502bSHong Zhang 1699b54502bSHong Zhang #undef __FUNCT__ 1709b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 1719b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 1729b54502bSHong Zhang { 1739b54502bSHong Zhang PetscErrorCode ierr; 1749b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1759b54502bSHong Zhang 1769b54502bSHong Zhang PetscFunctionBegin; 1779b54502bSHong Zhang ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 1789b54502bSHong Zhang PetscFunctionReturn(0); 1799b54502bSHong Zhang } 1809b54502bSHong Zhang 1819b54502bSHong Zhang #undef __FUNCT__ 1829b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ICC" 1839b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 1849b54502bSHong Zhang { 1859b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1869b54502bSHong Zhang 1879b54502bSHong Zhang PetscFunctionBegin; 1889b54502bSHong Zhang *mat = icc->fact; 1899b54502bSHong Zhang PetscFunctionReturn(0); 1909b54502bSHong Zhang } 1919b54502bSHong Zhang 1929b54502bSHong Zhang #undef __FUNCT__ 1939b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 1949b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 1959b54502bSHong Zhang { 1969b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1979b54502bSHong Zhang char tname[256]; 1989b54502bSHong Zhang PetscTruth flg; 1999b54502bSHong Zhang PetscErrorCode ierr; 2009b54502bSHong Zhang PetscFList ordlist; 2019b54502bSHong Zhang 2029b54502bSHong Zhang PetscFunctionBegin; 2039b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 2049b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 2052401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 20655ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 2079b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 208e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCFactorSetMatOrderingType",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 2099b54502bSHong Zhang if (flg) { 210e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 2119b54502bSHong Zhang } 2129f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 2139b54502bSHong Zhang if (flg) { 214afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 2159b54502bSHong Zhang } 2169f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",icc->info.shiftnz,&icc->info.shiftnz,0);CHKERRQ(ierr); 2172401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShift",&flg);CHKERRQ(ierr); 2189b54502bSHong Zhang if (flg) { 219afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 2209b54502bSHong Zhang } else { 221afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 2229b54502bSHong Zhang } 223ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 2249b54502bSHong Zhang 2259b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2269b54502bSHong Zhang PetscFunctionReturn(0); 2279b54502bSHong Zhang } 2289b54502bSHong Zhang 2299b54502bSHong Zhang #undef __FUNCT__ 2309b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 2319b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 2329b54502bSHong Zhang { 2339b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2349b54502bSHong Zhang PetscErrorCode ierr; 2359b54502bSHong Zhang PetscTruth isstring,iascii; 2369b54502bSHong Zhang 2379b54502bSHong Zhang PetscFunctionBegin; 2389b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 2399b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 2409b54502bSHong Zhang if (iascii) { 2419b54502bSHong Zhang if (icc->info.levels == 1) { 2429b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2439b54502bSHong Zhang } else { 2449b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2459b54502bSHong Zhang } 246f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: factor fill ratio allocated %G\n",icc->info.fill);CHKERRQ(ierr); 2470a29876aSHong Zhang if (icc->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 248f3a39becSBarry Smith if (icc->fact) { 249f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr); 250f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 251f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 252f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 253f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 254f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 255f3a39becSBarry Smith ierr = MatView(icc->fact,viewer);CHKERRQ(ierr); 256f3a39becSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 257f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 258f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 259f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 260f3a39becSBarry Smith } 2619b54502bSHong Zhang } else if (isstring) { 2629b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 2639b54502bSHong Zhang } else { 2649b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 2659b54502bSHong Zhang } 2669b54502bSHong Zhang PetscFunctionReturn(0); 2679b54502bSHong Zhang } 2689b54502bSHong Zhang 2699b54502bSHong Zhang /*MC 2709b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 2719b54502bSHong Zhang 2729b54502bSHong Zhang Options Database Keys: 2732401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 2742401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 2759b54502bSHong Zhang its factorization (overwrites original matrix) 27655ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2772401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 278f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 279f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 280f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 2819b54502bSHong Zhang 2829b54502bSHong Zhang Level: beginner 2839b54502bSHong Zhang 2849b54502bSHong Zhang Concepts: incomplete Cholesky factorization 2859b54502bSHong Zhang 2869b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 2879b54502bSHong Zhang 2889b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 2899b54502bSHong Zhang 2909b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 2919b54502bSHong Zhang 2922401956bSBarry Smith By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE); 2939b54502bSHong Zhang to turn off the shift. 2949b54502bSHong Zhang 2959b54502bSHong Zhang 2969b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 297ee45ca4aSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 298e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 2992401956bSBarry Smith PCFactorSetLevels(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), 3009b54502bSHong Zhang 3019b54502bSHong Zhang M*/ 3029b54502bSHong Zhang 3039b54502bSHong Zhang EXTERN_C_BEGIN 3049b54502bSHong Zhang #undef __FUNCT__ 3059b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 306dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 3079b54502bSHong Zhang { 3089b54502bSHong Zhang PetscErrorCode ierr; 3099b54502bSHong Zhang PC_ICC *icc; 3109b54502bSHong Zhang 3119b54502bSHong Zhang PetscFunctionBegin; 312*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 3139b54502bSHong Zhang 3149b54502bSHong Zhang icc->fact = 0; 3159b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 3169b54502bSHong Zhang ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 3179b54502bSHong Zhang icc->info.levels = 0; 3189b54502bSHong Zhang icc->info.fill = 1.0; 3199b54502bSHong Zhang icc->implctx = 0; 3209b54502bSHong Zhang 3219b54502bSHong Zhang icc->info.dtcol = PETSC_DEFAULT; 3220a29876aSHong Zhang icc->info.shiftnz = 0.0; 323fbf22428SSatish Balay icc->info.shiftpd = 1.0; /* true */ 3249b54502bSHong Zhang icc->info.shift_fraction = 0.0; 3259b54502bSHong Zhang icc->info.zeropivot = 1.e-12; 3269b54502bSHong Zhang pc->data = (void*)icc; 3279b54502bSHong Zhang 3289b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 3299b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 3309b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 3319b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 3329b54502bSHong Zhang pc->ops->view = PCView_ICC; 3339b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 3349b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 3359b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 3369b54502bSHong Zhang 337afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ICC", 338afaefe49SHong Zhang PCFactorSetZeroPivot_ICC);CHKERRQ(ierr); 339afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ICC", 340afaefe49SHong Zhang PCFactorSetShiftNonzero_ICC);CHKERRQ(ierr); 341afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ICC", 342afaefe49SHong Zhang PCFactorSetShiftPd_ICC);CHKERRQ(ierr); 343afaefe49SHong Zhang 3442401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ICC", 3452401956bSBarry Smith PCFactorSetLevels_ICC);CHKERRQ(ierr); 34655ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ICC", 34755ba2a51SBarry Smith PCFactorSetFill_ICC);CHKERRQ(ierr); 348e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ICC", 349e5a9bf91SBarry Smith PCFactorSetMatOrderingType_ICC);CHKERRQ(ierr); 3509b54502bSHong Zhang PetscFunctionReturn(0); 3519b54502bSHong Zhang } 3529b54502bSHong Zhang EXTERN_C_END 3539b54502bSHong Zhang 3549b54502bSHong Zhang 355