xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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) {
2248a46eb9SPierre Jolivet     if (!((PC_Factor *)icc)->fact) PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)icc)->solvertype, MAT_FACTOR_ICC, &((PC_Factor *)icc)->fact));
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) {
30f73b0415SBarry Smith     PetscBool canuseordering;
319566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
329566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(pc->pmat, ((PC_Factor *)icc)->solvertype, MAT_FACTOR_ICC, &((PC_Factor *)icc)->fact));
339566063dSJacob Faibussowitsch     PetscCall(MatFactorGetCanUseOrdering(((PC_Factor *)icc)->fact, &canuseordering));
34f73b0415SBarry Smith     if (canuseordering) {
359566063dSJacob Faibussowitsch       PetscCall(PCFactorSetDefaultOrdering_Factor(pc));
369566063dSJacob Faibussowitsch       PetscCall(MatGetOrdering(pc->pmat, ((PC_Factor *)icc)->ordering, &perm, &cperm));
372c7c0729SBarry Smith     }
389566063dSJacob Faibussowitsch     PetscCall(MatICCFactorSymbolic(((PC_Factor *)icc)->fact, pc->pmat, perm, &((PC_Factor *)icc)->info));
399b54502bSHong Zhang   }
409566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(((PC_Factor *)icc)->fact, MAT_LOCAL, &info));
413d1c1ea0SBarry Smith   icc->hdr.actualfill = info.fill_ratio_needed;
42f3a39becSBarry Smith 
439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cperm));
449566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
456baea169SHong Zhang 
469566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor *)icc)->fact, &err));
4700e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
4800e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
49*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
506baea169SHong Zhang   }
516baea169SHong Zhang 
529566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(((PC_Factor *)icc)->fact, pc->pmat, &((PC_Factor *)icc)->info));
539566063dSJacob Faibussowitsch   PetscCall(MatFactorGetError(((PC_Factor *)icc)->fact, &err));
5400e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
5500e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
566baea169SHong Zhang   }
5700c67f3bSHong Zhang 
589566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatSolverType(pc, &stype));
5900c67f3bSHong Zhang   if (!stype) {
60ea799195SBarry Smith     MatSolverType solverpackage;
619566063dSJacob Faibussowitsch     PetscCall(MatFactorGetSolverType(((PC_Factor *)icc)->fact, &solverpackage));
629566063dSJacob Faibussowitsch     PetscCall(PCFactorSetMatSolverType(pc, solverpackage));
6300c67f3bSHong Zhang   }
64*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
659b54502bSHong Zhang }
669b54502bSHong Zhang 
67d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ICC(PC pc)
68d71ae5a4SJacob Faibussowitsch {
69574deadeSBarry Smith   PC_ICC *icc = (PC_ICC *)pc->data;
70574deadeSBarry Smith 
71574deadeSBarry Smith   PetscFunctionBegin;
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&((PC_Factor *)icc)->fact));
73*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74574deadeSBarry Smith }
75574deadeSBarry Smith 
76d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ICC(PC pc)
77d71ae5a4SJacob Faibussowitsch {
789b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
799b54502bSHong Zhang 
809b54502bSHong Zhang   PetscFunctionBegin;
819566063dSJacob Faibussowitsch   PetscCall(PCReset_ICC(pc));
829566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)icc)->ordering));
839566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PC_Factor *)icc)->solvertype));
842e956fe4SStefano Zampini   PetscCall(PCFactorClearComposedFunctions(pc));
859566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
86*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
879b54502bSHong Zhang }
889b54502bSHong Zhang 
89d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ICC(PC pc, Vec x, Vec y)
90d71ae5a4SJacob Faibussowitsch {
919b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
929b54502bSHong Zhang 
939b54502bSHong Zhang   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(MatSolve(((PC_Factor *)icc)->fact, x, y));
95*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
969b54502bSHong Zhang }
979b54502bSHong Zhang 
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ICC(PC pc, Mat X, Mat Y)
99d71ae5a4SJacob Faibussowitsch {
1007b6e2003SPierre Jolivet   PC_ICC *icc = (PC_ICC *)pc->data;
1017b6e2003SPierre Jolivet 
1027b6e2003SPierre Jolivet   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(MatMatSolve(((PC_Factor *)icc)->fact, X, Y));
104*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1057b6e2003SPierre Jolivet }
1067b6e2003SPierre Jolivet 
107d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc, Vec x, Vec y)
108d71ae5a4SJacob Faibussowitsch {
1099b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1109b54502bSHong Zhang 
1119b54502bSHong Zhang   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(MatForwardSolve(((PC_Factor *)icc)->fact, x, y));
113*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1149b54502bSHong Zhang }
1159b54502bSHong Zhang 
116d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplySymmetricRight_ICC(PC pc, Vec x, Vec y)
117d71ae5a4SJacob Faibussowitsch {
1189b54502bSHong Zhang   PC_ICC *icc = (PC_ICC *)pc->data;
1199b54502bSHong Zhang 
1209b54502bSHong Zhang   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(MatBackwardSolve(((PC_Factor *)icc)->fact, x, y));
122*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1239b54502bSHong Zhang }
1249b54502bSHong Zhang 
125d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ICC(PC pc, PetscOptionItems *PetscOptionsObject)
126d71ae5a4SJacob Faibussowitsch {
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");
133dbbe0bcdSBarry Smith   PetscCall(PCSetFromOptions_Factor(pc, PetscOptionsObject));
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();
146*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
166f1580f4eSBarry Smith    For `MATSEQBAIJ` matrices this implements a point block ICC.
167f1580f4eSBarry Smith 
168f1580f4eSBarry Smith    By default, the Manteuffel is applied (for matrices with block size 1). Call `PCFactorSetShiftType`(pc,`MAT_SHIFT_POSITIVE_DEFINITE`);
169f1580f4eSBarry Smith    to turn off the shift.
1709b54502bSHong Zhang 
1719b54502bSHong Zhang    The Manteuffel shift is only implemented for matrices with block size 1
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 
178f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSOR`, `MatOrderingType`, `PCILU`, `PCLU`, `PCCHOLESKY`,
179db781477SPatrick Sanan           `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()`,
180db781477SPatrick Sanan           `PCFactorSetFill()`, `PCFactorSetMatOrderingType()`, `PCFactorSetReuseOrdering()`,
181db781477SPatrick Sanan           `PCFactorSetLevels()`
1829b54502bSHong Zhang M*/
1839b54502bSHong Zhang 
184d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
185d71ae5a4SJacob Faibussowitsch {
1869b54502bSHong Zhang   PC_ICC *icc;
1879b54502bSHong Zhang 
1889b54502bSHong Zhang   PetscFunctionBegin;
1894dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&icc));
1903d1c1ea0SBarry Smith   pc->data = (void *)icc;
1919566063dSJacob Faibussowitsch   PetscCall(PCFactorInitialize(pc, MAT_FACTOR_ICC));
1922fa5cd67SKarl Rupp 
193075768bcSBarry Smith   ((PC_Factor *)icc)->info.fill      = 1.0;
194075768bcSBarry Smith   ((PC_Factor *)icc)->info.dtcol     = PETSC_DEFAULT;
195f4db908eSBarry Smith   ((PC_Factor *)icc)->info.shifttype = (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE;
1969b54502bSHong Zhang 
1979b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
1987b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
199d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2004cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
201574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2029b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2039b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
20492e08861SBarry Smith   pc->ops->view                = PCView_Factor;
2059b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2069b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
207*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2089b54502bSHong Zhang }
209