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__ 109b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering_ICC" 119b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering_ICC(PC pc,MatOrderingType ordering) 129b54502bSHong Zhang { 139b54502bSHong Zhang PC_ICC *dir = (PC_ICC*)pc->data; 149b54502bSHong Zhang PetscErrorCode ierr; 159b54502bSHong Zhang 169b54502bSHong Zhang PetscFunctionBegin; 179b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 189b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 199b54502bSHong Zhang PetscFunctionReturn(0); 209b54502bSHong Zhang } 219b54502bSHong Zhang EXTERN_C_END 229b54502bSHong Zhang 239b54502bSHong Zhang EXTERN_C_BEGIN 249b54502bSHong Zhang #undef __FUNCT__ 259b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill_ICC" 269b54502bSHong Zhang PetscErrorCode PCICCSetFill_ICC(PC pc,PetscReal fill) 279b54502bSHong Zhang { 289b54502bSHong Zhang PC_ICC *dir; 299b54502bSHong Zhang 309b54502bSHong Zhang PetscFunctionBegin; 319b54502bSHong Zhang dir = (PC_ICC*)pc->data; 329b54502bSHong Zhang dir->info.fill = fill; 339b54502bSHong Zhang PetscFunctionReturn(0); 349b54502bSHong Zhang } 359b54502bSHong Zhang EXTERN_C_END 369b54502bSHong Zhang 379b54502bSHong Zhang EXTERN_C_BEGIN 389b54502bSHong Zhang #undef __FUNCT__ 399b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels_ICC" 409b54502bSHong Zhang PetscErrorCode PCICCSetLevels_ICC(PC pc,PetscInt levels) 419b54502bSHong Zhang { 429b54502bSHong Zhang PC_ICC *icc; 439b54502bSHong Zhang 449b54502bSHong Zhang PetscFunctionBegin; 459b54502bSHong Zhang icc = (PC_ICC*)pc->data; 469b54502bSHong Zhang icc->info.levels = levels; 479b54502bSHong Zhang PetscFunctionReturn(0); 489b54502bSHong Zhang } 499b54502bSHong Zhang EXTERN_C_END 509b54502bSHong Zhang 519b54502bSHong Zhang #undef __FUNCT__ 529b54502bSHong Zhang #define __FUNCT__ "PCICCSetMatOrdering" 539b54502bSHong Zhang /*@ 549b54502bSHong Zhang PCICCSetMatOrdering - Sets the ordering routine (to reduce fill) to 559b54502bSHong Zhang be used it the ICC factorization. 569b54502bSHong Zhang 579b54502bSHong Zhang Collective on PC 589b54502bSHong Zhang 599b54502bSHong Zhang Input Parameters: 609b54502bSHong Zhang + pc - the preconditioner context 619b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 629b54502bSHong Zhang 639b54502bSHong Zhang Options Database Key: 649b54502bSHong Zhang . -pc_icc_mat_ordering_type <nd,rcm,...> - Sets ordering routine 659b54502bSHong Zhang 669b54502bSHong Zhang Level: intermediate 679b54502bSHong Zhang 689b54502bSHong Zhang .seealso: PCLUSetMatOrdering() 699b54502bSHong Zhang 709b54502bSHong Zhang .keywords: PC, ICC, set, matrix, reordering 719b54502bSHong Zhang 729b54502bSHong Zhang @*/ 739b54502bSHong Zhang PetscErrorCode PCICCSetMatOrdering(PC pc,MatOrderingType ordering) 749b54502bSHong Zhang { 759b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 769b54502bSHong Zhang 779b54502bSHong Zhang PetscFunctionBegin; 789b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 799b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 809b54502bSHong Zhang if (f) { 819b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 829b54502bSHong Zhang } 839b54502bSHong Zhang PetscFunctionReturn(0); 849b54502bSHong Zhang } 859b54502bSHong Zhang 869b54502bSHong Zhang #undef __FUNCT__ 879b54502bSHong Zhang #define __FUNCT__ "PCICCSetLevels" 889b54502bSHong Zhang /*@ 899b54502bSHong Zhang PCICCSetLevels - Sets the number of levels of fill to use. 909b54502bSHong Zhang 919b54502bSHong Zhang Collective on PC 929b54502bSHong Zhang 939b54502bSHong Zhang Input Parameters: 949b54502bSHong Zhang + pc - the preconditioner context 959b54502bSHong Zhang - levels - number of levels of fill 969b54502bSHong Zhang 979b54502bSHong Zhang Options Database Key: 989b54502bSHong Zhang . -pc_icc_levels <levels> - Sets fill level 999b54502bSHong Zhang 1009b54502bSHong Zhang Level: intermediate 1019b54502bSHong Zhang 1029b54502bSHong Zhang Concepts: ICC^setting levels of fill 1039b54502bSHong Zhang 1049b54502bSHong Zhang @*/ 1059b54502bSHong Zhang PetscErrorCode PCICCSetLevels(PC pc,PetscInt levels) 1069b54502bSHong Zhang { 1079b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 1089b54502bSHong Zhang 1099b54502bSHong Zhang PetscFunctionBegin; 1109b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1119b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 1129b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 1139b54502bSHong Zhang if (f) { 1149b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 1159b54502bSHong Zhang } 1169b54502bSHong Zhang PetscFunctionReturn(0); 1179b54502bSHong Zhang } 1189b54502bSHong Zhang 1199b54502bSHong Zhang #undef __FUNCT__ 1209b54502bSHong Zhang #define __FUNCT__ "PCICCSetFill" 1219b54502bSHong Zhang /*@ 1229b54502bSHong Zhang PCICCSetFill - Indicate the amount of fill you expect in the factored matrix, 1239b54502bSHong Zhang where fill = number nonzeros in factor/number nonzeros in original matrix. 1249b54502bSHong Zhang 1259b54502bSHong Zhang Collective on PC 1269b54502bSHong Zhang 1279b54502bSHong Zhang Input Parameters: 1289b54502bSHong Zhang + pc - the preconditioner context 1299b54502bSHong Zhang - fill - amount of expected fill 1309b54502bSHong Zhang 1319b54502bSHong Zhang Options Database Key: 1329b54502bSHong Zhang $ -pc_icc_fill <fill> 1339b54502bSHong Zhang 1349b54502bSHong Zhang Note: 1359b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 1369b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 1379b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 1389b54502bSHong Zhang future runs. But default PETSc uses a value of 1.0 1399b54502bSHong Zhang 1409b54502bSHong Zhang Level: intermediate 1419b54502bSHong Zhang 1429b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 1439b54502bSHong Zhang 1449b54502bSHong Zhang .seealso: PCLUSetFill() 1459b54502bSHong Zhang @*/ 1469b54502bSHong Zhang PetscErrorCode PCICCSetFill(PC pc,PetscReal fill) 1479b54502bSHong Zhang { 1489b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 1499b54502bSHong Zhang 1509b54502bSHong Zhang PetscFunctionBegin; 1519b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1529b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 1539b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCICCSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 1549b54502bSHong Zhang if (f) { 1559b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 1569b54502bSHong Zhang } 1579b54502bSHong Zhang PetscFunctionReturn(0); 1589b54502bSHong Zhang } 1599b54502bSHong Zhang 1609b54502bSHong Zhang #undef __FUNCT__ 1619b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 1629b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 1639b54502bSHong Zhang { 1649b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1659b54502bSHong Zhang IS perm,cperm; 1669b54502bSHong Zhang PetscErrorCode ierr; 1679b54502bSHong Zhang 1689b54502bSHong Zhang PetscFunctionBegin; 1699b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,icc->ordering,&perm,&cperm);CHKERRQ(ierr); 1709b54502bSHong Zhang 1719b54502bSHong Zhang if (!pc->setupcalled) { 1729b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1739b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1749b54502bSHong Zhang ierr = MatDestroy(icc->fact);CHKERRQ(ierr); 1759b54502bSHong Zhang ierr = MatICCFactorSymbolic(pc->pmat,perm,&icc->info,&icc->fact);CHKERRQ(ierr); 1769b54502bSHong Zhang } 1779b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 1789b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 1799b54502bSHong Zhang ierr = MatCholeskyFactorNumeric(pc->pmat,&icc->info,&icc->fact);CHKERRQ(ierr); 1809b54502bSHong Zhang PetscFunctionReturn(0); 1819b54502bSHong Zhang } 1829b54502bSHong Zhang 1839b54502bSHong Zhang #undef __FUNCT__ 1849b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 1859b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 1869b54502bSHong Zhang { 1879b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1889b54502bSHong Zhang PetscErrorCode ierr; 1899b54502bSHong Zhang 1909b54502bSHong Zhang PetscFunctionBegin; 1919b54502bSHong Zhang if (icc->fact) {ierr = MatDestroy(icc->fact);CHKERRQ(ierr);} 1929b54502bSHong Zhang ierr = PetscStrfree(icc->ordering);CHKERRQ(ierr); 1939b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 1949b54502bSHong Zhang PetscFunctionReturn(0); 1959b54502bSHong Zhang } 1969b54502bSHong Zhang 1979b54502bSHong Zhang #undef __FUNCT__ 1989b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 1999b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 2009b54502bSHong Zhang { 2019b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2029b54502bSHong Zhang PetscErrorCode ierr; 2039b54502bSHong Zhang 2049b54502bSHong Zhang PetscFunctionBegin; 2059b54502bSHong Zhang ierr = MatSolve(icc->fact,x,y);CHKERRQ(ierr); 2069b54502bSHong Zhang PetscFunctionReturn(0); 2079b54502bSHong Zhang } 2089b54502bSHong Zhang 2099b54502bSHong Zhang #undef __FUNCT__ 2109b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 2119b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 2129b54502bSHong Zhang { 2139b54502bSHong Zhang PetscErrorCode ierr; 2149b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2159b54502bSHong Zhang 2169b54502bSHong Zhang PetscFunctionBegin; 2179b54502bSHong Zhang ierr = MatForwardSolve(icc->fact,x,y);CHKERRQ(ierr); 2189b54502bSHong Zhang PetscFunctionReturn(0); 2199b54502bSHong Zhang } 2209b54502bSHong Zhang 2219b54502bSHong Zhang #undef __FUNCT__ 2229b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 2239b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 2249b54502bSHong Zhang { 2259b54502bSHong Zhang PetscErrorCode ierr; 2269b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2279b54502bSHong Zhang 2289b54502bSHong Zhang PetscFunctionBegin; 2299b54502bSHong Zhang ierr = MatBackwardSolve(icc->fact,x,y);CHKERRQ(ierr); 2309b54502bSHong Zhang PetscFunctionReturn(0); 2319b54502bSHong Zhang } 2329b54502bSHong Zhang 2339b54502bSHong Zhang #undef __FUNCT__ 2349b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ICC" 2359b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ICC(PC pc,Mat *mat) 2369b54502bSHong Zhang { 2379b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2389b54502bSHong Zhang 2399b54502bSHong Zhang PetscFunctionBegin; 2409b54502bSHong Zhang *mat = icc->fact; 2419b54502bSHong Zhang PetscFunctionReturn(0); 2429b54502bSHong Zhang } 2439b54502bSHong Zhang 2449b54502bSHong Zhang #undef __FUNCT__ 2459b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 2469b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 2479b54502bSHong Zhang { 2489b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2499b54502bSHong Zhang char tname[256]; 2509b54502bSHong Zhang PetscTruth flg; 2519b54502bSHong Zhang PetscErrorCode ierr; 2529b54502bSHong Zhang PetscFList ordlist; 2539b54502bSHong Zhang 2549b54502bSHong Zhang PetscFunctionBegin; 2559b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 2569b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 2579b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_levels","levels of fill","PCICCSetLevels",icc->info.levels,&icc->info.levels,&flg);CHKERRQ(ierr); 2589b54502bSHong Zhang ierr = PetscOptionsReal("-pc_icc_fill","Expected fill in factorization","PCICCSetFill",icc->info.fill,&icc->info.fill,&flg);CHKERRQ(ierr); 2599b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 2609b54502bSHong Zhang ierr = PetscOptionsList("-pc_icc_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCICCSetMatOrdering",ordlist,icc->ordering,tname,256,&flg);CHKERRQ(ierr); 2619b54502bSHong Zhang if (flg) { 2629b54502bSHong Zhang ierr = PCICCSetMatOrdering(pc,tname);CHKERRQ(ierr); 2639b54502bSHong Zhang } 2640a29876aSHong Zhang ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 2659b54502bSHong Zhang if (flg) { 2660a29876aSHong Zhang ierr = PCFactorSetShiftNonzero((PetscReal) PETSC_DECIDE,&icc->info);CHKERRQ(ierr); 2679b54502bSHong Zhang } 2680a29876aSHong Zhang ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",icc->info.shiftnz,&icc->info.shiftnz,0);CHKERRQ(ierr); 2690a29876aSHong Zhang ierr = PetscOptionsName("-pc_factor_shiftpd","Manteuffel shift applied to diagonal","PCICCSetShift",&flg);CHKERRQ(ierr); 2709b54502bSHong Zhang if (flg) { 271*ee45ca4aSHong Zhang ierr = PCFactorSetShiftPd(PETSC_TRUE,&icc->info);CHKERRQ(ierr); 2729b54502bSHong Zhang } else { 273*ee45ca4aSHong Zhang ierr = PCFactorSetShiftPd(PETSC_FALSE,&icc->info);CHKERRQ(ierr); 2749b54502bSHong Zhang } 275*ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",icc->info.zeropivot,&icc->info.zeropivot,0);CHKERRQ(ierr); 2769b54502bSHong Zhang 2779b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2789b54502bSHong Zhang PetscFunctionReturn(0); 2799b54502bSHong Zhang } 2809b54502bSHong Zhang 2819b54502bSHong Zhang #undef __FUNCT__ 2829b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 2839b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 2849b54502bSHong Zhang { 2859b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 2869b54502bSHong Zhang PetscErrorCode ierr; 2879b54502bSHong Zhang PetscTruth isstring,iascii; 2889b54502bSHong Zhang 2899b54502bSHong Zhang PetscFunctionBegin; 2909b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 2919b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 2929b54502bSHong Zhang if (iascii) { 2939b54502bSHong Zhang if (icc->info.levels == 1) { 2949b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2959b54502bSHong Zhang } else { 2969b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)icc->info.levels);CHKERRQ(ierr); 2979b54502bSHong Zhang } 2989b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICC: max fill ratio allocated %g\n",icc->info.fill);CHKERRQ(ierr); 2990a29876aSHong Zhang if (icc->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 3009b54502bSHong Zhang } else if (isstring) { 3019b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)icc->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 3029b54502bSHong Zhang } else { 3039b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 3049b54502bSHong Zhang } 3059b54502bSHong Zhang PetscFunctionReturn(0); 3069b54502bSHong Zhang } 3079b54502bSHong Zhang 3089b54502bSHong Zhang /*MC 3099b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 3109b54502bSHong Zhang 3119b54502bSHong Zhang Options Database Keys: 3129b54502bSHong Zhang + -pc_icc_levels <k> - number of levels of fill for ICC(k) 3139b54502bSHong Zhang . -pc_icc_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 3149b54502bSHong Zhang its factorization (overwrites original matrix) 3159b54502bSHong Zhang . -pc_icc_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 3169b54502bSHong Zhang - -pc_icc_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 3179b54502bSHong Zhang 3189b54502bSHong Zhang Level: beginner 3199b54502bSHong Zhang 3209b54502bSHong Zhang Concepts: incomplete Cholesky factorization 3219b54502bSHong Zhang 3229b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 3239b54502bSHong Zhang 3249b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 3259b54502bSHong Zhang 3269b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 3279b54502bSHong Zhang 3289b54502bSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCICCSetShift(pc,PETSC_FALSE); 3299b54502bSHong Zhang to turn off the shift. 3309b54502bSHong Zhang 3319b54502bSHong Zhang 3329b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 333*ee45ca4aSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), 3349b54502bSHong Zhang PCICCSetFill(), PCICCSetMatOrdering(), PCICCSetReuseOrdering(), 3359b54502bSHong Zhang PCICCSetLevels() 3369b54502bSHong Zhang 3379b54502bSHong Zhang M*/ 3389b54502bSHong Zhang 3399b54502bSHong Zhang EXTERN_C_BEGIN 3409b54502bSHong Zhang #undef __FUNCT__ 3419b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 3429b54502bSHong Zhang PetscErrorCode PCCreate_ICC(PC pc) 3439b54502bSHong Zhang { 3449b54502bSHong Zhang PetscErrorCode ierr; 3459b54502bSHong Zhang PC_ICC *icc; 3469b54502bSHong Zhang 3479b54502bSHong Zhang PetscFunctionBegin; 3489b54502bSHong Zhang ierr = PetscNew(PC_ICC,&icc);CHKERRQ(ierr); 3499b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_ICC)); 3509b54502bSHong Zhang 3519b54502bSHong Zhang icc->fact = 0; 3529b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&icc->ordering);CHKERRQ(ierr); 3539b54502bSHong Zhang ierr = MatFactorInfoInitialize(&icc->info);CHKERRQ(ierr); 3549b54502bSHong Zhang icc->info.levels = 0; 3559b54502bSHong Zhang icc->info.fill = 1.0; 3569b54502bSHong Zhang icc->implctx = 0; 3579b54502bSHong Zhang 3589b54502bSHong Zhang icc->info.dtcol = PETSC_DEFAULT; 3590a29876aSHong Zhang icc->info.shiftnz = 0.0; 3600a29876aSHong Zhang icc->info.shiftpd = PETSC_TRUE; 3619b54502bSHong Zhang icc->info.shift_fraction = 0.0; 3629b54502bSHong Zhang icc->info.zeropivot = 1.e-12; 3639b54502bSHong Zhang pc->data = (void*)icc; 3649b54502bSHong Zhang 3659b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 3669b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 3679b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 3689b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 3699b54502bSHong Zhang pc->ops->view = PCView_ICC; 3709b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ICC; 3719b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 3729b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 3739b54502bSHong Zhang 3749b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetLevels_C","PCICCSetLevels_ICC", 3759b54502bSHong Zhang PCICCSetLevels_ICC);CHKERRQ(ierr); 3769b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetFill_C","PCICCSetFill_ICC", 3779b54502bSHong Zhang PCICCSetFill_ICC);CHKERRQ(ierr); 3789b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCICCSetMatOrdering_C","PCICCSetMatOrdering_ICC", 3799b54502bSHong Zhang PCICCSetMatOrdering_ICC);CHKERRQ(ierr); 3809b54502bSHong Zhang PetscFunctionReturn(0); 3819b54502bSHong Zhang } 3829b54502bSHong Zhang EXTERN_C_END 3839b54502bSHong Zhang 3849b54502bSHong Zhang 385