xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 00e125f844f1e8da6dd97b92f2bb419dde0a009f)
1dba47a55SKris Buschelman 
2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/icc/icc.h>   /*I "petscpc.h" I*/
39b54502bSHong Zhang 
49b54502bSHong Zhang #undef __FUNCT__
54cccfbddSHong Zhang #define __FUNCT__ "PCSetUp_ICC"
64cccfbddSHong 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;
1200c67f3bSHong Zhang   const MatSolverPackage stype;
13*00e125f8SBarry Smith   MatFactorError         err;
149b54502bSHong Zhang 
159b54502bSHong Zhang   PetscFunctionBegin;
16075768bcSBarry Smith   ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr);
179b54502bSHong Zhang 
1884d44b13SHong Zhang   ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr);
199b54502bSHong Zhang   if (!pc->setupcalled) {
20d09a07f4SBarry Smith     if (!((PC_Factor*)icc)->fact) {
21174acd27SBarry Smith       ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
22a1f19f5aSHong Zhang     }
23075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
249b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
256bf464f9SBarry Smith     ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
2614798fb4SJed Brown     ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
27075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
289b54502bSHong Zhang   }
29075768bcSBarry Smith   ierr            = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
30f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
31f3a39becSBarry Smith 
32fcfd50ebSBarry Smith   ierr = ISDestroy(&cperm);CHKERRQ(ierr);
33fcfd50ebSBarry Smith   ierr = ISDestroy(&perm);CHKERRQ(ierr);
346baea169SHong Zhang 
35*00e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
36*00e125f8SBarry Smith   if (err) { /* FactorSymbolic() fails */
37*00e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
386baea169SHong Zhang     PetscFunctionReturn(0);
396baea169SHong Zhang   }
406baea169SHong Zhang 
41075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
42*00e125f8SBarry Smith   ierr = MatFactorGetError(((PC_Factor*)icc)->fact,&err);CHKERRQ(ierr);
43*00e125f8SBarry Smith   if (err) { /* FactorNumeric() fails */
44*00e125f8SBarry Smith     pc->failedreason = (PCFailedReason)err;
456baea169SHong Zhang   }
4600c67f3bSHong Zhang 
4700c67f3bSHong Zhang   ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
4800c67f3bSHong Zhang   if (!stype) {
49*00e125f8SBarry Smith     const MatSolverPackage solverpackage;
50*00e125f8SBarry Smith     ierr = MatFactorGetSolverPackage(((PC_Factor*)icc)->fact,&solverpackage);CHKERRQ(ierr);
51*00e125f8SBarry Smith     ierr = PCFactorSetMatSolverPackage(pc,solverpackage);CHKERRQ(ierr);
5200c67f3bSHong Zhang   }
539b54502bSHong Zhang   PetscFunctionReturn(0);
549b54502bSHong Zhang }
559b54502bSHong Zhang 
569b54502bSHong Zhang #undef __FUNCT__
57574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC"
58574deadeSBarry Smith static PetscErrorCode PCReset_ICC(PC pc)
59574deadeSBarry Smith {
60574deadeSBarry Smith   PC_ICC         *icc = (PC_ICC*)pc->data;
61574deadeSBarry Smith   PetscErrorCode ierr;
62574deadeSBarry Smith 
63574deadeSBarry Smith   PetscFunctionBegin;
646bf464f9SBarry Smith   ierr = MatDestroy(&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
65574deadeSBarry Smith   PetscFunctionReturn(0);
66574deadeSBarry Smith }
67574deadeSBarry Smith 
68574deadeSBarry Smith #undef __FUNCT__
699b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
709b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
719b54502bSHong Zhang {
729b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
739b54502bSHong Zhang   PetscErrorCode ierr;
749b54502bSHong Zhang 
759b54502bSHong Zhang   PetscFunctionBegin;
76574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
77503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
78503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
79c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
809b54502bSHong Zhang   PetscFunctionReturn(0);
819b54502bSHong Zhang }
829b54502bSHong Zhang 
839b54502bSHong Zhang #undef __FUNCT__
849b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
859b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
869b54502bSHong Zhang {
879b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
889b54502bSHong Zhang   PetscErrorCode ierr;
899b54502bSHong Zhang 
909b54502bSHong Zhang   PetscFunctionBegin;
91075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
929b54502bSHong Zhang   PetscFunctionReturn(0);
939b54502bSHong Zhang }
949b54502bSHong Zhang 
959b54502bSHong Zhang #undef __FUNCT__
969b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
979b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
989b54502bSHong Zhang {
999b54502bSHong Zhang   PetscErrorCode ierr;
1009b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1019b54502bSHong Zhang 
1029b54502bSHong Zhang   PetscFunctionBegin;
103075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1049b54502bSHong Zhang   PetscFunctionReturn(0);
1059b54502bSHong Zhang }
1069b54502bSHong Zhang 
1079b54502bSHong Zhang #undef __FUNCT__
1089b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
1099b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_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 = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1169b54502bSHong Zhang   PetscFunctionReturn(0);
1179b54502bSHong Zhang }
1189b54502bSHong Zhang 
1199b54502bSHong Zhang #undef __FUNCT__
1209b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1214416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1229b54502bSHong Zhang {
1239b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
124ace3abfcSBarry Smith   PetscBool      flg;
1259b54502bSHong Zhang   PetscErrorCode ierr;
1262fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1279b54502bSHong Zhang 
1289b54502bSHong Zhang   PetscFunctionBegin;
129e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
130e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1319b54502bSHong Zhang 
1328ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1335a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1344c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1354c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
13678fc6b22SHong Zhang   PetscInt       dtmax = 3;
137b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1384c9036c7SBarry Smith   if (flg) {
139b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1404c9036c7SBarry Smith   }
14178fc6b22SHong Zhang   */
1429b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1439b54502bSHong Zhang   PetscFunctionReturn(0);
1449b54502bSHong Zhang }
1459b54502bSHong Zhang 
1469b54502bSHong Zhang #undef __FUNCT__
1479b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1489b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1499b54502bSHong Zhang {
1509b54502bSHong Zhang   PetscErrorCode ierr;
1519b54502bSHong Zhang 
1529b54502bSHong Zhang   PetscFunctionBegin;
153914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1549b54502bSHong Zhang   PetscFunctionReturn(0);
1559b54502bSHong Zhang }
1569b54502bSHong Zhang 
1577087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1584c9036c7SBarry Smith 
159c60c7ad4SBarry Smith #undef __FUNCT__
160c60c7ad4SBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_ICC"
161c60c7ad4SBarry Smith PetscErrorCode  PCFactorGetUseInPlace_ICC(PC pc,PetscBool *flg)
162c60c7ad4SBarry Smith {
163c60c7ad4SBarry Smith   PetscFunctionBegin;
164c60c7ad4SBarry Smith   *flg = PETSC_FALSE;
165c60c7ad4SBarry Smith   PetscFunctionReturn(0);
166c60c7ad4SBarry Smith }
167c60c7ad4SBarry Smith 
1689b54502bSHong Zhang /*MC
1699b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1709b54502bSHong Zhang 
1719b54502bSHong Zhang    Options Database Keys:
1722401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1732401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1749b54502bSHong Zhang                       its factorization (overwrites original matrix)
17555ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
176145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1779b54502bSHong Zhang 
1789b54502bSHong Zhang    Level: beginner
1799b54502bSHong Zhang 
1809b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1819b54502bSHong Zhang 
182d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel.
1839b54502bSHong Zhang 
1849b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1859b54502bSHong Zhang 
1869b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1879b54502bSHong Zhang 
18814d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1899b54502bSHong Zhang           to turn off the shift.
1909b54502bSHong Zhang 
191c582cd25SBarry Smith    References:
19296a0c994SBarry Smith .  1. - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
19396a0c994SBarry Smith       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
19496a0c994SBarry Smith       Science and Engineering, Kluwer.
195c582cd25SBarry Smith 
1969b54502bSHong Zhang 
1979b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
19881f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
199e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
2008ff23777SHong Zhang            PCFactorSetLevels()
2019b54502bSHong Zhang 
2029b54502bSHong Zhang M*/
2039b54502bSHong Zhang 
2049b54502bSHong Zhang #undef __FUNCT__
2059b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
2068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
2079b54502bSHong Zhang {
2089b54502bSHong Zhang   PetscErrorCode ierr;
2099b54502bSHong Zhang   PC_ICC         *icc;
2109b54502bSHong Zhang 
2119b54502bSHong Zhang   PetscFunctionBegin;
212b00a9115SJed Brown   ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr);
2139b54502bSHong Zhang 
214075768bcSBarry Smith   ((PC_Factor*)icc)->fact = 0;
2152fa5cd67SKarl Rupp 
21619fd82e9SBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
217075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
2182fa5cd67SKarl Rupp 
219879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype  = MAT_FACTOR_ICC;
22075567043SBarry Smith   ((PC_Factor*)icc)->info.levels = 0.;
221075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill   = 1.0;
2229b54502bSHong Zhang   icc->implctx                   = 0;
2239b54502bSHong Zhang 
224075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
225f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2260ed735ceSHong Zhang   ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
2270ed735ceSHong Zhang   ((PC_Factor*)icc)->info.zeropivot   = 100.0*PETSC_MACHINE_EPSILON;
2289b54502bSHong Zhang 
229915743fcSHong Zhang   pc->data                     = (void*)icc;
2309b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
231d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2324cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
233574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2349b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2359b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
2369b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
23785317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
2389b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2399b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2409b54502bSHong Zhang 
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
244bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
245bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
246bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
2472591b318SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
248bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
249bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
250bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
251bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
252c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ICC);CHKERRQ(ierr);
2539b54502bSHong Zhang   PetscFunctionReturn(0);
2549b54502bSHong Zhang }
2559b54502bSHong Zhang 
2569b54502bSHong Zhang 
257