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; 8f3a39becSBarry Smith MatInfo info; 9ea799195SBarry Smith MatSolverType stype; 1000e125f8SBarry Smith MatFactorError err; 114ac6704cSBarry Smith PetscBool canuseordering; 129b54502bSHong Zhang 139b54502bSHong Zhang PetscFunctionBegin; 14c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 159b54502bSHong Zhang 169566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure)); 179b54502bSHong Zhang if (!pc->setupcalled) { 18d09a07f4SBarry Smith if (!((PC_Factor*)icc)->fact) { 199566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact)); 204ac6704cSBarry Smith } 219566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering)); 22f73b0415SBarry Smith if (canuseordering) { 239566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 249566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm)); 252c7c0729SBarry Smith } 269566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info)); 279b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 28f73b0415SBarry Smith PetscBool canuseordering; 299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)icc)->fact)); 309566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact)); 319566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering)); 32f73b0415SBarry Smith if (canuseordering) { 339566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 349566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm)); 352c7c0729SBarry Smith } 369566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info)); 379b54502bSHong Zhang } 389566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info)); 393d1c1ea0SBarry Smith icc->hdr.actualfill = info.fill_ratio_needed; 40f3a39becSBarry Smith 419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 436baea169SHong Zhang 449566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err)); 4500e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 4600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 476baea169SHong Zhang PetscFunctionReturn(0); 486baea169SHong Zhang } 496baea169SHong Zhang 509566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info)); 519566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err)); 5200e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 5300e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 546baea169SHong Zhang } 5500c67f3bSHong Zhang 569566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 5700c67f3bSHong Zhang if (!stype) { 58ea799195SBarry Smith MatSolverType solverpackage; 599566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage)); 609566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 6100c67f3bSHong Zhang } 629b54502bSHong Zhang PetscFunctionReturn(0); 639b54502bSHong Zhang } 649b54502bSHong Zhang 65574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc) 66574deadeSBarry Smith { 67574deadeSBarry Smith PC_ICC *icc = (PC_ICC*)pc->data; 68574deadeSBarry Smith 69574deadeSBarry Smith PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)icc)->fact)); 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 789b54502bSHong Zhang PetscFunctionBegin; 799566063dSJacob Faibussowitsch PetscCall(PCReset_ICC(pc)); 809566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)icc)->ordering)); 819566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)icc)->solvertype)); 829566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 839b54502bSHong Zhang PetscFunctionReturn(0); 849b54502bSHong Zhang } 859b54502bSHong Zhang 869b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 879b54502bSHong Zhang { 889b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 899b54502bSHong Zhang 909b54502bSHong Zhang PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor*)icc)->fact,x,y)); 929b54502bSHong Zhang PetscFunctionReturn(0); 939b54502bSHong Zhang } 949b54502bSHong Zhang 957b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y) 967b6e2003SPierre Jolivet { 977b6e2003SPierre Jolivet PC_ICC *icc = (PC_ICC*)pc->data; 987b6e2003SPierre Jolivet 997b6e2003SPierre Jolivet PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor*)icc)->fact,X,Y)); 1017b6e2003SPierre Jolivet PetscFunctionReturn(0); 1027b6e2003SPierre Jolivet } 1037b6e2003SPierre Jolivet 1049b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 1059b54502bSHong Zhang { 1069b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1079b54502bSHong Zhang 1089b54502bSHong Zhang PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y)); 1109b54502bSHong Zhang PetscFunctionReturn(0); 1119b54502bSHong Zhang } 1129b54502bSHong Zhang 1139b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 1149b54502bSHong Zhang { 1159b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1169b54502bSHong Zhang 1179b54502bSHong Zhang PetscFunctionBegin; 1189566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y)); 1199b54502bSHong Zhang PetscFunctionReturn(0); 1209b54502bSHong Zhang } 1219b54502bSHong Zhang 1224416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc) 1239b54502bSHong Zhang { 1249b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 125ace3abfcSBarry Smith PetscBool flg; 1262fa5cd67SKarl Rupp /* PetscReal dt[3];*/ 1279b54502bSHong Zhang 1289b54502bSHong Zhang PetscFunctionBegin; 129*d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"ICC Options"); 1309566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc)); 1319b54502bSHong Zhang 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg)); 1335a586d82SBarry Smith /*dt[0] = ((PC_Factor*)icc)->info.dt; 1344c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 1354c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 13678fc6b22SHong Zhang PetscInt dtmax = 3; 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg)); 1384c9036c7SBarry Smith if (flg) { 1399566063dSJacob Faibussowitsch PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2])); 1404c9036c7SBarry Smith } 14178fc6b22SHong Zhang */ 142*d0609cedSBarry Smith PetscOptionsHeadEnd(); 1439b54502bSHong Zhang PetscFunctionReturn(0); 1449b54502bSHong Zhang } 1459b54502bSHong Zhang 1467087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 1474c9036c7SBarry Smith 1489b54502bSHong Zhang /*MC 1499b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1509b54502bSHong Zhang 1519b54502bSHong Zhang Options Database Keys: 1522401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1532401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1549b54502bSHong Zhang its factorization (overwrites original matrix) 15555ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 156145b9266SHong Zhang - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 1579b54502bSHong Zhang 1589b54502bSHong Zhang Level: beginner 1599b54502bSHong Zhang 16095452b02SPatrick Sanan Notes: 16195452b02SPatrick Sanan Only implemented for some matrix formats. Not implemented in parallel. 1629b54502bSHong Zhang 1639b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1649b54502bSHong Zhang 1659b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1669b54502bSHong Zhang 16714d2772eSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 1689b54502bSHong Zhang to turn off the shift. 1699b54502bSHong Zhang 170c582cd25SBarry Smith References: 171606c0280SSatish Balay . * - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 17296a0c994SBarry Smith Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 17396a0c994SBarry Smith Science and Engineering, Kluwer. 174c582cd25SBarry Smith 1759b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 17681f96dccSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 177e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1788ff23777SHong Zhang PCFactorSetLevels() 1799b54502bSHong Zhang 1809b54502bSHong Zhang M*/ 1819b54502bSHong Zhang 1828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc) 1839b54502bSHong Zhang { 1849b54502bSHong Zhang PC_ICC *icc; 1859b54502bSHong Zhang 1869b54502bSHong Zhang PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&icc)); 1883d1c1ea0SBarry Smith pc->data = (void*)icc; 1899566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC)); 1902fa5cd67SKarl Rupp 191075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 192075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 193f4db908eSBarry Smith ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 1949b54502bSHong Zhang 1959b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 1967b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ICC; 197d674540eSJed Brown pc->ops->applytranspose = PCApply_ICC; 1984cccfbddSHong Zhang pc->ops->setup = PCSetUp_ICC; 199574deadeSBarry Smith pc->ops->reset = PCReset_ICC; 2009b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 2019b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 20292e08861SBarry Smith pc->ops->view = PCView_Factor; 2039b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 2049b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 2059b54502bSHong Zhang PetscFunctionReturn(0); 2069b54502bSHong Zhang } 207