1dba47a55SKris Buschelman 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h> /*I "petscpc.h" I*/ 39b54502bSHong Zhang 49b54502bSHong Zhang #undef __FUNCT__ 59b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC" 69b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc) 79b54502bSHong Zhang { 89b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 99b54502bSHong Zhang IS perm,cperm; 109b54502bSHong Zhang PetscErrorCode ierr; 11f3a39becSBarry Smith MatInfo info; 129b54502bSHong Zhang 139b54502bSHong Zhang PetscFunctionBegin; 14075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr); 159b54502bSHong Zhang 169b54502bSHong Zhang if (!pc->setupcalled) { 17d09a07f4SBarry Smith if (!((PC_Factor*)icc)->fact) { 18174acd27SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 19a1f19f5aSHong Zhang } 20075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 219b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 226bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 2314798fb4SJed Brown ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 24075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 259b54502bSHong Zhang } 26075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 27f3a39becSBarry Smith icc->actualfill = info.fill_ratio_needed; 28f3a39becSBarry Smith 29fcfd50ebSBarry Smith ierr = ISDestroy(&cperm);CHKERRQ(ierr); 30fcfd50ebSBarry Smith ierr = ISDestroy(&perm);CHKERRQ(ierr); 31075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 329b54502bSHong Zhang PetscFunctionReturn(0); 339b54502bSHong Zhang } 349b54502bSHong Zhang 359b54502bSHong Zhang #undef __FUNCT__ 36574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC" 37574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc) 38574deadeSBarry Smith { 39574deadeSBarry Smith PC_ICC *icc = (PC_ICC*)pc->data; 40574deadeSBarry Smith PetscErrorCode ierr; 41574deadeSBarry Smith 42574deadeSBarry Smith PetscFunctionBegin; 436bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 44574deadeSBarry Smith PetscFunctionReturn(0); 45574deadeSBarry Smith } 46574deadeSBarry Smith 47574deadeSBarry Smith #undef __FUNCT__ 489b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC" 499b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 509b54502bSHong Zhang { 519b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 529b54502bSHong Zhang PetscErrorCode ierr; 539b54502bSHong Zhang 549b54502bSHong Zhang PetscFunctionBegin; 55574deadeSBarry Smith ierr = PCReset_ICC(pc);CHKERRQ(ierr); 56503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 57503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 58c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 599b54502bSHong Zhang PetscFunctionReturn(0); 609b54502bSHong Zhang } 619b54502bSHong Zhang 629b54502bSHong Zhang #undef __FUNCT__ 639b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC" 649b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 659b54502bSHong Zhang { 669b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 679b54502bSHong Zhang PetscErrorCode ierr; 689b54502bSHong Zhang 699b54502bSHong Zhang PetscFunctionBegin; 70075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 719b54502bSHong Zhang PetscFunctionReturn(0); 729b54502bSHong Zhang } 739b54502bSHong Zhang 749b54502bSHong Zhang #undef __FUNCT__ 759b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC" 769b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 779b54502bSHong Zhang { 789b54502bSHong Zhang PetscErrorCode ierr; 799b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 809b54502bSHong Zhang 819b54502bSHong Zhang PetscFunctionBegin; 82075768bcSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 839b54502bSHong Zhang PetscFunctionReturn(0); 849b54502bSHong Zhang } 859b54502bSHong Zhang 869b54502bSHong Zhang #undef __FUNCT__ 879b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC" 889b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 899b54502bSHong Zhang { 909b54502bSHong Zhang PetscErrorCode ierr; 919b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 929b54502bSHong Zhang 939b54502bSHong Zhang PetscFunctionBegin; 94075768bcSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 959b54502bSHong Zhang PetscFunctionReturn(0); 969b54502bSHong Zhang } 979b54502bSHong Zhang 989b54502bSHong Zhang #undef __FUNCT__ 999b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC" 1009b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc) 1019b54502bSHong Zhang { 1029b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 103ace3abfcSBarry Smith PetscBool flg; 1049b54502bSHong Zhang PetscErrorCode ierr; 105*2fa5cd67SKarl Rupp /* PetscReal dt[3];*/ 1069b54502bSHong Zhang 1079b54502bSHong Zhang PetscFunctionBegin; 1089b54502bSHong Zhang ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 1098ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 1109b54502bSHong Zhang 1118ff23777SHong Zhang ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr); 1125a586d82SBarry Smith /*dt[0] = ((PC_Factor*)icc)->info.dt; 1134c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 1144c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 11578fc6b22SHong Zhang PetscInt dtmax = 3; 116b7c853c4SBarry Smith ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 1174c9036c7SBarry Smith if (flg) { 118b7c853c4SBarry Smith ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 1194c9036c7SBarry Smith } 12078fc6b22SHong Zhang */ 1219b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1229b54502bSHong Zhang PetscFunctionReturn(0); 1239b54502bSHong Zhang } 1249b54502bSHong Zhang 1259b54502bSHong Zhang #undef __FUNCT__ 1269b54502bSHong Zhang #define __FUNCT__ "PCView_ICC" 1279b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 1289b54502bSHong Zhang { 1299b54502bSHong Zhang PetscErrorCode ierr; 1309b54502bSHong Zhang 1319b54502bSHong Zhang PetscFunctionBegin; 132914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 1339b54502bSHong Zhang PetscFunctionReturn(0); 1349b54502bSHong Zhang } 1359b54502bSHong Zhang 136295a483bSBarry Smith EXTERN_C_BEGIN 1377087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 138295a483bSBarry Smith EXTERN_C_END 1394c9036c7SBarry Smith 1409b54502bSHong Zhang /*MC 1419b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1429b54502bSHong Zhang 1439b54502bSHong Zhang Options Database Keys: 1442401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1452401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1469b54502bSHong Zhang its factorization (overwrites original matrix) 14755ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 148145b9266SHong Zhang - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 1499b54502bSHong Zhang 1509b54502bSHong Zhang Level: beginner 1519b54502bSHong Zhang 1529b54502bSHong Zhang Concepts: incomplete Cholesky factorization 1539b54502bSHong Zhang 154d9ba8547SBarry Smith Notes: Only implemented for some matrix formats. Not implemented in parallel. 1559b54502bSHong Zhang 1569b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1579b54502bSHong Zhang 1589b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1599b54502bSHong Zhang 16014d2772eSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 1619b54502bSHong Zhang to turn off the shift. 1629b54502bSHong Zhang 163c582cd25SBarry Smith References: 164c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 165c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 166c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 167c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 168c582cd25SBarry Smith 1699b54502bSHong Zhang 1709b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 17181f96dccSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 172e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1738ff23777SHong Zhang PCFactorSetLevels() 1749b54502bSHong Zhang 1759b54502bSHong Zhang M*/ 1769b54502bSHong Zhang 1779b54502bSHong Zhang EXTERN_C_BEGIN 1789b54502bSHong Zhang #undef __FUNCT__ 1799b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC" 1807087cfbeSBarry Smith PetscErrorCode PCCreate_ICC(PC pc) 1819b54502bSHong Zhang { 1829b54502bSHong Zhang PetscErrorCode ierr; 1839b54502bSHong Zhang PC_ICC *icc; 1849b54502bSHong Zhang 1859b54502bSHong Zhang PetscFunctionBegin; 18638f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 1879b54502bSHong Zhang 188075768bcSBarry Smith ((PC_Factor*)icc)->fact = 0; 189*2fa5cd67SKarl Rupp 19019fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 1912692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 192075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr); 193*2fa5cd67SKarl Rupp 194879e8a4dSBarry Smith ((PC_Factor*)icc)->factortype = MAT_FACTOR_ICC; 19575567043SBarry Smith ((PC_Factor*)icc)->info.levels = 0.; 196075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 1979b54502bSHong Zhang icc->implctx = 0; 1989b54502bSHong Zhang 199075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 200f4db908eSBarry Smith ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 2010ed735ceSHong Zhang ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 2020ed735ceSHong Zhang ((PC_Factor*)icc)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 2039b54502bSHong Zhang 204915743fcSHong Zhang pc->data = (void*)icc; 2059b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 206d674540eSJed Brown pc->ops->applytranspose = PCApply_ICC; 2079b54502bSHong Zhang pc->ops->setup = PCSetup_ICC; 208574deadeSBarry Smith pc->ops->reset = PCReset_ICC; 2099b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 2109b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 2119b54502bSHong Zhang pc->ops->view = PCView_ICC; 21285317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 2139b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 2149b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 2159b54502bSHong Zhang 216f8260c8fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 217f8260c8fSBarry Smith PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 2187112b564SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 2197112b564SBarry Smith PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 22085317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 22185317021SBarry Smith PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 222d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 223d90ac83dSHong Zhang PCFactorSetShiftType_Factor);CHKERRQ(ierr); 224d90ac83dSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 225d90ac83dSHong Zhang PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 22685317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 22785317021SBarry Smith PCFactorSetLevels_Factor);CHKERRQ(ierr); 22885317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 22985317021SBarry Smith PCFactorSetFill_Factor);CHKERRQ(ierr); 23085317021SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 23185317021SBarry Smith PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 232174acd27SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 233174acd27SBarry Smith PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 234b7c853c4SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU", 235b7c853c4SBarry Smith PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 2369b54502bSHong Zhang PetscFunctionReturn(0); 2379b54502bSHong Zhang } 2389b54502bSHong Zhang EXTERN_C_END 2399b54502bSHong Zhang 2409b54502bSHong Zhang 241