1dba47a55SKris Buschelman 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h> /*I "petscpc.h" I*/ 39b54502bSHong Zhang 44cccfbddSHong Zhang static PetscErrorCode PCSetUp_ICC(PC pc) 59b54502bSHong Zhang { 69b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 72c7c0729SBarry Smith IS perm = NULL,cperm = NULL; 89b54502bSHong Zhang PetscErrorCode ierr; 9f3a39becSBarry Smith MatInfo info; 10ea799195SBarry Smith MatSolverType stype; 1100e125f8SBarry Smith MatFactorError err; 129b54502bSHong Zhang 139b54502bSHong Zhang PetscFunctionBegin; 14c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 159b54502bSHong Zhang 1684d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 179b54502bSHong Zhang if (!pc->setupcalled) { 18d09a07f4SBarry Smith if (!((PC_Factor*)icc)->fact) { 19*f73b0415SBarry Smith PetscBool canuseordering; 20174acd27SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 21*f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering);CHKERRQ(ierr); 22*f73b0415SBarry Smith if (canuseordering) { 232c7c0729SBarry Smith ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr); 242c7c0729SBarry Smith } 25a1f19f5aSHong Zhang } 26075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 279b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 28*f73b0415SBarry Smith PetscBool canuseordering; 296bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 3014798fb4SJed Brown ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 31*f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering);CHKERRQ(ierr); 32*f73b0415SBarry Smith if (canuseordering) { 332c7c0729SBarry Smith ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr); 342c7c0729SBarry Smith } 35075768bcSBarry Smith ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 369b54502bSHong Zhang } 37075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 383d1c1ea0SBarry Smith icc->hdr.actualfill = info.fill_ratio_needed; 39f3a39becSBarry Smith 40fcfd50ebSBarry Smith ierr = ISDestroy(&cperm);CHKERRQ(ierr); 41fcfd50ebSBarry Smith ierr = ISDestroy(&perm);CHKERRQ(ierr); 426baea169SHong Zhang 4300e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr); 4400e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 4500e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 466baea169SHong Zhang PetscFunctionReturn(0); 476baea169SHong Zhang } 486baea169SHong Zhang 49075768bcSBarry Smith ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 5000e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr); 5100e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 5200e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 536baea169SHong Zhang } 5400c67f3bSHong Zhang 553ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 5600c67f3bSHong Zhang if (!stype) { 57ea799195SBarry Smith MatSolverType solverpackage; 583ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage);CHKERRQ(ierr); 593ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 6000c67f3bSHong Zhang } 619b54502bSHong Zhang PetscFunctionReturn(0); 629b54502bSHong Zhang } 639b54502bSHong Zhang 64574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc) 65574deadeSBarry Smith { 66574deadeSBarry Smith PC_ICC *icc = (PC_ICC*)pc->data; 67574deadeSBarry Smith PetscErrorCode ierr; 68574deadeSBarry Smith 69574deadeSBarry Smith PetscFunctionBegin; 706bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 71574deadeSBarry Smith PetscFunctionReturn(0); 72574deadeSBarry Smith } 73574deadeSBarry Smith 749b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 759b54502bSHong Zhang { 769b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 779b54502bSHong Zhang PetscErrorCode ierr; 789b54502bSHong Zhang 799b54502bSHong Zhang PetscFunctionBegin; 80574deadeSBarry Smith ierr = PCReset_ICC(pc);CHKERRQ(ierr); 81503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 82503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 83c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 849b54502bSHong Zhang PetscFunctionReturn(0); 859b54502bSHong Zhang } 869b54502bSHong Zhang 879b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 889b54502bSHong Zhang { 899b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 909b54502bSHong Zhang PetscErrorCode ierr; 919b54502bSHong Zhang 929b54502bSHong Zhang PetscFunctionBegin; 93075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 949b54502bSHong Zhang PetscFunctionReturn(0); 959b54502bSHong Zhang } 969b54502bSHong Zhang 977b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y) 987b6e2003SPierre Jolivet { 997b6e2003SPierre Jolivet PC_ICC *icc = (PC_ICC*)pc->data; 1007b6e2003SPierre Jolivet PetscErrorCode ierr; 1017b6e2003SPierre Jolivet 1027b6e2003SPierre Jolivet PetscFunctionBegin; 1037b6e2003SPierre Jolivet ierr = MatMatSolve(((PC_Factor*)icc)->fact,X,Y);CHKERRQ(ierr); 1047b6e2003SPierre Jolivet PetscFunctionReturn(0); 1057b6e2003SPierre Jolivet } 1067b6e2003SPierre Jolivet 1079b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 1089b54502bSHong Zhang { 1099b54502bSHong Zhang PetscErrorCode ierr; 1109b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1119b54502bSHong Zhang 1129b54502bSHong Zhang PetscFunctionBegin; 113075768bcSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 1149b54502bSHong Zhang PetscFunctionReturn(0); 1159b54502bSHong Zhang } 1169b54502bSHong Zhang 1179b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 1189b54502bSHong Zhang { 1199b54502bSHong Zhang PetscErrorCode ierr; 1209b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1219b54502bSHong Zhang 1229b54502bSHong Zhang PetscFunctionBegin; 123075768bcSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 1249b54502bSHong Zhang PetscFunctionReturn(0); 1259b54502bSHong Zhang } 1269b54502bSHong Zhang 1274416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc) 1289b54502bSHong Zhang { 1299b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 130ace3abfcSBarry Smith PetscBool flg; 1319b54502bSHong Zhang PetscErrorCode ierr; 1322fa5cd67SKarl Rupp /* PetscReal dt[3];*/ 1339b54502bSHong Zhang 1349b54502bSHong Zhang PetscFunctionBegin; 135e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr); 136e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 1379b54502bSHong Zhang 1388ff23777SHong Zhang ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr); 1395a586d82SBarry Smith /*dt[0] = ((PC_Factor*)icc)->info.dt; 1404c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 1414c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 14278fc6b22SHong Zhang PetscInt dtmax = 3; 143b7c853c4SBarry Smith ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 1444c9036c7SBarry Smith if (flg) { 145b7c853c4SBarry Smith ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 1464c9036c7SBarry Smith } 14778fc6b22SHong Zhang */ 1489b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1499b54502bSHong Zhang PetscFunctionReturn(0); 1509b54502bSHong Zhang } 1519b54502bSHong Zhang 1527087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 1534c9036c7SBarry Smith 1549b54502bSHong Zhang /*MC 1559b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1569b54502bSHong Zhang 1579b54502bSHong Zhang Options Database Keys: 1582401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1592401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1609b54502bSHong Zhang its factorization (overwrites original matrix) 16155ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 162145b9266SHong Zhang - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 1639b54502bSHong Zhang 1649b54502bSHong Zhang Level: beginner 1659b54502bSHong Zhang 16695452b02SPatrick Sanan Notes: 16795452b02SPatrick Sanan Only implemented for some matrix formats. Not implemented in parallel. 1689b54502bSHong Zhang 1699b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1709b54502bSHong Zhang 1719b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1729b54502bSHong Zhang 17314d2772eSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 1749b54502bSHong Zhang to turn off the shift. 1759b54502bSHong Zhang 176c582cd25SBarry Smith References: 17796a0c994SBarry Smith . 1. - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 17896a0c994SBarry Smith Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 17996a0c994SBarry Smith Science and Engineering, Kluwer. 180c582cd25SBarry Smith 1819b54502bSHong Zhang 1829b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 18381f96dccSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 184e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1858ff23777SHong Zhang PCFactorSetLevels() 1869b54502bSHong Zhang 1879b54502bSHong Zhang M*/ 1889b54502bSHong Zhang 1898cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc) 1909b54502bSHong Zhang { 1919b54502bSHong Zhang PetscErrorCode ierr; 1929b54502bSHong Zhang PC_ICC *icc; 1939b54502bSHong Zhang 1949b54502bSHong Zhang PetscFunctionBegin; 195b00a9115SJed Brown ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr); 1963d1c1ea0SBarry Smith pc->data = (void*)icc; 1973d1c1ea0SBarry Smith ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 19819fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 1992fa5cd67SKarl Rupp 200879e8a4dSBarry Smith ((PC_Factor*)icc)->factortype = MAT_FACTOR_ICC; 201075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 202075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 203f4db908eSBarry Smith ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 2049b54502bSHong Zhang 2059b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 2067b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ICC; 207d674540eSJed Brown pc->ops->applytranspose = PCApply_ICC; 2084cccfbddSHong Zhang pc->ops->setup = PCSetUp_ICC; 209574deadeSBarry Smith pc->ops->reset = PCReset_ICC; 2109b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 2119b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 21292e08861SBarry Smith pc->ops->view = PCView_Factor; 2139b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 2149b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 2159b54502bSHong Zhang PetscFunctionReturn(0); 2169b54502bSHong Zhang } 217