xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 03e5aca47aa9bd9d5a48c628e4311d43ab2119ff)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h> /*I "petscpc.h" I*/
39b54502bSHong Zhang 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ICC(PC pc)
5d71ae5a4SJacob Faibussowitsch {
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;
15c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
169b54502bSHong Zhang 
1726cc229bSBarry Smith   PetscCall(PCGetOptionsPrefix(pc, &prefix));
1826cc229bSBarry Smith   PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix));
1926cc229bSBarry Smith 
209566063dSJacob Faibussowitsch   PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure));
219b54502bSHong Zhang   if (!pc->setupcalled) {
22*03e5aca4SStefano Zampini     PetscCall(PCFactorSetUpMatSolverType(pc));
239566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)icc)->fact, &canuseordering));
24f73b0415SBarry Smith     if (canuseordering) {
259566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
269566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)icc)->ordering, &perm, &cperm));
272c7c0729SBarry Smith     }
289566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor *)icc)->fact, pc->pmat, perm, &((PC_Factor *)icc)->info));
299b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
309566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
31*03e5aca4SStefano Zampini     PetscCall(PCFactorSetUpMatSolverType(pc));
329566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)icc)->fact, &canuseordering));
33f73b0415SBarry Smith     if (canuseordering) {
349566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
359566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)icc)->ordering, &perm, &cperm));
362c7c0729SBarry Smith     }
379566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor *)icc)->fact, pc->pmat, perm, &((PC_Factor *)icc)->info));
389b54502bSHong Zhang   }
399566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(((PC_Factor *)icc)->fact, MAT_LOCAL, &info));
403d1c1ea0SBarry Smith   icc->hdr.actualfill = info.fill_ratio_needed;
41f3a39becSBarry Smith 
429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
446baea169SHong Zhang 
459566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor *)icc)->fact, &err));
4600e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
4700e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
483ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
496baea169SHong Zhang   }
506baea169SHong Zhang 
519566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)icc)->fact, pc->pmat, &((PC_Factor *)icc)->info));
529566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor *)icc)->fact, &err));
5300e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
5400e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
556baea169SHong Zhang   }
5600c67f3bSHong Zhang 
579566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
5800c67f3bSHong Zhang   if (!stype) {
59ea799195SBarry Smith     MatSolverType solverpackage;
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 
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 
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 
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 
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 
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 
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 
124d71ae5a4SJacob Faibussowitsch 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));
1404c9036c7SBarry Smith   if (flg) {
1419566063dSJacob Faibussowitsch     PetscCall(PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]));
1424c9036c7SBarry Smith   }
14378fc6b22SHong Zhang   */
144d0609cedSBarry Smith   PetscOptionsHeadEnd();
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1469b54502bSHong Zhang }
1479b54502bSHong Zhang 
1487087cfbeSBarry Smith extern PetscErrorCode PCFactorSetDropTolerance_ILU(PC, PetscReal, PetscReal, PetscInt);
1494c9036c7SBarry Smith 
1509b54502bSHong Zhang /*MC
1519b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1529b54502bSHong Zhang 
1539b54502bSHong Zhang    Options Database Keys:
1542401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1552401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1569b54502bSHong Zhang                       its factorization (overwrites original matrix)
15755ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
158145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1599b54502bSHong Zhang 
1609b54502bSHong Zhang    Level: beginner
1619b54502bSHong Zhang 
16295452b02SPatrick Sanan    Notes:
16395452b02SPatrick Sanan    Only implemented for some matrix formats. Not implemented in parallel.
1649b54502bSHong Zhang 
165f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ICC.
166f1580f4eSBarry Smith 
167f1580f4eSBarry Smith    By default, the Manteuffel is applied (for matrices with block size 1). Call `PCFactorSetShiftType`(pc,`MAT_SHIFT_POSITIVE_DEFINITE`);
168f1580f4eSBarry Smith    to turn off the shift.
1699b54502bSHong Zhang 
1709b54502bSHong Zhang    The Manteuffel shift is only implemented for matrices with block size 1
1719b54502bSHong Zhang 
172c582cd25SBarry Smith    References:
173606c0280SSatish Balay .  * - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
17496a0c994SBarry Smith       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
17596a0c994SBarry Smith       Science and Engineering, Kluwer.
176c582cd25SBarry Smith 
177f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCILU`, `PCLU`, `PCCHOLESKY`,
178db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
179db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
180db781477SPatrick Sanan           `PCFactorSetLevels()`
1819b54502bSHong Zhang M*/
1829b54502bSHong Zhang 
183d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
184d71ae5a4SJacob Faibussowitsch {
1859b54502bSHong Zhang   PC_ICC *icc;
1869b54502bSHong Zhang 
1879b54502bSHong Zhang   PetscFunctionBegin;
1884dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&icc));
1893d1c1ea0SBarry Smith   pc->data = (void *)icc;
1909566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
1912fa5cd67SKarl Rupp 
192075768bcSBarry Smith   ((PC_Factor *)icc)->info.fill      = 1.0;
193075768bcSBarry Smith   ((PC_Factor *)icc)->info.dtcol     = PETSC_DEFAULT;
194f4db908eSBarry Smith   ((PC_Factor *)icc)->info.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
1959b54502bSHong Zhang 
1969b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
1977b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
198d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
1994cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
200574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2019b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2029b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
20392e08861SBarry Smith   pc->ops->view                = PCView_Factor;
2049b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2059b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2079b54502bSHong Zhang }
208