xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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;
79b54502bSHong Zhang   IS                     perm,cperm;
89b54502bSHong Zhang   PetscErrorCode         ierr;
9f3a39becSBarry Smith   MatInfo                info;
10ea799195SBarry Smith   MatSolverType          stype;
1100e125f8SBarry Smith   MatFactorError         err;
129b54502bSHong Zhang 
139b54502bSHong Zhang   PetscFunctionBegin;
14c6e4fdc6SHong Zhang   pc->failedreason = PC_NOERROR;
15075768bcSBarry Smith   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
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);
21a1f19f5aSHong Zhang     }
22075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
239b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
246bf464f9SBarry Smith     ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
2514798fb4SJed Brown     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
26075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
279b54502bSHong Zhang   }
28075768bcSBarry Smith   ierr                = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
293d1c1ea0SBarry Smith   icc->hdr.actualfill = info.fill_ratio_needed;
30f3a39becSBarry Smith 
31fcfd50ebSBarry Smith   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
32fcfd50ebSBarry Smith   ierr = ISDestroy(&perm);CHKERRQ(ierr);
336baea169SHong Zhang 
3400e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
3500e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
3600e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
376baea169SHong Zhang     PetscFunctionReturn(0);
386baea169SHong Zhang   }
396baea169SHong Zhang 
40075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
4100e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
4200e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
4300e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
446baea169SHong Zhang   }
4500c67f3bSHong Zhang 
463ca39a21SBarry Smith   ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr);
4700c67f3bSHong Zhang   if (!stype) {
48ea799195SBarry Smith     MatSolverType solverpackage;
493ca39a21SBarry Smith     ierr = MatFactorGetSolverType(((PC_Factor*)icc)->fact,&solverpackage);CHKERRQ(ierr);
503ca39a21SBarry Smith     ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr);
5100c67f3bSHong Zhang   }
529b54502bSHong Zhang   PetscFunctionReturn(0);
539b54502bSHong Zhang }
549b54502bSHong Zhang 
55574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
56574deadeSBarry Smith {
57574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
58574deadeSBarry Smith   PetscErrorCode ierr;
59574deadeSBarry Smith 
60574deadeSBarry Smith   PetscFunctionBegin;
616bf464f9SBarry Smith   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
62574deadeSBarry Smith   PetscFunctionReturn(0);
63574deadeSBarry Smith }
64574deadeSBarry Smith 
659b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
669b54502bSHong Zhang {
679b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
689b54502bSHong Zhang   PetscErrorCode ierr;
699b54502bSHong Zhang 
709b54502bSHong Zhang   PetscFunctionBegin;
71574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
72503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
73503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
74c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
759b54502bSHong Zhang   PetscFunctionReturn(0);
769b54502bSHong Zhang }
779b54502bSHong Zhang 
789b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
799b54502bSHong Zhang {
809b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
819b54502bSHong Zhang   PetscErrorCode ierr;
829b54502bSHong Zhang 
839b54502bSHong Zhang   PetscFunctionBegin;
84075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
859b54502bSHong Zhang   PetscFunctionReturn(0);
869b54502bSHong Zhang }
879b54502bSHong Zhang 
889b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
899b54502bSHong Zhang {
909b54502bSHong Zhang   PetscErrorCode ierr;
919b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
929b54502bSHong Zhang 
939b54502bSHong Zhang   PetscFunctionBegin;
94075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
959b54502bSHong Zhang   PetscFunctionReturn(0);
969b54502bSHong Zhang }
979b54502bSHong Zhang 
989b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
999b54502bSHong Zhang {
1009b54502bSHong Zhang   PetscErrorCode ierr;
1019b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1029b54502bSHong Zhang 
1039b54502bSHong Zhang   PetscFunctionBegin;
104075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1059b54502bSHong Zhang   PetscFunctionReturn(0);
1069b54502bSHong Zhang }
1079b54502bSHong Zhang 
1084416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1099b54502bSHong Zhang {
1109b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
111ace3abfcSBarry Smith   PetscBool      flg;
1129b54502bSHong Zhang   PetscErrorCode ierr;
1132fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1149b54502bSHong Zhang 
1159b54502bSHong Zhang   PetscFunctionBegin;
116e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
117e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1189b54502bSHong Zhang 
1198ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1205a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1214c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1224c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
12378fc6b22SHong Zhang   PetscInt       dtmax = 3;
124b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1254c9036c7SBarry Smith   if (flg) {
126b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1274c9036c7SBarry Smith   }
12878fc6b22SHong Zhang   */
1299b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1309b54502bSHong Zhang   PetscFunctionReturn(0);
1319b54502bSHong Zhang }
1329b54502bSHong Zhang 
1339b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1349b54502bSHong Zhang {
1359b54502bSHong Zhang   PetscErrorCode ierr;
1369b54502bSHong Zhang 
1379b54502bSHong Zhang   PetscFunctionBegin;
138914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1399b54502bSHong Zhang   PetscFunctionReturn(0);
1409b54502bSHong Zhang }
1419b54502bSHong Zhang 
1427087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1434c9036c7SBarry Smith 
1449b54502bSHong Zhang /*MC
1459b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1469b54502bSHong Zhang 
1479b54502bSHong Zhang    Options Database Keys:
1482401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1492401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1509b54502bSHong Zhang                       its factorization (overwrites original matrix)
15155ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
152145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1539b54502bSHong Zhang 
1549b54502bSHong Zhang    Level: beginner
1559b54502bSHong Zhang 
1569b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1579b54502bSHong Zhang 
158*95452b02SPatrick Sanan    Notes:
159*95452b02SPatrick Sanan     Only implemented for some matrix formats. Not implemented in parallel.
1609b54502bSHong Zhang 
1619b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1629b54502bSHong Zhang 
1639b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1649b54502bSHong Zhang 
16514d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1669b54502bSHong Zhang           to turn off the shift.
1679b54502bSHong Zhang 
168c582cd25SBarry Smith    References:
16996a0c994SBarry Smith .  1. - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
17096a0c994SBarry Smith       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
17196a0c994SBarry Smith       Science and Engineering, Kluwer.
172c582cd25SBarry Smith 
1739b54502bSHong Zhang 
1749b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
17581f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
176e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1778ff23777SHong Zhang            PCFactorSetLevels()
1789b54502bSHong Zhang 
1799b54502bSHong Zhang M*/
1809b54502bSHong Zhang 
1818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
1829b54502bSHong Zhang {
1839b54502bSHong Zhang   PetscErrorCode ierr;
1849b54502bSHong Zhang   PC_ICC         *icc;
1859b54502bSHong Zhang 
1869b54502bSHong Zhang   PetscFunctionBegin;
187b00a9115SJed Brown   ierr     = PetscNewLog(pc,&icc);CHKERRQ(ierr);
1883d1c1ea0SBarry Smith   pc->data = (void*)icc;
1893d1c1ea0SBarry Smith   ierr     = PCFactorInitialize(pc);CHKERRQ(ierr);
19019fd82e9SBarry Smith   ierr     = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
1912fa5cd67SKarl Rupp 
192879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype     = MAT_FACTOR_ICC;
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;
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;
2039b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
2049b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2059b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2069b54502bSHong Zhang   PetscFunctionReturn(0);
2079b54502bSHong Zhang }
2089b54502bSHong Zhang 
2099b54502bSHong Zhang 
210