xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 6baea16979a9c8ee8132e9efba43c7ba974ff0b3)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h>   /*I "petscpc.h" I*/
39b54502bSHong Zhang 
49b54502bSHong Zhang #undef __FUNCT__
59b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC"
69b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc)
79b54502bSHong Zhang {
89b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
99b54502bSHong Zhang   IS             perm,cperm;
109b54502bSHong Zhang   PetscErrorCode ierr;
11f3a39becSBarry Smith   MatInfo        info;
129b54502bSHong Zhang 
139b54502bSHong Zhang   PetscFunctionBegin;
14075768bcSBarry Smith   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
159b54502bSHong Zhang 
169b54502bSHong Zhang   if (!pc->setupcalled) {
17d09a07f4SBarry Smith     if (!((PC_Factor*)icc)->fact) {
18174acd27SBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
19a1f19f5aSHong Zhang     }
20075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
219b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
226bf464f9SBarry Smith     ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
2314798fb4SJed Brown     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
24075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
259b54502bSHong Zhang   }
26075768bcSBarry Smith   ierr            = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
27f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
28f3a39becSBarry Smith 
29fcfd50ebSBarry Smith   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
30fcfd50ebSBarry Smith   ierr = ISDestroy(&perm);CHKERRQ(ierr);
31*6baea169SHong Zhang 
32*6baea169SHong Zhang   if (((PC_Factor*)icc)->info.errortype) { /* FactorSymbolic() fails */
33*6baea169SHong Zhang     MatFactorInfo factinfo=((PC_Factor*)icc)->info;
34*6baea169SHong Zhang     pc->failedreason = (PCFailedReason)factinfo.errortype;
35*6baea169SHong Zhang     PetscFunctionReturn(0);
36*6baea169SHong Zhang   }
37*6baea169SHong Zhang 
38075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
39*6baea169SHong Zhang   if (((PC_Factor*)icc)->info.errortype) { /* FactorNumeric() fails */
40*6baea169SHong Zhang     MatFactorInfo factinfo=((PC_Factor*)icc)->info;
41*6baea169SHong Zhang     pc->failedreason = (PCFailedReason)factinfo.errortype;
42*6baea169SHong Zhang   }
439b54502bSHong Zhang   PetscFunctionReturn(0);
449b54502bSHong Zhang }
459b54502bSHong Zhang 
469b54502bSHong Zhang #undef __FUNCT__
47574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC"
48574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
49574deadeSBarry Smith {
50574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
51574deadeSBarry Smith   PetscErrorCode ierr;
52574deadeSBarry Smith 
53574deadeSBarry Smith   PetscFunctionBegin;
546bf464f9SBarry Smith   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
55574deadeSBarry Smith   PetscFunctionReturn(0);
56574deadeSBarry Smith }
57574deadeSBarry Smith 
58574deadeSBarry Smith #undef __FUNCT__
599b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
609b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
619b54502bSHong Zhang {
629b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
639b54502bSHong Zhang   PetscErrorCode ierr;
649b54502bSHong Zhang 
659b54502bSHong Zhang   PetscFunctionBegin;
66574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
67503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
68503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
69c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
709b54502bSHong Zhang   PetscFunctionReturn(0);
719b54502bSHong Zhang }
729b54502bSHong Zhang 
739b54502bSHong Zhang #undef __FUNCT__
749b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
759b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
769b54502bSHong Zhang {
779b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
789b54502bSHong Zhang   PetscErrorCode ierr;
799b54502bSHong Zhang 
809b54502bSHong Zhang   PetscFunctionBegin;
81075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
829b54502bSHong Zhang   PetscFunctionReturn(0);
839b54502bSHong Zhang }
849b54502bSHong Zhang 
859b54502bSHong Zhang #undef __FUNCT__
869b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
879b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
889b54502bSHong Zhang {
899b54502bSHong Zhang   PetscErrorCode ierr;
909b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
919b54502bSHong Zhang 
929b54502bSHong Zhang   PetscFunctionBegin;
93075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
949b54502bSHong Zhang   PetscFunctionReturn(0);
959b54502bSHong Zhang }
969b54502bSHong Zhang 
979b54502bSHong Zhang #undef __FUNCT__
989b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
999b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1009b54502bSHong Zhang {
1019b54502bSHong Zhang   PetscErrorCode ierr;
1029b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1039b54502bSHong Zhang 
1049b54502bSHong Zhang   PetscFunctionBegin;
105075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1069b54502bSHong Zhang   PetscFunctionReturn(0);
1079b54502bSHong Zhang }
1089b54502bSHong Zhang 
1099b54502bSHong Zhang #undef __FUNCT__
1109b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1114416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1129b54502bSHong Zhang {
1139b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
114ace3abfcSBarry Smith   PetscBool      flg;
1159b54502bSHong Zhang   PetscErrorCode ierr;
1162fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1179b54502bSHong Zhang 
1189b54502bSHong Zhang   PetscFunctionBegin;
119e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
120e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1219b54502bSHong Zhang 
1228ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1235a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1244c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1254c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
12678fc6b22SHong Zhang   PetscInt       dtmax = 3;
127b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1284c9036c7SBarry Smith   if (flg) {
129b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1304c9036c7SBarry Smith   }
13178fc6b22SHong Zhang   */
1329b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1339b54502bSHong Zhang   PetscFunctionReturn(0);
1349b54502bSHong Zhang }
1359b54502bSHong Zhang 
1369b54502bSHong Zhang #undef __FUNCT__
1379b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1389b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1399b54502bSHong Zhang {
1409b54502bSHong Zhang   PetscErrorCode ierr;
1419b54502bSHong Zhang 
1429b54502bSHong Zhang   PetscFunctionBegin;
143914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1449b54502bSHong Zhang   PetscFunctionReturn(0);
1459b54502bSHong Zhang }
1469b54502bSHong Zhang 
1477087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1484c9036c7SBarry Smith 
149c60c7ad4SBarry Smith #undef __FUNCT__
150c60c7ad4SBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_ICC"
151c60c7ad4SBarry Smith PetscErrorCode  PCFactorGetUseInPlace_ICC(PC pc,PetscBool *flg)
152c60c7ad4SBarry Smith {
153c60c7ad4SBarry Smith   PetscFunctionBegin;
154c60c7ad4SBarry Smith   *flg = PETSC_FALSE;
155c60c7ad4SBarry Smith   PetscFunctionReturn(0);
156c60c7ad4SBarry Smith }
157c60c7ad4SBarry Smith 
1589b54502bSHong Zhang /*MC
1599b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1609b54502bSHong Zhang 
1619b54502bSHong Zhang    Options Database Keys:
1622401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1632401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1649b54502bSHong Zhang                       its factorization (overwrites original matrix)
16555ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
166145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1679b54502bSHong Zhang 
1689b54502bSHong Zhang    Level: beginner
1699b54502bSHong Zhang 
1709b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1719b54502bSHong Zhang 
172d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel.
1739b54502bSHong Zhang 
1749b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1759b54502bSHong Zhang 
1769b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1779b54502bSHong Zhang 
17814d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1799b54502bSHong Zhang           to turn off the shift.
1809b54502bSHong Zhang 
181c582cd25SBarry Smith    References:
182c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
183c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
184c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
185c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
186c582cd25SBarry Smith 
1879b54502bSHong Zhang 
1889b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
18981f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
190e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1918ff23777SHong Zhang            PCFactorSetLevels()
1929b54502bSHong Zhang 
1939b54502bSHong Zhang M*/
1949b54502bSHong Zhang 
1959b54502bSHong Zhang #undef __FUNCT__
1969b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
1978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
1989b54502bSHong Zhang {
1999b54502bSHong Zhang   PetscErrorCode ierr;
2009b54502bSHong Zhang   PC_ICC         *icc;
2019b54502bSHong Zhang 
2029b54502bSHong Zhang   PetscFunctionBegin;
203b00a9115SJed Brown   ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr);
2049b54502bSHong Zhang 
205075768bcSBarry Smith   ((PC_Factor*)icc)->fact = 0;
2062fa5cd67SKarl Rupp 
20719fd82e9SBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
2082692d6eeSBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
209075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
2102fa5cd67SKarl Rupp 
211879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype  = MAT_FACTOR_ICC;
21275567043SBarry Smith   ((PC_Factor*)icc)->info.levels = 0.;
213075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill   = 1.0;
2149b54502bSHong Zhang   icc->implctx                   = 0;
2159b54502bSHong Zhang 
216075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
217f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2180ed735ceSHong Zhang   ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
2190ed735ceSHong Zhang   ((PC_Factor*)icc)->info.zeropivot   = 100.0*PETSC_MACHINE_EPSILON;
2209b54502bSHong Zhang 
221915743fcSHong Zhang   pc->data                     = (void*)icc;
2229b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
223d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2249b54502bSHong Zhang   pc->ops->setup               = PCSetup_ICC;
225574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2269b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2279b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
2289b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
22985317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
2309b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2319b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2329b54502bSHong Zhang 
233bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
2392591b318SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
240bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
244c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ICC);CHKERRQ(ierr);
2459b54502bSHong Zhang   PetscFunctionReturn(0);
2469b54502bSHong Zhang }
2479b54502bSHong Zhang 
2489b54502bSHong Zhang 
249