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); 222692d6eeSBarry Smith ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,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);} 43503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 44503cfb0cSBarry Smith ierr = PetscFree(((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*5a586d82SBarry Smith /* PetscReal dt[3];*/ 929b54502bSHong Zhang PetscErrorCode ierr; 939b54502bSHong Zhang 949b54502bSHong Zhang PetscFunctionBegin; 959b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 968ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 979b54502bSHong Zhang 988ff23777SHong Zhang ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr); 99*5a586d82SBarry Smith /*dt[0] = ((PC_Factor*)icc)->info.dt; 1004c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 1014c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 10278fc6b22SHong Zhang PetscInt dtmax = 3; 103b7c853c4SBarry Smith ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 1044c9036c7SBarry Smith if (flg) { 105b7c853c4SBarry Smith ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 1064c9036c7SBarry Smith } 10778fc6b22SHong Zhang */ 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 PetscErrorCode ierr; 1179b54502bSHong Zhang 1189b54502bSHong Zhang PetscFunctionBegin; 119914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 1209b54502bSHong Zhang PetscFunctionReturn(0); 1219b54502bSHong Zhang } 1229b54502bSHong Zhang 123295a483bSBarry Smith EXTERN_C_BEGIN 124b7c853c4SBarry Smith extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 125295a483bSBarry Smith EXTERN_C_END 1264c9036c7SBarry Smith 1279b54502bSHong Zhang /*MC 1289b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1299b54502bSHong Zhang 1309b54502bSHong Zhang Options Database Keys: 1312401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1322401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1339b54502bSHong Zhang its factorization (overwrites original matrix) 13455ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 135145b9266SHong Zhang - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 1369b54502bSHong Zhang 1379b54502bSHong Zhang Level: beginner 1389b54502bSHong Zhang 1399b54502bSHong Zhang Concepts: incomplete Cholesky factorization 1409b54502bSHong Zhang 141d9ba8547SBarry Smith Notes: Only implemented for some matrix formats. Not implemented in parallel. 1429b54502bSHong Zhang 1439b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1449b54502bSHong Zhang 1459b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1469b54502bSHong Zhang 14714d2772eSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 1489b54502bSHong Zhang to turn off the shift. 1499b54502bSHong Zhang 150c582cd25SBarry Smith References: 151c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 152c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 153c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 154c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 155c582cd25SBarry Smith 1569b54502bSHong Zhang 1579b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 15881f96dccSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 159e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1608ff23777SHong Zhang PCFactorSetLevels() 1619b54502bSHong Zhang 1629b54502bSHong Zhang M*/ 1639b54502bSHong Zhang 1649b54502bSHong Zhang EXTERN_C_BEGIN 1659b54502bSHong Zhang #undef __FUNCT__ 1669b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 167dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 1689b54502bSHong Zhang { 1699b54502bSHong Zhang PetscErrorCode ierr; 1709b54502bSHong Zhang PC_ICC *icc; 1719b54502bSHong Zhang 1729b54502bSHong Zhang PetscFunctionBegin; 17338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 1749b54502bSHong Zhang 175075768bcSBarry Smith ((PC_Factor*)icc)->fact = 0; 1762692d6eeSBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 1772692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 178075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr); 179879e8a4dSBarry Smith ((PC_Factor*)icc)->factortype = MAT_FACTOR_ICC; 18075567043SBarry Smith ((PC_Factor*)icc)->info.levels = 0.; 181075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 1829b54502bSHong Zhang icc->implctx = 0; 1839b54502bSHong Zhang 184075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 185f4db908eSBarry Smith ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 186daa17b54SHong Zhang ((PC_Factor*)icc)->info.shiftamount = 1.e-12; 187075768bcSBarry Smith ((PC_Factor*)icc)->info.zeropivot = 1.e-12; 1889b54502bSHong Zhang 189915743fcSHong Zhang pc->data = (void*)icc; 1909b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 1919b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 1929b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 1939b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 1949b54502bSHong Zhang pc->ops->view = PCView_ICC; 19585317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 1969b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 1979b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 1989b54502bSHong Zhang 1997112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 2007112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 20185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 20285317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 203d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 204d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 205d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 206d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 20785317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 20885317021SBarry Smith PCFactorSetLevels_Factor);CHKERRQ(ierr); 20985317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 21085317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 21185317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 21285317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 213174acd27SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 214174acd27SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 215b7c853c4SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU", 216b7c853c4SBarry Smith PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 2179b54502bSHong Zhang PetscFunctionReturn(0); 2189b54502bSHong Zhang } 2199b54502bSHong Zhang EXTERN_C_END 2209b54502bSHong Zhang 2219b54502bSHong Zhang 222