xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 4dbf25a8fa98e38799e7b47dcb2d8a9309975f41)
1c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h> /*I "petscpc.h" I*/
29b54502bSHong Zhang 
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;
599566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)icc)->fact, &solverpackage));
609566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
6100c67f3bSHong Zhang   }
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
639b54502bSHong Zhang }
649b54502bSHong Zhang 
65d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ICC(PC pc)
66d71ae5a4SJacob Faibussowitsch {
67574deadeSBarry Smith   PC_ICC *icc = (PC_ICC *)pc->data;
68574deadeSBarry Smith 
69574deadeSBarry Smith   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72574deadeSBarry Smith }
73574deadeSBarry Smith 
74d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ICC(PC pc)
75d71ae5a4SJacob Faibussowitsch {
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));
822e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
839566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
859b54502bSHong Zhang }
869b54502bSHong Zhang 
87d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ICC(PC pc, Vec x, Vec y)
88d71ae5a4SJacob Faibussowitsch {
899b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
909b54502bSHong Zhang 
919b54502bSHong Zhang   PetscFunctionBegin;
929566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)icc)->fact, x, y));
933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
949b54502bSHong Zhang }
959b54502bSHong Zhang 
96d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ICC(PC pc, Mat X, Mat Y)
97d71ae5a4SJacob Faibussowitsch {
987b6e2003SPierre Jolivet   PC_ICC *icc = (PC_ICC *)pc->data;
997b6e2003SPierre Jolivet 
1007b6e2003SPierre Jolivet   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)icc)->fact, X, Y));
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1037b6e2003SPierre Jolivet }
1047b6e2003SPierre Jolivet 
105d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc, Vec x, Vec y)
106d71ae5a4SJacob Faibussowitsch {
1079b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1089b54502bSHong Zhang 
1099b54502bSHong Zhang   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1129b54502bSHong Zhang }
1139b54502bSHong Zhang 
114d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ICC(PC pc, Vec x, Vec y)
115d71ae5a4SJacob Faibussowitsch {
1169b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1179b54502bSHong Zhang 
1189b54502bSHong Zhang   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1219b54502bSHong Zhang }
1229b54502bSHong Zhang 
123ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PC pc, PetscOptionItems PetscOptionsObject)
124d71ae5a4SJacob Faibussowitsch {
1259b54502bSHong Zhang   PC_ICC   *icc = (PC_ICC *)pc->data;
126ace3abfcSBarry Smith   PetscBool flg;
1272fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1289b54502bSHong Zhang 
1299b54502bSHong Zhang   PetscFunctionBegin;
130d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ICC Options");
131dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
1329b54502bSHong Zhang 
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_factor_levels", "levels of fill", "PCFactorSetLevels", ((PC_Factor *)icc)->info.levels, &((PC_Factor *)icc)->info.levels, &flg));
1345a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1354c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1364c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
13778fc6b22SHong Zhang   PetscInt       dtmax = 3;
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg));
1394c9036c7SBarry Smith   if (flg) {
1409566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]));
1414c9036c7SBarry Smith   }
14278fc6b22SHong Zhang   */
143d0609cedSBarry Smith   PetscOptionsHeadEnd();
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1459b54502bSHong Zhang }
1469b54502bSHong Zhang 
1477087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC, PetscReal, PetscReal, PetscInt);
1484c9036c7SBarry Smith 
1499b54502bSHong Zhang /*MC
1501d27aa22SBarry Smith      PCICC - Incomplete Cholesky factorization preconditioners {cite}`chan1997approximate`
1519b54502bSHong Zhang 
1529b54502bSHong Zhang    Options Database Keys:
1532401956bSBarry Smith +  -pc_factor_levels <k>                                 - number of levels of fill for ICC(k)
1542401956bSBarry Smith .  -pc_factor_in_place                                   - only for ICC(0) with natural ordering, reuses the space of the matrix for
1559b54502bSHong Zhang                                                          its factorization (overwrites original matrix)
15655ba2a51SBarry Smith .  -pc_factor_fill <nfill>                               - expected amount of fill in factored matrix compared to original matrix, nfill > 1
157145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1589b54502bSHong Zhang 
1599b54502bSHong Zhang    Level: beginner
1609b54502bSHong Zhang 
16195452b02SPatrick Sanan    Notes:
16295452b02SPatrick Sanan    Only implemented for some matrix formats. Not implemented in parallel.
1639b54502bSHong Zhang 
164f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ICC.
165f1580f4eSBarry Smith 
1661d27aa22SBarry Smith    By default, the Manteuffel shift {cite}`manteuffel1979shifted` is applied, for matrices with block size 1 only. Call `PCFactorSetShiftType`(pc,`MAT_SHIFT_POSITIVE_DEFINITE`);
167f1580f4eSBarry Smith    to turn off the shift.
1689b54502bSHong Zhang 
169562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCILU`, `PCLU`, `PCCHOLESKY`,
170db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
171db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
172db781477SPatrick Sanan           `PCFactorSetLevels()`
1739b54502bSHong Zhang M*/
1749b54502bSHong Zhang 
175d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
176d71ae5a4SJacob Faibussowitsch {
1779b54502bSHong Zhang   PC_ICC *icc;
1789b54502bSHong Zhang 
1799b54502bSHong Zhang   PetscFunctionBegin;
1804dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&icc));
1813d1c1ea0SBarry Smith   pc->data = (void *)icc;
1829566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
1832fa5cd67SKarl Rupp 
184075768bcSBarry Smith   ((PC_Factor *)icc)->info.fill      = 1.0;
185075768bcSBarry Smith   ((PC_Factor *)icc)->info.dtcol     = PETSC_DEFAULT;
186f4db908eSBarry Smith   ((PC_Factor *)icc)->info.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
1879b54502bSHong Zhang 
1889b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
1897b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
190d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
191*4dbf25a8SPierre Jolivet   pc->ops->matapplytranspose   = PCMatApply_ICC;
1924cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
193574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
1949b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
1959b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
19692e08861SBarry Smith   pc->ops->view                = PCView_Factor;
1979b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
1989b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2009b54502bSHong Zhang }
201