xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 606c028001d6169ccbb8dc4a5aa61aa3bda31a09)
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;
89b54502bSHong Zhang   PetscErrorCode         ierr;
9f3a39becSBarry Smith   MatInfo                info;
10ea799195SBarry Smith   MatSolverType          stype;
1100e125f8SBarry Smith   MatFactorError         err;
124ac6704cSBarry Smith   PetscBool              canuseordering;
139b54502bSHong Zhang 
149b54502bSHong Zhang   PetscFunctionBegin;
15c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
169b54502bSHong Zhang 
1784d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
189b54502bSHong Zhang   if (!pc->setupcalled) {
19d09a07f4SBarry Smith     if (!((PC_Factor*)icc)->fact) {
20174acd27SBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
214ac6704cSBarry Smith     }
22f73b0415SBarry Smith     ierr = MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering);CHKERRQ(ierr);
23f73b0415SBarry Smith     if (canuseordering) {
244ac6704cSBarry Smith       ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr);
252c7c0729SBarry Smith       ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
262c7c0729SBarry Smith     }
27075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
289b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
29f73b0415SBarry Smith     PetscBool canuseordering;
306bf464f9SBarry Smith     ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
3114798fb4SJed Brown     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
32f73b0415SBarry Smith     ierr = MatFactorGetCanUseOrdering(((PC_Factor*)icc)->fact,&canuseordering);CHKERRQ(ierr);
33f73b0415SBarry Smith     if (canuseordering) {
344ac6704cSBarry Smith       ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr);
352c7c0729SBarry Smith       ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
362c7c0729SBarry Smith     }
37075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
389b54502bSHong Zhang   }
39075768bcSBarry Smith   ierr                = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
403d1c1ea0SBarry Smith   icc->hdr.actualfill = info.fill_ratio_needed;
41f3a39becSBarry Smith 
42fcfd50ebSBarry Smith   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
43fcfd50ebSBarry Smith   ierr = ISDestroy(&perm);CHKERRQ(ierr);
446baea169SHong Zhang 
4500e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
4600e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
4700e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
486baea169SHong Zhang     PetscFunctionReturn(0);
496baea169SHong Zhang   }
506baea169SHong Zhang 
51075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
5200e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
5300e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
5400e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
556baea169SHong Zhang   }
5600c67f3bSHong Zhang 
573ca39a21SBarry Smith   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
5800c67f3bSHong Zhang   if (!stype) {
59ea799195SBarry Smith     MatSolverType solverpackage;
603ca39a21SBarry Smith     ierr = MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage);CHKERRQ(ierr);
613ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
6200c67f3bSHong Zhang   }
639b54502bSHong Zhang   PetscFunctionReturn(0);
649b54502bSHong Zhang }
659b54502bSHong Zhang 
66574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
67574deadeSBarry Smith {
68574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
69574deadeSBarry Smith   PetscErrorCode ierr;
70574deadeSBarry Smith 
71574deadeSBarry Smith   PetscFunctionBegin;
726bf464f9SBarry Smith   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
73574deadeSBarry Smith   PetscFunctionReturn(0);
74574deadeSBarry Smith }
75574deadeSBarry Smith 
769b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
779b54502bSHong Zhang {
789b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
799b54502bSHong Zhang   PetscErrorCode ierr;
809b54502bSHong Zhang 
819b54502bSHong Zhang   PetscFunctionBegin;
82574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
83503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
84503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
85c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
869b54502bSHong Zhang   PetscFunctionReturn(0);
879b54502bSHong Zhang }
889b54502bSHong Zhang 
899b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
909b54502bSHong Zhang {
919b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
929b54502bSHong Zhang   PetscErrorCode ierr;
939b54502bSHong Zhang 
949b54502bSHong Zhang   PetscFunctionBegin;
95075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
969b54502bSHong Zhang   PetscFunctionReturn(0);
979b54502bSHong Zhang }
989b54502bSHong Zhang 
997b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ICC(PC pc,Mat X,Mat Y)
1007b6e2003SPierre Jolivet {
1017b6e2003SPierre Jolivet   PC_ICC         *icc = (PC_ICC*)pc->data;
1027b6e2003SPierre Jolivet   PetscErrorCode ierr;
1037b6e2003SPierre Jolivet 
1047b6e2003SPierre Jolivet   PetscFunctionBegin;
1057b6e2003SPierre Jolivet   ierr = MatMatSolve(((PC_Factor*)icc)->fact,X,Y);CHKERRQ(ierr);
1067b6e2003SPierre Jolivet   PetscFunctionReturn(0);
1077b6e2003SPierre Jolivet }
1087b6e2003SPierre Jolivet 
1099b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
1109b54502bSHong Zhang {
1119b54502bSHong Zhang   PetscErrorCode ierr;
1129b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1139b54502bSHong Zhang 
1149b54502bSHong Zhang   PetscFunctionBegin;
115075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1169b54502bSHong Zhang   PetscFunctionReturn(0);
1179b54502bSHong Zhang }
1189b54502bSHong Zhang 
1199b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1209b54502bSHong Zhang {
1219b54502bSHong Zhang   PetscErrorCode ierr;
1229b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1239b54502bSHong Zhang 
1249b54502bSHong Zhang   PetscFunctionBegin;
125075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1269b54502bSHong Zhang   PetscFunctionReturn(0);
1279b54502bSHong Zhang }
1289b54502bSHong Zhang 
1294416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1309b54502bSHong Zhang {
1319b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
132ace3abfcSBarry Smith   PetscBool      flg;
1339b54502bSHong Zhang   PetscErrorCode ierr;
1342fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1359b54502bSHong Zhang 
1369b54502bSHong Zhang   PetscFunctionBegin;
137e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
138e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1399b54502bSHong Zhang 
1408ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1415a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1424c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1434c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
14478fc6b22SHong Zhang   PetscInt       dtmax = 3;
145b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1464c9036c7SBarry Smith   if (flg) {
147b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1484c9036c7SBarry Smith   }
14978fc6b22SHong Zhang   */
1509b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1519b54502bSHong Zhang   PetscFunctionReturn(0);
1529b54502bSHong Zhang }
1539b54502bSHong Zhang 
1547087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1554c9036c7SBarry Smith 
1569b54502bSHong Zhang /*MC
1579b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1589b54502bSHong Zhang 
1599b54502bSHong Zhang    Options Database Keys:
1602401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1612401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1629b54502bSHong Zhang                       its factorization (overwrites original matrix)
16355ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
164145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1659b54502bSHong Zhang 
1669b54502bSHong Zhang    Level: beginner
1679b54502bSHong Zhang 
16895452b02SPatrick Sanan    Notes:
16995452b02SPatrick Sanan     Only implemented for some matrix formats. Not implemented in parallel.
1709b54502bSHong Zhang 
1719b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1729b54502bSHong Zhang 
1739b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1749b54502bSHong Zhang 
17514d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1769b54502bSHong Zhang           to turn off the shift.
1779b54502bSHong Zhang 
178c582cd25SBarry Smith    References:
179*606c0280SSatish Balay .  * - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
18096a0c994SBarry Smith       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
18196a0c994SBarry Smith       Science and Engineering, Kluwer.
182c582cd25SBarry Smith 
1839b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
18481f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
185e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1868ff23777SHong Zhang            PCFactorSetLevels()
1879b54502bSHong Zhang 
1889b54502bSHong Zhang M*/
1899b54502bSHong Zhang 
1908cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
1919b54502bSHong Zhang {
1929b54502bSHong Zhang   PetscErrorCode ierr;
1939b54502bSHong Zhang   PC_ICC         *icc;
1949b54502bSHong Zhang 
1959b54502bSHong Zhang   PetscFunctionBegin;
196b00a9115SJed Brown   ierr     = PetscNewLog(pc,&icc);CHKERRQ(ierr);
1973d1c1ea0SBarry Smith   pc->data = (void*)icc;
1984ac6704cSBarry Smith   ierr     = PCFactorInitialize(pc, MAT_FACTOR_ICC);CHKERRQ(ierr);
1992fa5cd67SKarl Rupp 
200075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill      = 1.0;
201075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol     = PETSC_DEFAULT;
202f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2039b54502bSHong Zhang 
2049b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
2057b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_ICC;
206d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2074cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
208574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2099b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2109b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
21192e08861SBarry Smith   pc->ops->view                = PCView_Factor;
2129b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2139b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2149b54502bSHong Zhang   PetscFunctionReturn(0);
2159b54502bSHong Zhang }
216