xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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;
12f023e1d5SPierre Jolivet   const char             *prefix;
139b54502bSHong Zhang 
149b54502bSHong Zhang   PetscFunctionBegin;
15a7f29c02SPierre Jolivet   if (!((PetscObject)pc->pmat)->prefix) {
16f023e1d5SPierre Jolivet     PetscCall(PCGetOptionsPrefix(pc,&prefix));
17f023e1d5SPierre Jolivet     PetscCall(MatSetOptionsPrefix(pc->pmat,prefix));
18a7f29c02SPierre Jolivet   }
19c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
209b54502bSHong Zhang 
219566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat,pc->erroriffailure));
229b54502bSHong Zhang   if (!pc->setupcalled) {
23d09a07f4SBarry Smith     if (!((PC_Factor*)icc)->fact) {
249566063dSJacob Faibussowitsch       PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact));
254ac6704cSBarry Smith     }
269566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering));
27f73b0415SBarry Smith     if (canuseordering) {
289566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
299566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm));
302c7c0729SBarry Smith     }
319566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info));
329b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
33f73b0415SBarry Smith     PetscBool canuseordering;
349566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&((PC_Factor*)icc)->fact));
359566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact));
369566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering));
37f73b0415SBarry Smith     if (canuseordering) {
389566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
399566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm));
402c7c0729SBarry Smith     }
419566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info));
429b54502bSHong Zhang   }
439566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info));
443d1c1ea0SBarry Smith   icc->hdr.actualfill = info.fill_ratio_needed;
45f3a39becSBarry Smith 
469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
479566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
486baea169SHong Zhang 
499566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err));
5000e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
5100e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
526baea169SHong Zhang     PetscFunctionReturn(0);
536baea169SHong Zhang   }
546baea169SHong Zhang 
559566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info));
569566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor*)icc)->fact,&err));
5700e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
5800e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
596baea169SHong Zhang   }
6000c67f3bSHong Zhang 
619566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc,&stype));
6200c67f3bSHong Zhang   if (!stype) {
63ea799195SBarry Smith     MatSolverType solverpackage;
649566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage));
659566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc,solverpackage));
6600c67f3bSHong Zhang   }
679b54502bSHong Zhang   PetscFunctionReturn(0);
689b54502bSHong Zhang }
699b54502bSHong Zhang 
70574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
71574deadeSBarry Smith {
72574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
73574deadeSBarry Smith 
74574deadeSBarry Smith   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&((PC_Factor*)icc)->fact));
76574deadeSBarry Smith   PetscFunctionReturn(0);
77574deadeSBarry Smith }
78574deadeSBarry Smith 
799b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
809b54502bSHong Zhang {
819b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
829b54502bSHong Zhang 
839b54502bSHong Zhang   PetscFunctionBegin;
849566063dSJacob Faibussowitsch   PetscCall(PCReset_ICC(pc));
859566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)icc)->ordering));
869566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor*)icc)->solvertype));
87*2e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
889566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
899b54502bSHong Zhang   PetscFunctionReturn(0);
909b54502bSHong Zhang }
919b54502bSHong Zhang 
929b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
939b54502bSHong Zhang {
949b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
959b54502bSHong Zhang 
969b54502bSHong Zhang   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor*)icc)->fact,x,y));
989b54502bSHong Zhang   PetscFunctionReturn(0);
999b54502bSHong Zhang }
1009b54502bSHong Zhang 
1017b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y)
1027b6e2003SPierre Jolivet {
1037b6e2003SPierre Jolivet   PC_ICC         *icc = (PC_ICC*)pc->data;
1047b6e2003SPierre Jolivet 
1057b6e2003SPierre Jolivet   PetscFunctionBegin;
1069566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor*)icc)->fact,X,Y));
1077b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1087b6e2003SPierre Jolivet }
1097b6e2003SPierre Jolivet 
1109b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
1119b54502bSHong Zhang {
1129b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1139b54502bSHong Zhang 
1149b54502bSHong Zhang   PetscFunctionBegin;
1159566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor*)icc)->fact,x,y));
1169b54502bSHong Zhang   PetscFunctionReturn(0);
1179b54502bSHong Zhang }
1189b54502bSHong Zhang 
1199b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1209b54502bSHong Zhang {
1219b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1229b54502bSHong Zhang 
1239b54502bSHong Zhang   PetscFunctionBegin;
1249566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor*)icc)->fact,x,y));
1259b54502bSHong Zhang   PetscFunctionReturn(0);
1269b54502bSHong Zhang }
1279b54502bSHong Zhang 
1284416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1299b54502bSHong Zhang {
1309b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
131ace3abfcSBarry Smith   PetscBool      flg;
1322fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1339b54502bSHong Zhang 
1349b54502bSHong Zhang   PetscFunctionBegin;
135d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"ICC Options");
1369566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions_Factor(PetscOptionsObject,pc));
1379b54502bSHong Zhang 
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg));
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;
1439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg));
1444c9036c7SBarry Smith   if (flg) {
1459566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]));
1464c9036c7SBarry Smith   }
14778fc6b22SHong Zhang   */
148d0609cedSBarry Smith   PetscOptionsHeadEnd();
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:
177606c0280SSatish Balay .  * - 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 
181db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`,
182db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
183db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
184db781477SPatrick Sanan           `PCFactorSetLevels()`
1859b54502bSHong Zhang 
1869b54502bSHong Zhang M*/
1879b54502bSHong Zhang 
1888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
1899b54502bSHong Zhang {
1909b54502bSHong Zhang   PC_ICC         *icc;
1919b54502bSHong Zhang 
1929b54502bSHong Zhang   PetscFunctionBegin;
1939566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&icc));
1943d1c1ea0SBarry Smith   pc->data = (void*)icc;
1959566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
1962fa5cd67SKarl Rupp 
197075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill      = 1.0;
198075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol     = PETSC_DEFAULT;
199f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2009b54502bSHong Zhang 
2019b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
2027b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
203d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2044cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
205574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2069b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2079b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
20892e08861SBarry Smith   pc->ops->view                = PCView_Factor;
2099b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2109b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2119b54502bSHong Zhang   PetscFunctionReturn(0);
2129b54502bSHong Zhang }
213