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 { 10afaefe49SHong Zhang PC_ICC *icc; 11afaefe49SHong Zhang 12afaefe49SHong Zhang PetscFunctionBegin; 13afaefe49SHong Zhang icc = (PC_ICC*)pc->data; 14afaefe49SHong Zhang icc->info.zeropivot = z; 15afaefe49SHong Zhang PetscFunctionReturn(0); 16afaefe49SHong Zhang } 17afaefe49SHong Zhang EXTERN_C_END 18afaefe49SHong Zhang 19afaefe49SHong Zhang EXTERN_C_BEGIN 20afaefe49SHong Zhang #undef __FUNCT__ 21afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ICC" 22dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ICC(PC pc,PetscReal shift) 23afaefe49SHong Zhang { 24afaefe49SHong Zhang PC_ICC *dir; 25afaefe49SHong Zhang 26afaefe49SHong Zhang PetscFunctionBegin; 27afaefe49SHong Zhang dir = (PC_ICC*)pc->data; 28afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 29afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 30afaefe49SHong Zhang } else { 31afaefe49SHong Zhang dir->info.shiftnz = shift; 32afaefe49SHong Zhang } 33afaefe49SHong Zhang PetscFunctionReturn(0); 34afaefe49SHong Zhang } 35afaefe49SHong Zhang EXTERN_C_END 36afaefe49SHong Zhang 37afaefe49SHong Zhang EXTERN_C_BEGIN 38afaefe49SHong Zhang #undef __FUNCT__ 39afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ICC" 40dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ICC(PC pc,PetscTruth shift) 41afaefe49SHong Zhang { 42afaefe49SHong Zhang PC_ICC *dir; 43afaefe49SHong Zhang 44afaefe49SHong Zhang PetscFunctionBegin; 45afaefe49SHong Zhang dir = (PC_ICC*)pc->data; 46fbf22428SSatish Balay if (shift) { 47fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 48fbf22428SSatish Balay dir->info.shiftpd = 1.0; 49fbf22428SSatish Balay } else { 50fbf22428SSatish Balay dir->info.shiftpd = 0.0; 51fbf22428SSatish Balay } 52afaefe49SHong Zhang PetscFunctionReturn(0); 53afaefe49SHong Zhang } 54afaefe49SHong Zhang EXTERN_C_END 55afaefe49SHong Zhang 56afaefe49SHong Zhang EXTERN_C_BEGIN 57afaefe49SHong Zhang #undef __FUNCT__ 58e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ICC" 59e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ICC(PC pc,MatOrderingType ordering) 609b54502bSHong Zhang { 619b54502bSHong Zhang PC_ICC *dir = (PC_ICC*)pc->data; 629b54502bSHong Zhang PetscErrorCode ierr; 639b54502bSHong Zhang 649b54502bSHong Zhang PetscFunctionBegin; 65*05dd20ceSSatish Balay ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 66*05dd20ceSSatish Balay ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 679b54502bSHong Zhang PetscFunctionReturn(0); 689b54502bSHong Zhang } 699b54502bSHong Zhang EXTERN_C_END 709b54502bSHong Zhang 719b54502bSHong Zhang EXTERN_C_BEGIN 729b54502bSHong Zhang #undef __FUNCT__ 7355ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ICC" 7455ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ICC(PC pc,PetscReal fill) 759b54502bSHong Zhang { 769b54502bSHong Zhang PC_ICC *dir; 779b54502bSHong Zhang 789b54502bSHong Zhang PetscFunctionBegin; 799b54502bSHong Zhang dir = (PC_ICC*)pc->data; 809b54502bSHong Zhang dir->info.fill = fill; 819b54502bSHong Zhang PetscFunctionReturn(0); 829b54502bSHong Zhang } 839b54502bSHong Zhang EXTERN_C_END 849b54502bSHong Zhang 859b54502bSHong Zhang EXTERN_C_BEGIN 869b54502bSHong Zhang #undef __FUNCT__ 872401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels_ICC" 882401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ICC(PC pc,PetscInt levels) 899b54502bSHong Zhang { 909b54502bSHong Zhang PC_ICC *icc; 919b54502bSHong Zhang 929b54502bSHong Zhang PetscFunctionBegin; 939b54502bSHong Zhang icc = (PC_ICC*)pc->data; 949b54502bSHong Zhang icc->info.levels = levels; 959b54502bSHong Zhang PetscFunctionReturn(0); 969b54502bSHong Zhang } 979b54502bSHong Zhang EXTERN_C_END 989b54502bSHong Zhang 999b54502bSHong Zhang #undef __FUNCT__ 1009b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 1019b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 1029b54502bSHong Zhang { 1039b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1049b54502bSHong Zhang IS perm,cperm; 1059b54502bSHong Zhang PetscErrorCode ierr; 106f3a39becSBarry Smith MatInfo info; 1079b54502bSHong Zhang 1089b54502bSHong Zhang PetscFunctionBegin; 1099b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 1109b54502bSHong Zhang 1119b54502bSHong Zhang if (!pc->setupcalled) { 1129b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1139b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1149b54502bSHong Zhang ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 1159b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1169b54502bSHong Zhang } 117f3a39becSBarry Smith ierr = MatGetInfo(icc->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 118f3a39becSBarry Smith icc->actualfill = info.fill_ratio_needed; 119f3a39becSBarry Smith 1209b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 1219b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 1229b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 1239b54502bSHong Zhang PetscFunctionReturn(0); 1249b54502bSHong Zhang } 1259b54502bSHong Zhang 1269b54502bSHong Zhang #undef __FUNCT__ 1279b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 1289b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 1299b54502bSHong Zhang { 1309b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1319b54502bSHong Zhang PetscErrorCode ierr; 1329b54502bSHong Zhang 1339b54502bSHong Zhang PetscFunctionBegin; 1349b54502bSHong Zhang if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 135*05dd20ceSSatish Balay ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 1369b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 1379b54502bSHong Zhang PetscFunctionReturn(0); 1389b54502bSHong Zhang } 1399b54502bSHong Zhang 1409b54502bSHong Zhang #undef __FUNCT__ 1419b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 1429b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 1439b54502bSHong Zhang { 1449b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1459b54502bSHong Zhang PetscErrorCode ierr; 1469b54502bSHong Zhang 1479b54502bSHong Zhang PetscFunctionBegin; 1489b54502bSHong Zhang ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 1499b54502bSHong Zhang PetscFunctionReturn(0); 1509b54502bSHong Zhang } 1519b54502bSHong Zhang 1529b54502bSHong Zhang #undef __FUNCT__ 1539b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 1549b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 1559b54502bSHong Zhang { 1569b54502bSHong Zhang PetscErrorCode ierr; 1579b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1589b54502bSHong Zhang 1599b54502bSHong Zhang PetscFunctionBegin; 1609b54502bSHong Zhang ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 1619b54502bSHong Zhang PetscFunctionReturn(0); 1629b54502bSHong Zhang } 1639b54502bSHong Zhang 1649b54502bSHong Zhang #undef __FUNCT__ 1659b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 1669b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 1679b54502bSHong Zhang { 1689b54502bSHong Zhang PetscErrorCode ierr; 1699b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1709b54502bSHong Zhang 1719b54502bSHong Zhang PetscFunctionBegin; 1729b54502bSHong Zhang ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 1739b54502bSHong Zhang PetscFunctionReturn(0); 1749b54502bSHong Zhang } 1759b54502bSHong Zhang 1769b54502bSHong Zhang #undef __FUNCT__ 177a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_ICC" 178a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_ICC(PC pc,Mat *mat) 1799b54502bSHong Zhang { 1809b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1819b54502bSHong Zhang 1829b54502bSHong Zhang PetscFunctionBegin; 1839b54502bSHong Zhang *mat = icc->fact; 1849b54502bSHong Zhang PetscFunctionReturn(0); 1859b54502bSHong Zhang } 1869b54502bSHong Zhang 1879b54502bSHong Zhang #undef __FUNCT__ 1889b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 1899b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 1909b54502bSHong Zhang { 1919b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1929b54502bSHong Zhang char tname[256]; 1939b54502bSHong Zhang PetscTruth flg; 1949b54502bSHong Zhang PetscErrorCode ierr; 1959b54502bSHong Zhang PetscFList ordlist; 1969b54502bSHong Zhang 1979b54502bSHong Zhang PetscFunctionBegin; 1989b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 1999b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 2002401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 20155ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 2029b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 203e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCFactorSetMatOrderingType",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 2049b54502bSHong Zhang if (flg) { 205e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 2069b54502bSHong Zhang } 2079f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 2089b54502bSHong Zhang if (flg) { 209afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 2109b54502bSHong Zhang } 2119f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",icc->info.shiftnz,&icc->info.shiftnz,0);CHKERRQ(ierr); 21262bba022SBarry Smith flg = (icc->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE; 21362bba022SBarry Smith ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 21462bba022SBarry Smith ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr); 215ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 2169b54502bSHong Zhang 2179b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2189b54502bSHong Zhang PetscFunctionReturn(0); 2199b54502bSHong Zhang } 2209b54502bSHong Zhang 2219b54502bSHong Zhang #undef __FUNCT__ 2229b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 2239b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 2249b54502bSHong Zhang { 2259b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2269b54502bSHong Zhang PetscErrorCode ierr; 2279b54502bSHong Zhang PetscTruth isstring,iascii; 2289b54502bSHong Zhang 2299b54502bSHong Zhang PetscFunctionBegin; 2309b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 2319b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 2329b54502bSHong Zhang if (iascii) { 2339b54502bSHong Zhang if (icc->info.levels == 1) { 2349b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2359b54502bSHong Zhang } else { 2369b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2379b54502bSHong Zhang } 238f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: factor fill ratio allocated %G\n",icc->info.fill);CHKERRQ(ierr); 2390a29876aSHong Zhang if (icc->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 24062bba022SBarry Smith if (icc->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);} 241f3a39becSBarry Smith if (icc->fact) { 242f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr); 243f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 244f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 245f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 246f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 247f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 248f3a39becSBarry Smith ierr = MatView(icc->fact,viewer);CHKERRQ(ierr); 249f3a39becSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 250f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 251f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 253f3a39becSBarry Smith } 2549b54502bSHong Zhang } else if (isstring) { 2559b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 2569b54502bSHong Zhang } else { 2579b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 2589b54502bSHong Zhang } 2599b54502bSHong Zhang PetscFunctionReturn(0); 2609b54502bSHong Zhang } 2619b54502bSHong Zhang 2629b54502bSHong Zhang /*MC 2639b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 2649b54502bSHong Zhang 2659b54502bSHong Zhang Options Database Keys: 2662401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 2672401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 2689b54502bSHong Zhang its factorization (overwrites original matrix) 26955ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2702401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 271f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 272f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 273f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 2749b54502bSHong Zhang 2759b54502bSHong Zhang Level: beginner 2769b54502bSHong Zhang 2779b54502bSHong Zhang Concepts: incomplete Cholesky factorization 2789b54502bSHong Zhang 279c582cd25SBarry Smith Notes: Only implemented for some matrix formats. Not implemented in parallel (for parallel use you 280c582cd25SBarry Smith must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ICC(0) and this is not recommended 281c582cd25SBarry Smith unless you really want a parallel ICC). 2829b54502bSHong Zhang 2839b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 2849b54502bSHong Zhang 2859b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 2869b54502bSHong Zhang 2872401956bSBarry Smith By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE); 2889b54502bSHong Zhang to turn off the shift. 2899b54502bSHong Zhang 290c582cd25SBarry Smith References: 291c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 292c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 293c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 294c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 295c582cd25SBarry Smith 2969b54502bSHong Zhang 2979b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 298ee45ca4aSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 299e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 3002401956bSBarry Smith PCFactorSetLevels(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(), 3019b54502bSHong Zhang 3029b54502bSHong Zhang M*/ 3039b54502bSHong Zhang 3049b54502bSHong Zhang EXTERN_C_BEGIN 3059b54502bSHong Zhang #undef __FUNCT__ 3069b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 307dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 3089b54502bSHong Zhang { 3099b54502bSHong Zhang PetscErrorCode ierr; 3109b54502bSHong Zhang PC_ICC *icc; 3119b54502bSHong Zhang 3129b54502bSHong Zhang PetscFunctionBegin; 31338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 3149b54502bSHong Zhang 3159b54502bSHong Zhang icc->fact = 0; 316*05dd20ceSSatish Balay ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 3179b54502bSHong Zhang ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 3189b54502bSHong Zhang icc->info.levels = 0; 3199b54502bSHong Zhang icc->info.fill = 1.0; 3209b54502bSHong Zhang icc->implctx = 0; 3219b54502bSHong Zhang 3229b54502bSHong Zhang icc->info.dtcol = PETSC_DEFAULT; 3230a29876aSHong Zhang icc->info.shiftnz = 0.0; 324fbf22428SSatish Balay icc->info.shiftpd = 1.0; /* true */ 3259b54502bSHong Zhang icc->info.shift_fraction = 0.0; 3269b54502bSHong Zhang icc->info.zeropivot = 1.e-12; 3279b54502bSHong Zhang pc->data = (void*)icc; 3289b54502bSHong Zhang 3299b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 3309b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 3319b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 3329b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 3339b54502bSHong Zhang pc->ops->view = PCView_ICC; 334a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_ICC; 3359b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 3369b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 3379b54502bSHong Zhang 338afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ICC", 339afaefe49SHong Zhang PCFactorSetZeroPivot_ICC);CHKERRQ(ierr); 340afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ICC", 341afaefe49SHong Zhang PCFactorSetShiftNonzero_ICC);CHKERRQ(ierr); 342afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ICC", 343afaefe49SHong Zhang PCFactorSetShiftPd_ICC);CHKERRQ(ierr); 344afaefe49SHong Zhang 3452401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ICC", 3462401956bSBarry Smith PCFactorSetLevels_ICC);CHKERRQ(ierr); 34755ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ICC", 34855ba2a51SBarry Smith PCFactorSetFill_ICC);CHKERRQ(ierr); 349e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ICC", 350e5a9bf91SBarry Smith PCFactorSetMatOrderingType_ICC);CHKERRQ(ierr); 3519b54502bSHong Zhang PetscFunctionReturn(0); 3529b54502bSHong Zhang } 3539b54502bSHong Zhang EXTERN_C_END 3549b54502bSHong Zhang 3559b54502bSHong Zhang 356