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; 12*f023e1d5SPierre Jolivet const char *prefix; 139b54502bSHong Zhang 149b54502bSHong Zhang PetscFunctionBegin; 15*f023e1d5SPierre Jolivet PetscCall(PCGetOptionsPrefix(pc,&prefix)); 16*f023e1d5SPierre Jolivet PetscCall(MatSetOptionsPrefix(pc->pmat,prefix)); 17c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 189b54502bSHong Zhang 199566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure)); 209b54502bSHong Zhang if (!pc->setupcalled) { 21d09a07f4SBarry Smith if (!((PC_Factor*)icc)->fact) { 229566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact)); 234ac6704cSBarry Smith } 249566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering)); 25f73b0415SBarry Smith if (canuseordering) { 269566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 279566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm)); 282c7c0729SBarry Smith } 299566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info)); 309b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 31f73b0415SBarry Smith PetscBool canuseordering; 329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)icc)->fact)); 339566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact)); 349566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering)); 35f73b0415SBarry Smith if (canuseordering) { 369566063dSJacob Faibussowitsch PetscCall(PCFactorSetDefaultOrdering_Factor(pc)); 379566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm)); 382c7c0729SBarry Smith } 399566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info)); 409b54502bSHong Zhang } 419566063dSJacob Faibussowitsch PetscCall(MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info)); 423d1c1ea0SBarry Smith icc->hdr.actualfill = info.fill_ratio_needed; 43f3a39becSBarry Smith 449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cperm)); 459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 466baea169SHong Zhang 479566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err)); 4800e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 4900e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 506baea169SHong Zhang PetscFunctionReturn(0); 516baea169SHong Zhang } 526baea169SHong Zhang 539566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info)); 549566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err)); 5500e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 5600e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 576baea169SHong Zhang } 5800c67f3bSHong Zhang 599566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatSolverType(pc,&stype)); 6000c67f3bSHong Zhang if (!stype) { 61ea799195SBarry Smith MatSolverType solverpackage; 629566063dSJacob Faibussowitsch PetscCall(MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage)); 639566063dSJacob Faibussowitsch PetscCall(PCFactorSetMatSolverType(pc,solverpackage)); 6400c67f3bSHong Zhang } 659b54502bSHong Zhang PetscFunctionReturn(0); 669b54502bSHong Zhang } 679b54502bSHong Zhang 68574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc) 69574deadeSBarry Smith { 70574deadeSBarry Smith PC_ICC *icc = (PC_ICC*)pc->data; 71574deadeSBarry Smith 72574deadeSBarry Smith PetscFunctionBegin; 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((PC_Factor*)icc)->fact)); 74574deadeSBarry Smith PetscFunctionReturn(0); 75574deadeSBarry Smith } 76574deadeSBarry Smith 779b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc) 789b54502bSHong Zhang { 799b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 809b54502bSHong Zhang 819b54502bSHong Zhang PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(PCReset_ICC(pc)); 839566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)icc)->ordering)); 849566063dSJacob Faibussowitsch PetscCall(PetscFree(((PC_Factor*)icc)->solvertype)); 859566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 869b54502bSHong Zhang PetscFunctionReturn(0); 879b54502bSHong Zhang } 889b54502bSHong Zhang 899b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 909b54502bSHong Zhang { 919b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 929b54502bSHong Zhang 939b54502bSHong Zhang PetscFunctionBegin; 949566063dSJacob Faibussowitsch PetscCall(MatSolve(((PC_Factor*)icc)->fact,x,y)); 959b54502bSHong Zhang PetscFunctionReturn(0); 969b54502bSHong Zhang } 979b54502bSHong Zhang 987b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y) 997b6e2003SPierre Jolivet { 1007b6e2003SPierre Jolivet PC_ICC *icc = (PC_ICC*)pc->data; 1017b6e2003SPierre Jolivet 1027b6e2003SPierre Jolivet PetscFunctionBegin; 1039566063dSJacob Faibussowitsch PetscCall(MatMatSolve(((PC_Factor*)icc)->fact,X,Y)); 1047b6e2003SPierre Jolivet PetscFunctionReturn(0); 1057b6e2003SPierre Jolivet } 1067b6e2003SPierre Jolivet 1079b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 1089b54502bSHong Zhang { 1099b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1109b54502bSHong Zhang 1119b54502bSHong Zhang PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y)); 1139b54502bSHong Zhang PetscFunctionReturn(0); 1149b54502bSHong Zhang } 1159b54502bSHong Zhang 1169b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 1179b54502bSHong Zhang { 1189b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 1199b54502bSHong Zhang 1209b54502bSHong Zhang PetscFunctionBegin; 1219566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y)); 1229b54502bSHong Zhang PetscFunctionReturn(0); 1239b54502bSHong Zhang } 1249b54502bSHong Zhang 1254416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc) 1269b54502bSHong Zhang { 1279b54502bSHong Zhang PC_ICC *icc = (PC_ICC*)pc->data; 128ace3abfcSBarry Smith PetscBool flg; 1292fa5cd67SKarl Rupp /* PetscReal dt[3];*/ 1309b54502bSHong Zhang 1319b54502bSHong Zhang PetscFunctionBegin; 132d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"ICC Options"); 1339566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc)); 1349b54502bSHong Zhang 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg)); 1365a586d82SBarry Smith /*dt[0] = ((PC_Factor*)icc)->info.dt; 1374c9036c7SBarry Smith dt[1] = ((PC_Factor*)icc)->info.dtcol; 1384c9036c7SBarry Smith dt[2] = ((PC_Factor*)icc)->info.dtcount; 13978fc6b22SHong Zhang PetscInt dtmax = 3; 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg)); 1414c9036c7SBarry Smith if (flg) { 1429566063dSJacob Faibussowitsch PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2])); 1434c9036c7SBarry Smith } 14478fc6b22SHong Zhang */ 145d0609cedSBarry Smith PetscOptionsHeadEnd(); 1469b54502bSHong Zhang PetscFunctionReturn(0); 1479b54502bSHong Zhang } 1489b54502bSHong Zhang 1497087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 1504c9036c7SBarry Smith 1519b54502bSHong Zhang /*MC 1529b54502bSHong Zhang PCICC - Incomplete Cholesky factorization preconditioners. 1539b54502bSHong Zhang 1549b54502bSHong Zhang Options Database Keys: 1552401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ICC(k) 1562401956bSBarry Smith . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 1579b54502bSHong Zhang its factorization (overwrites original matrix) 15855ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 159145b9266SHong Zhang - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 1609b54502bSHong Zhang 1619b54502bSHong Zhang Level: beginner 1629b54502bSHong Zhang 16395452b02SPatrick Sanan Notes: 16495452b02SPatrick Sanan Only implemented for some matrix formats. Not implemented in parallel. 1659b54502bSHong Zhang 1669b54502bSHong Zhang For BAIJ matrices this implements a point block ICC. 1679b54502bSHong Zhang 1689b54502bSHong Zhang The Manteuffel shift is only implemented for matrices with block size 1 1699b54502bSHong Zhang 17014d2772eSHong Zhang By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 1719b54502bSHong Zhang to turn off the shift. 1729b54502bSHong Zhang 173c582cd25SBarry Smith References: 174606c0280SSatish Balay . * - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 17596a0c994SBarry Smith Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 17696a0c994SBarry Smith Science and Engineering, Kluwer. 177c582cd25SBarry Smith 1789b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 17981f96dccSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 180e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 1818ff23777SHong Zhang PCFactorSetLevels() 1829b54502bSHong Zhang 1839b54502bSHong Zhang M*/ 1849b54502bSHong Zhang 1858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc) 1869b54502bSHong Zhang { 1879b54502bSHong Zhang PC_ICC *icc; 1889b54502bSHong Zhang 1899b54502bSHong Zhang PetscFunctionBegin; 1909566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&icc)); 1913d1c1ea0SBarry Smith pc->data = (void*)icc; 1929566063dSJacob Faibussowitsch PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC)); 1932fa5cd67SKarl Rupp 194075768bcSBarry Smith ((PC_Factor*)icc)->info.fill = 1.0; 195075768bcSBarry Smith ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 196f4db908eSBarry Smith ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 1979b54502bSHong Zhang 1989b54502bSHong Zhang pc->ops->apply = PCApply_ICC; 1997b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ICC; 200d674540eSJed Brown pc->ops->applytranspose = PCApply_ICC; 2014cccfbddSHong Zhang pc->ops->setup = PCSetUp_ICC; 202574deadeSBarry Smith pc->ops->reset = PCReset_ICC; 2039b54502bSHong Zhang pc->ops->destroy = PCDestroy_ICC; 2049b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ICC; 20592e08861SBarry Smith pc->ops->view = PCView_Factor; 2069b54502bSHong Zhang pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 2079b54502bSHong Zhang pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 2089b54502bSHong Zhang PetscFunctionReturn(0); 2099b54502bSHong Zhang } 210