xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 4cccfbdd2a20c44cba925148732fe999d2dfbe44)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h>   /*I "petscpc.h" I*/
39b54502bSHong Zhang 
49b54502bSHong Zhang #undef __FUNCT__
5*4cccfbddSHong Zhang #define __FUNCT__ "PCSetUp_ICC"
6*4cccfbddSHong 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;
12*4cccfbddSHong Zhang   Mat            F;
139b54502bSHong Zhang 
149b54502bSHong Zhang   PetscFunctionBegin;
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);
29f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
30f3a39becSBarry Smith 
31fcfd50ebSBarry Smith   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
32fcfd50ebSBarry Smith   ierr = ISDestroy(&perm);CHKERRQ(ierr);
336baea169SHong Zhang 
34*4cccfbddSHong Zhang   F = ((PC_Factor*)icc)->fact;
35*4cccfbddSHong Zhang   if (F->errortype) { /* FactorSymbolic() fails */
36*4cccfbddSHong Zhang     pc->failedreason = (PCFailedReason)F->errortype;
376baea169SHong Zhang     PetscFunctionReturn(0);
386baea169SHong Zhang   }
396baea169SHong Zhang 
40075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
41*4cccfbddSHong Zhang   if (F->errortype) { /* FactorNumeric() fails */
42*4cccfbddSHong Zhang     pc->failedreason = (PCFailedReason)F->errortype;
436baea169SHong Zhang   }
449b54502bSHong Zhang   PetscFunctionReturn(0);
459b54502bSHong Zhang }
469b54502bSHong Zhang 
479b54502bSHong Zhang #undef __FUNCT__
48574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC"
49574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
50574deadeSBarry Smith {
51574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
52574deadeSBarry Smith   PetscErrorCode ierr;
53574deadeSBarry Smith 
54574deadeSBarry Smith   PetscFunctionBegin;
556bf464f9SBarry Smith   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
56574deadeSBarry Smith   PetscFunctionReturn(0);
57574deadeSBarry Smith }
58574deadeSBarry Smith 
59574deadeSBarry Smith #undef __FUNCT__
609b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
619b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
629b54502bSHong Zhang {
639b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
649b54502bSHong Zhang   PetscErrorCode ierr;
659b54502bSHong Zhang 
669b54502bSHong Zhang   PetscFunctionBegin;
67574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
68503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
69503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
70c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
719b54502bSHong Zhang   PetscFunctionReturn(0);
729b54502bSHong Zhang }
739b54502bSHong Zhang 
749b54502bSHong Zhang #undef __FUNCT__
759b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
769b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
779b54502bSHong Zhang {
789b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
799b54502bSHong Zhang   PetscErrorCode ierr;
809b54502bSHong Zhang 
819b54502bSHong Zhang   PetscFunctionBegin;
82075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
839b54502bSHong Zhang   PetscFunctionReturn(0);
849b54502bSHong Zhang }
859b54502bSHong Zhang 
869b54502bSHong Zhang #undef __FUNCT__
879b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
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 #undef __FUNCT__
999b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
1009b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1019b54502bSHong Zhang {
1029b54502bSHong Zhang   PetscErrorCode ierr;
1039b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1049b54502bSHong Zhang 
1059b54502bSHong Zhang   PetscFunctionBegin;
106075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1079b54502bSHong Zhang   PetscFunctionReturn(0);
1089b54502bSHong Zhang }
1099b54502bSHong Zhang 
1109b54502bSHong Zhang #undef __FUNCT__
1119b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1124416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1139b54502bSHong Zhang {
1149b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
115ace3abfcSBarry Smith   PetscBool      flg;
1169b54502bSHong Zhang   PetscErrorCode ierr;
1172fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1189b54502bSHong Zhang 
1199b54502bSHong Zhang   PetscFunctionBegin;
120e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
121e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1229b54502bSHong Zhang 
1238ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1245a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1254c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1264c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
12778fc6b22SHong Zhang   PetscInt       dtmax = 3;
128b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1294c9036c7SBarry Smith   if (flg) {
130b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1314c9036c7SBarry Smith   }
13278fc6b22SHong Zhang   */
1339b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1349b54502bSHong Zhang   PetscFunctionReturn(0);
1359b54502bSHong Zhang }
1369b54502bSHong Zhang 
1379b54502bSHong Zhang #undef __FUNCT__
1389b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1399b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1409b54502bSHong Zhang {
1419b54502bSHong Zhang   PetscErrorCode ierr;
1429b54502bSHong Zhang 
1439b54502bSHong Zhang   PetscFunctionBegin;
144914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1459b54502bSHong Zhang   PetscFunctionReturn(0);
1469b54502bSHong Zhang }
1479b54502bSHong Zhang 
1487087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1494c9036c7SBarry Smith 
150c60c7ad4SBarry Smith #undef __FUNCT__
151c60c7ad4SBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_ICC"
152c60c7ad4SBarry Smith PetscErrorCode  PCFactorGetUseInPlace_ICC(PC pc,PetscBool *flg)
153c60c7ad4SBarry Smith {
154c60c7ad4SBarry Smith   PetscFunctionBegin;
155c60c7ad4SBarry Smith   *flg = PETSC_FALSE;
156c60c7ad4SBarry Smith   PetscFunctionReturn(0);
157c60c7ad4SBarry Smith }
158c60c7ad4SBarry Smith 
1599b54502bSHong Zhang /*MC
1609b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1619b54502bSHong Zhang 
1629b54502bSHong Zhang    Options Database Keys:
1632401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1642401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1659b54502bSHong Zhang                       its factorization (overwrites original matrix)
16655ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
167145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1689b54502bSHong Zhang 
1699b54502bSHong Zhang    Level: beginner
1709b54502bSHong Zhang 
1719b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1729b54502bSHong Zhang 
173d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel.
1749b54502bSHong Zhang 
1759b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1769b54502bSHong Zhang 
1779b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1789b54502bSHong Zhang 
17914d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1809b54502bSHong Zhang           to turn off the shift.
1819b54502bSHong Zhang 
182c582cd25SBarry Smith    References:
183c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
184c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
185c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
186c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
187c582cd25SBarry Smith 
1889b54502bSHong Zhang 
1899b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
19081f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
191e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1928ff23777SHong Zhang            PCFactorSetLevels()
1939b54502bSHong Zhang 
1949b54502bSHong Zhang M*/
1959b54502bSHong Zhang 
1969b54502bSHong Zhang #undef __FUNCT__
1979b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
1988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
1999b54502bSHong Zhang {
2009b54502bSHong Zhang   PetscErrorCode ierr;
2019b54502bSHong Zhang   PC_ICC         *icc;
2029b54502bSHong Zhang 
2039b54502bSHong Zhang   PetscFunctionBegin;
204b00a9115SJed Brown   ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr);
2059b54502bSHong Zhang 
206075768bcSBarry Smith   ((PC_Factor*)icc)->fact = 0;
2072fa5cd67SKarl Rupp 
20819fd82e9SBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
2092692d6eeSBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
210075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
2112fa5cd67SKarl Rupp 
212879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype  = MAT_FACTOR_ICC;
21375567043SBarry Smith   ((PC_Factor*)icc)->info.levels = 0.;
214075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill   = 1.0;
2159b54502bSHong Zhang   icc->implctx                   = 0;
2169b54502bSHong Zhang 
217075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
218f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2190ed735ceSHong Zhang   ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
2200ed735ceSHong Zhang   ((PC_Factor*)icc)->info.zeropivot   = 100.0*PETSC_MACHINE_EPSILON;
2219b54502bSHong Zhang 
222915743fcSHong Zhang   pc->data                     = (void*)icc;
2239b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
224d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
225*4cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
226574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2279b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2289b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
2299b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
23085317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
2319b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2329b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2339b54502bSHong Zhang 
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
2402591b318SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
244bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
245c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ICC);CHKERRQ(ierr);
2469b54502bSHong Zhang   PetscFunctionReturn(0);
2479b54502bSHong Zhang }
2489b54502bSHong Zhang 
2499b54502bSHong Zhang 
250