xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 574deadea3b5f81f544fd27eaca8164a96d1c060)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
37c4f633dSBarry Smith #include "../src/ksp/pc/impls/factor/icc/icc.h"   /*I "petscpc.h" I*/
49b54502bSHong Zhang 
59b54502bSHong Zhang #undef __FUNCT__
69b54502bSHong Zhang #define __FUNCT__ "PCSetup_ICC"
79b54502bSHong Zhang static PetscErrorCode PCSetup_ICC(PC pc)
89b54502bSHong Zhang {
99b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
109b54502bSHong Zhang   IS             perm,cperm;
119b54502bSHong Zhang   PetscErrorCode ierr;
12f3a39becSBarry Smith   MatInfo        info;
139b54502bSHong Zhang 
149b54502bSHong Zhang   PetscFunctionBegin;
15075768bcSBarry Smith   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
169b54502bSHong Zhang 
179b54502bSHong Zhang   if (!pc->setupcalled) {
18d09a07f4SBarry Smith     if (!((PC_Factor*)icc)->fact){
19174acd27SBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr);
20a1f19f5aSHong Zhang     }
21075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
229b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
23075768bcSBarry Smith     ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);
242692d6eeSBarry Smith     ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
25075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
269b54502bSHong Zhang   }
27075768bcSBarry Smith   ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
28f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
29f3a39becSBarry Smith 
309b54502bSHong Zhang   ierr = ISDestroy(cperm);CHKERRQ(ierr);
319b54502bSHong Zhang   ierr = ISDestroy(perm);CHKERRQ(ierr);
32075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
339b54502bSHong Zhang   PetscFunctionReturn(0);
349b54502bSHong Zhang }
359b54502bSHong Zhang 
369b54502bSHong Zhang #undef __FUNCT__
37*574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC"
38*574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
39*574deadeSBarry Smith {
40*574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
41*574deadeSBarry Smith   PetscErrorCode ierr;
42*574deadeSBarry Smith 
43*574deadeSBarry Smith   PetscFunctionBegin;
44*574deadeSBarry Smith   if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);}
45*574deadeSBarry Smith   PetscFunctionReturn(0);
46*574deadeSBarry Smith }
47*574deadeSBarry Smith 
48*574deadeSBarry Smith #undef __FUNCT__
499b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
509b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
519b54502bSHong Zhang {
529b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
539b54502bSHong Zhang   PetscErrorCode ierr;
549b54502bSHong Zhang 
559b54502bSHong Zhang   PetscFunctionBegin;
56*574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
57503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
58503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
59c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
609b54502bSHong Zhang   PetscFunctionReturn(0);
619b54502bSHong Zhang }
629b54502bSHong Zhang 
639b54502bSHong Zhang #undef __FUNCT__
649b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
659b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
669b54502bSHong Zhang {
679b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
689b54502bSHong Zhang   PetscErrorCode ierr;
699b54502bSHong Zhang 
709b54502bSHong Zhang   PetscFunctionBegin;
71075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
729b54502bSHong Zhang   PetscFunctionReturn(0);
739b54502bSHong Zhang }
749b54502bSHong Zhang 
759b54502bSHong Zhang #undef __FUNCT__
769b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
779b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
789b54502bSHong Zhang {
799b54502bSHong Zhang   PetscErrorCode ierr;
809b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
819b54502bSHong Zhang 
829b54502bSHong Zhang   PetscFunctionBegin;
83075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
849b54502bSHong Zhang   PetscFunctionReturn(0);
859b54502bSHong Zhang }
869b54502bSHong Zhang 
879b54502bSHong Zhang #undef __FUNCT__
889b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
899b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
909b54502bSHong Zhang {
919b54502bSHong Zhang   PetscErrorCode ierr;
929b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
939b54502bSHong Zhang 
949b54502bSHong Zhang   PetscFunctionBegin;
95075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
969b54502bSHong Zhang   PetscFunctionReturn(0);
979b54502bSHong Zhang }
989b54502bSHong Zhang 
999b54502bSHong Zhang #undef __FUNCT__
1009b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1019b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc)
1029b54502bSHong Zhang {
1039b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
104ace3abfcSBarry Smith   PetscBool      flg;
1055a586d82SBarry Smith   /* PetscReal      dt[3];*/
1069b54502bSHong Zhang   PetscErrorCode ierr;
1079b54502bSHong Zhang 
1089b54502bSHong Zhang   PetscFunctionBegin;
1099b54502bSHong Zhang   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
1108ff23777SHong Zhang     ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr);
1119b54502bSHong Zhang 
1128ff23777SHong Zhang     ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1135a586d82SBarry Smith     /*dt[0] = ((PC_Factor*)icc)->info.dt;
1144c9036c7SBarry Smith     dt[1] = ((PC_Factor*)icc)->info.dtcol;
1154c9036c7SBarry Smith     dt[2] = ((PC_Factor*)icc)->info.dtcount;
11678fc6b22SHong Zhang     PetscInt       dtmax = 3;
117b7c853c4SBarry Smith     ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1184c9036c7SBarry Smith     if (flg) {
119b7c853c4SBarry Smith       ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1204c9036c7SBarry Smith     }
12178fc6b22SHong Zhang     */
1229b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1239b54502bSHong Zhang   PetscFunctionReturn(0);
1249b54502bSHong Zhang }
1259b54502bSHong Zhang 
1269b54502bSHong Zhang #undef __FUNCT__
1279b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1289b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1299b54502bSHong Zhang {
1309b54502bSHong Zhang   PetscErrorCode ierr;
1319b54502bSHong Zhang 
1329b54502bSHong Zhang   PetscFunctionBegin;
133914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1349b54502bSHong Zhang   PetscFunctionReturn(0);
1359b54502bSHong Zhang }
1369b54502bSHong Zhang 
137295a483bSBarry Smith EXTERN_C_BEGIN
1387087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
139295a483bSBarry Smith EXTERN_C_END
1404c9036c7SBarry Smith 
1419b54502bSHong Zhang /*MC
1429b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1439b54502bSHong Zhang 
1449b54502bSHong Zhang    Options Database Keys:
1452401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1462401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1479b54502bSHong Zhang                       its factorization (overwrites original matrix)
14855ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
149145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1509b54502bSHong Zhang 
1519b54502bSHong Zhang    Level: beginner
1529b54502bSHong Zhang 
1539b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1549b54502bSHong Zhang 
155d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel.
1569b54502bSHong Zhang 
1579b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1589b54502bSHong Zhang 
1599b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1609b54502bSHong Zhang 
16114d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1629b54502bSHong Zhang           to turn off the shift.
1639b54502bSHong Zhang 
164c582cd25SBarry Smith    References:
165c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
166c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
167c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
168c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
169c582cd25SBarry Smith 
1709b54502bSHong Zhang 
1719b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
17281f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
173e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1748ff23777SHong Zhang            PCFactorSetLevels()
1759b54502bSHong Zhang 
1769b54502bSHong Zhang M*/
1779b54502bSHong Zhang 
1789b54502bSHong Zhang EXTERN_C_BEGIN
1799b54502bSHong Zhang #undef __FUNCT__
1809b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
1817087cfbeSBarry Smith PetscErrorCode  PCCreate_ICC(PC pc)
1829b54502bSHong Zhang {
1839b54502bSHong Zhang   PetscErrorCode ierr;
1849b54502bSHong Zhang   PC_ICC         *icc;
1859b54502bSHong Zhang 
1869b54502bSHong Zhang   PetscFunctionBegin;
18738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr);
1889b54502bSHong Zhang 
189075768bcSBarry Smith   ((PC_Factor*)icc)->fact	          = 0;
1902692d6eeSBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
1912692d6eeSBarry Smith   ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
192075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
193879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype         = MAT_FACTOR_ICC;
19475567043SBarry Smith   ((PC_Factor*)icc)->info.levels	= 0.;
195075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill          = 1.0;
1969b54502bSHong Zhang   icc->implctx            = 0;
1979b54502bSHong Zhang 
198075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
199f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
200daa17b54SHong Zhang   ((PC_Factor*)icc)->info.shiftamount = 1.e-12;
201075768bcSBarry Smith   ((PC_Factor*)icc)->info.zeropivot   = 1.e-12;
2029b54502bSHong Zhang 
203915743fcSHong Zhang   pc->data	               = (void*)icc;
2049b54502bSHong Zhang   pc->ops->apply	       = PCApply_ICC;
2059b54502bSHong Zhang   pc->ops->setup               = PCSetup_ICC;
206*574deadeSBarry Smith   pc->ops->reset  	       = PCReset_ICC;
2079b54502bSHong Zhang   pc->ops->destroy	       = PCDestroy_ICC;
2089b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
2099b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
21085317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
2119b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2129b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2139b54502bSHong Zhang 
214f8260c8fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor",
215f8260c8fSBarry Smith                     PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
2167112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
2177112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
21885317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
21985317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
220d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor",
221d90ac83dSHong Zhang                     PCFactorSetShiftType_Factor);CHKERRQ(ierr);
222d90ac83dSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor",
223d90ac83dSHong Zhang                     PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
22485317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
22585317021SBarry Smith                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
22685317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
22785317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
22885317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
22985317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
230174acd27SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor",
231174acd27SBarry Smith                     PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
232b7c853c4SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU",
233b7c853c4SBarry Smith                     PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
2349b54502bSHong Zhang   PetscFunctionReturn(0);
2359b54502bSHong Zhang }
2369b54502bSHong Zhang EXTERN_C_END
2379b54502bSHong Zhang 
2389b54502bSHong Zhang 
239