xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
1c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h> /*I "petscpc.h" I*/
29b54502bSHong Zhang 
PCSetUp_ICC(PC pc)3d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ICC(PC pc)
4d71ae5a4SJacob Faibussowitsch {
59b54502bSHong Zhang   PC_ICC        *icc  = (PC_ICC *)pc->data;
62c7c0729SBarry Smith   IS             perm = NULL, cperm = NULL;
7f3a39becSBarry Smith   MatInfo        info;
8ea799195SBarry Smith   MatSolverType  stype;
900e125f8SBarry Smith   MatFactorError err;
104ac6704cSBarry Smith   PetscBool      canuseordering;
11f023e1d5SPierre Jolivet   const char    *prefix;
129b54502bSHong Zhang 
139b54502bSHong Zhang   PetscFunctionBegin;
14c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
159b54502bSHong Zhang 
1626cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
1726cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
1826cc229bSBarry Smith 
199566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
209b54502bSHong Zhang   if (!pc->setupcalled) {
2103e5aca4SStefano Zampini     PetscCall(PCFactorSetUpMatSolverType(pc));
229566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)icc)->fact, &canuseordering));
23f73b0415SBarry Smith     if (canuseordering) {
249566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
259566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)icc)->ordering, &perm, &cperm));
262c7c0729SBarry Smith     }
279566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor *)icc)->fact, pc->pmat, perm, &((PC_Factor *)icc)->info));
289b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
299566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
3003e5aca4SStefano Zampini     PetscCall(PCFactorSetUpMatSolverType(pc));
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;
473ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
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;
59*c294e005SBarry Smith 
609566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)icc)->fact, &solverpackage));
619566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
6200c67f3bSHong Zhang   }
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
649b54502bSHong Zhang }
659b54502bSHong Zhang 
PCReset_ICC(PC pc)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ICC(PC pc)
67d71ae5a4SJacob Faibussowitsch {
68574deadeSBarry Smith   PC_ICC *icc = (PC_ICC *)pc->data;
69574deadeSBarry Smith 
70574deadeSBarry Smith   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73574deadeSBarry Smith }
74574deadeSBarry Smith 
PCDestroy_ICC(PC pc)75d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ICC(PC pc)
76d71ae5a4SJacob Faibussowitsch {
779b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
789b54502bSHong Zhang 
799b54502bSHong Zhang   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(PCReset_ICC(pc));
819566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)icc)->ordering));
829566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)icc)->solvertype));
832e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
849566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
869b54502bSHong Zhang }
879b54502bSHong Zhang 
PCApply_ICC(PC pc,Vec x,Vec y)88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ICC(PC pc, Vec x, Vec y)
89d71ae5a4SJacob Faibussowitsch {
909b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
919b54502bSHong Zhang 
929b54502bSHong Zhang   PetscFunctionBegin;
939566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)icc)->fact, x, y));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
959b54502bSHong Zhang }
969b54502bSHong Zhang 
PCMatApply_ICC(PC pc,Mat X,Mat Y)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ICC(PC pc, Mat X, Mat Y)
98d71ae5a4SJacob Faibussowitsch {
997b6e2003SPierre Jolivet   PC_ICC *icc = (PC_ICC *)pc->data;
1007b6e2003SPierre Jolivet 
1017b6e2003SPierre Jolivet   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)icc)->fact, X, Y));
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1047b6e2003SPierre Jolivet }
1057b6e2003SPierre Jolivet 
PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)106d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc, Vec x, Vec y)
107d71ae5a4SJacob Faibussowitsch {
1089b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1099b54502bSHong Zhang 
1109b54502bSHong Zhang   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1139b54502bSHong Zhang }
1149b54502bSHong Zhang 
PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)115d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ICC(PC pc, Vec x, Vec y)
116d71ae5a4SJacob Faibussowitsch {
1179b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1189b54502bSHong Zhang 
1199b54502bSHong Zhang   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1229b54502bSHong Zhang }
1239b54502bSHong Zhang 
PCSetFromOptions_ICC(PC pc,PetscOptionItems PetscOptionsObject)124ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PC pc, PetscOptionItems PetscOptionsObject)
125d71ae5a4SJacob Faibussowitsch {
1269b54502bSHong Zhang   PC_ICC   *icc = (PC_ICC *)pc->data;
127ace3abfcSBarry Smith   PetscBool flg;
1282fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1299b54502bSHong Zhang 
1309b54502bSHong Zhang   PetscFunctionBegin;
131d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ICC Options");
132dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
1339b54502bSHong Zhang 
1349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_levels", "levels of fill", "PCFactorSetLevels", ((PC_Factor *)icc)->info.levels, &((PC_Factor *)icc)->info.levels, &flg));
1355a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1364c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1374c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
13878fc6b22SHong Zhang   PetscInt       dtmax = 3;
1399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg));
1403a7d0413SPierre Jolivet   if (flg) PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]));
14178fc6b22SHong Zhang   */
142d0609cedSBarry Smith   PetscOptionsHeadEnd();
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1449b54502bSHong Zhang }
1459b54502bSHong Zhang 
1469b54502bSHong Zhang /*MC
1471d27aa22SBarry Smith      PCICC - Incomplete Cholesky factorization preconditioners {cite}`chan1997approximate`
1489b54502bSHong Zhang 
1499b54502bSHong Zhang    Options Database Keys:
1502401956bSBarry Smith +  -pc_factor_levels <k>                                 - number of levels of fill for ICC(k)
1512401956bSBarry Smith .  -pc_factor_in_place                                   - only for ICC(0) with natural ordering, reuses the space of the matrix for
1529b54502bSHong Zhang                                                          its factorization (overwrites original matrix)
15355ba2a51SBarry Smith .  -pc_factor_fill <nfill>                               - expected amount of fill in factored matrix compared to original matrix, nfill > 1
154145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1559b54502bSHong Zhang 
1569b54502bSHong Zhang    Level: beginner
1579b54502bSHong Zhang 
15895452b02SPatrick Sanan    Notes:
15995452b02SPatrick Sanan    Only implemented for some matrix formats. Not implemented in parallel.
1609b54502bSHong Zhang 
161f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ICC.
162f1580f4eSBarry Smith 
1631d27aa22SBarry Smith    By default, the Manteuffel shift {cite}`manteuffel1979shifted` is applied, for matrices with block size 1 only. Call `PCFactorSetShiftType`(pc,`MAT_SHIFT_POSITIVE_DEFINITE`);
164f1580f4eSBarry Smith    to turn off the shift.
1659b54502bSHong Zhang 
166562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCILU`, `PCLU`, `PCCHOLESKY`,
167db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
168db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
169db781477SPatrick Sanan           `PCFactorSetLevels()`
1709b54502bSHong Zhang M*/
1719b54502bSHong Zhang 
PCCreate_ICC(PC pc)172d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
173d71ae5a4SJacob Faibussowitsch {
1749b54502bSHong Zhang   PC_ICC *icc;
1759b54502bSHong Zhang 
1769b54502bSHong Zhang   PetscFunctionBegin;
1774dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&icc));
1783d1c1ea0SBarry Smith   pc->data = (void *)icc;
1799566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
1802fa5cd67SKarl Rupp 
181075768bcSBarry Smith   ((PC_Factor *)icc)->info.fill      = 1.0;
182075768bcSBarry Smith   ((PC_Factor *)icc)->info.dtcol     = PETSC_DEFAULT;
183f4db908eSBarry Smith   ((PC_Factor *)icc)->info.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
1849b54502bSHong Zhang 
1859b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
1867b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
187d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
1884dbf25a8SPierre Jolivet   pc->ops->matapplytranspose   = PCMatApply_ICC;
1894cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
190574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
1919b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
1929b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
19392e08861SBarry Smith   pc->ops->view                = PCView_Factor;
1949b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
1959b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1979b54502bSHong Zhang }
198