1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 37c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 49b54502bSHong Zhang 59b54502bSHong Zhang #undef __FUNCT__ 69b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 79b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 89b54502bSHong Zhang { 99b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 109b54502bSHong Zhang IS perm,cperm; 119b54502bSHong Zhang PetscErrorCode ierr; 12f3a39becSBarry Smith MatInfo info; 139b54502bSHong Zhang 149b54502bSHong Zhang PetscFunctionBegin; 15075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr); 169b54502bSHong Zhang 179b54502bSHong Zhang if (!pc->setupcalled) { 18174acd27SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr); 19075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 209b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 21075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr); 22075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 23075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 249b54502bSHong Zhang } 25075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 26f3a39becSBarry Smith icc->actualfill = info.fill_ratio_needed; 27f3a39becSBarry Smith 289b54502bSHong Zhang ierr = ISDestroy(cperm);CHKERRQ(ierr); 299b54502bSHong Zhang ierr = ISDestroy(perm);CHKERRQ(ierr); 30075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 319b54502bSHong Zhang PetscFunctionReturn(0); 329b54502bSHong Zhang } 339b54502bSHong Zhang 349b54502bSHong Zhang #undef __FUNCT__ 359b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 369b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 379b54502bSHong Zhang { 389b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 399b54502bSHong Zhang PetscErrorCode ierr; 409b54502bSHong Zhang 419b54502bSHong Zhang PetscFunctionBegin; 42075768bcSBarry Smith if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);} 43075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 44174acd27SBarry Smith ierr = PetscStrfree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 459b54502bSHong Zhang ierr = PetscFree(icc);CHKERRQ(ierr); 469b54502bSHong Zhang PetscFunctionReturn(0); 479b54502bSHong Zhang } 489b54502bSHong Zhang 499b54502bSHong Zhang #undef __FUNCT__ 509b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 519b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 529b54502bSHong Zhang { 539b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 549b54502bSHong Zhang PetscErrorCode ierr; 559b54502bSHong Zhang 569b54502bSHong Zhang PetscFunctionBegin; 57075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 589b54502bSHong Zhang PetscFunctionReturn(0); 599b54502bSHong Zhang } 609b54502bSHong Zhang 619b54502bSHong Zhang #undef __FUNCT__ 629b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 639b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 649b54502bSHong Zhang { 659b54502bSHong Zhang PetscErrorCode ierr; 669b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 679b54502bSHong Zhang 689b54502bSHong Zhang PetscFunctionBegin; 69075768bcSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 709b54502bSHong Zhang PetscFunctionReturn(0); 719b54502bSHong Zhang } 729b54502bSHong Zhang 739b54502bSHong Zhang #undef __FUNCT__ 749b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 759b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 769b54502bSHong Zhang { 779b54502bSHong Zhang PetscErrorCode ierr; 789b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 799b54502bSHong Zhang 809b54502bSHong Zhang PetscFunctionBegin; 81075768bcSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 829b54502bSHong Zhang PetscFunctionReturn(0); 839b54502bSHong Zhang } 849b54502bSHong Zhang 859b54502bSHong Zhang #undef __FUNCT__ 869b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 879b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 889b54502bSHong Zhang { 899b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 909b54502bSHong Zhang PetscTruth flg; 91*4c9036c7SBarry Smith PetscInt dtmax = 3; 92*4c9036c7SBarry Smith PetscReal dt[3]; 939b54502bSHong Zhang PetscErrorCode ierr; 949b54502bSHong Zhang 959b54502bSHong Zhang PetscFunctionBegin; 969b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 978ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 989b54502bSHong Zhang 998ff23777SHong Zhang ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr); 100*4c9036c7SBarry Smith dt[0] = ((PC_Factor*)icc)->info.dt; 101*4c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 102*4c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 103*4c9036c7SBarry Smith ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 104*4c9036c7SBarry Smith if (flg) { 105*4c9036c7SBarry Smith ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 106*4c9036c7SBarry Smith } 107*4c9036c7SBarry Smith 1089b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1099b54502bSHong Zhang PetscFunctionReturn(0); 1109b54502bSHong Zhang } 1119b54502bSHong Zhang 1129b54502bSHong Zhang #undef __FUNCT__ 1139b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 1149b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 1159b54502bSHong Zhang { 1169b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1179b54502bSHong Zhang PetscErrorCode ierr; 1189b54502bSHong Zhang PetscTruth isstring,iascii; 1199b54502bSHong Zhang 1209b54502bSHong Zhang PetscFunctionBegin; 1219b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 1229b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1239b54502bSHong Zhang if (iascii) { 124075768bcSBarry Smith if (((PC_Factor*)icc)->info.levels == 1) { 125075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: %D level of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr); 1269b54502bSHong Zhang } else { 127075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: %D levels of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr); 1289b54502bSHong Zhang } 129075768bcSBarry 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); 130075768bcSBarry Smith if (((PC_Factor*)icc)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using Manteuffel shift\n");CHKERRQ(ierr);} 131075768bcSBarry Smith if (((PC_Factor*)icc)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer," ICC: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);} 132075768bcSBarry Smith if (((PC_Factor*)icc)->fact) { 133f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr); 134f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 135f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 136f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 137f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 138f3a39becSBarry Smith ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 139075768bcSBarry Smith ierr = MatView(((PC_Factor*)icc)->fact,viewer);CHKERRQ(ierr); 140f3a39becSBarry Smith ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 141f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 142f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 143f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 144f3a39becSBarry Smith } 1459b54502bSHong Zhang } else if (isstring) { 146075768bcSBarry Smith ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);CHKERRQ(ierr); 1479b54502bSHong Zhang } else { 1489b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name); 1499b54502bSHong Zhang } 1509b54502bSHong Zhang PetscFunctionReturn(0); 1519b54502bSHong Zhang } 1529b54502bSHong Zhang 153*4c9036c7SBarry Smith extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 154*4c9036c7SBarry Smith 1559b54502bSHong Zhang /*MC 1569b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1579b54502bSHong Zhang 1589b54502bSHong Zhang Options Database Keys: 1592401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1602401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1619b54502bSHong Zhang its factorization (overwrites original matrix) 16255ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 1632401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 164f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 165f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 166f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 1679b54502bSHong Zhang 1689b54502bSHong Zhang Level: beginner 1699b54502bSHong Zhang 1709b54502bSHong Zhang Concepts: incomplete Cholesky factorization 1719b54502bSHong Zhang 172d9ba8547SBarry Smith Notes: Only implemented for some matrix formats. Not implemented in parallel. 1739b54502bSHong Zhang 1749b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1759b54502bSHong Zhang 1769b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1779b54502bSHong Zhang 1782401956bSBarry Smith By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE); 1799b54502bSHong Zhang to turn off the shift. 1809b54502bSHong Zhang 181c582cd25SBarry Smith References: 182c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 183c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 184c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 185c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 186c582cd25SBarry Smith 1879b54502bSHong Zhang 1889b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 1898ff23777SHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),PCFactorSetShiftInBlocks(), 190e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1918ff23777SHong Zhang PCFactorSetLevels() 1929b54502bSHong Zhang 1939b54502bSHong Zhang M*/ 1949b54502bSHong Zhang 1959b54502bSHong Zhang EXTERN_C_BEGIN 1969b54502bSHong Zhang #undef __FUNCT__ 1979b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 198dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 1999b54502bSHong Zhang { 2009b54502bSHong Zhang PetscErrorCode ierr; 2019b54502bSHong Zhang PC_ICC *icc; 2029b54502bSHong Zhang 2039b54502bSHong Zhang PetscFunctionBegin; 20438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 2059b54502bSHong Zhang 206075768bcSBarry Smith ((PC_Factor*)icc)->fact = 0; 207075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 208174acd27SBarry Smith ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 209075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr); 21075567043SBarry Smith ((PC_Factor*)icc)->info.levels = 0.; 211075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 2129b54502bSHong Zhang icc->implctx = 0; 2139b54502bSHong Zhang 214075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 215075768bcSBarry Smith ((PC_Factor*)icc)->info.shiftnz = 0.0; 216075768bcSBarry Smith ((PC_Factor*)icc)->info.shiftpd = 1.0; /* true */ 217075768bcSBarry Smith ((PC_Factor*)icc)->info.zeropivot = 1.e-12; 2189b54502bSHong Zhang pc->data = (void*)icc; 2199b54502bSHong Zhang 2209b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 2219b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 2229b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 2239b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 2249b54502bSHong Zhang pc->ops->view = PCView_ICC; 22585317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 2269b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 2279b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 2289b54502bSHong Zhang 2297112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 2307112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 23185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 23285317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 23385317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor", 23485317021SBarry Smith PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr); 23585317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor", 23685317021SBarry Smith PCFactorSetShiftPd_Factor);CHKERRQ(ierr); 2378ff23777SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor", 2388ff23777SHong Zhang PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr); 239afaefe49SHong Zhang 24085317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 24185317021SBarry Smith PCFactorSetLevels_Factor);CHKERRQ(ierr); 24285317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 24385317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 24485317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 24585317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 246174acd27SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 247174acd27SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 248*4c9036c7SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 249*4c9036c7SBarry Smith PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 2509b54502bSHong Zhang PetscFunctionReturn(0); 2519b54502bSHong Zhang } 2529b54502bSHong Zhang EXTERN_C_END 2539b54502bSHong Zhang 2549b54502bSHong Zhang 255