xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 00c67f3b09c0bcda06af5ed306d845d9138e5003)
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;
124cccfbddSHong Zhang   Mat            F;
13*00c67f3bSHong Zhang   const MatSolverPackage stype;
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 
354cccfbddSHong Zhang   F = ((PC_Factor*)icc)->fact;
364cccfbddSHong Zhang   if (F->errortype) { /* FactorSymbolic() fails */
374cccfbddSHong Zhang     pc->failedreason = (PCFailedReason)F->errortype;
386baea169SHong Zhang     PetscFunctionReturn(0);
396baea169SHong Zhang   }
406baea169SHong Zhang 
41075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
424cccfbddSHong Zhang   if (F->errortype) { /* FactorNumeric() fails */
434cccfbddSHong Zhang     pc->failedreason = (PCFailedReason)F->errortype;
446baea169SHong Zhang   }
45*00c67f3bSHong Zhang 
46*00c67f3bSHong Zhang     ierr = PCFactorGetMatSolverPackage(pc,&stype);CHKERRQ(ierr);
47*00c67f3bSHong Zhang   if (!stype) {
48*00c67f3bSHong Zhang     ierr = PCFactorSetMatSolverPackage(pc,((PC_Factor*)icc)->fact->solvertype);CHKERRQ(ierr);
49*00c67f3bSHong Zhang   }
509b54502bSHong Zhang   PetscFunctionReturn(0);
519b54502bSHong Zhang }
529b54502bSHong Zhang 
539b54502bSHong Zhang #undef __FUNCT__
54574deadeSBarry Smith #define __FUNCT__ "PCReset_ICC"
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 
65574deadeSBarry Smith #undef __FUNCT__
669b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
679b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
689b54502bSHong Zhang {
699b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
709b54502bSHong Zhang   PetscErrorCode ierr;
719b54502bSHong Zhang 
729b54502bSHong Zhang   PetscFunctionBegin;
73574deadeSBarry Smith   ierr = PCReset_ICC(pc);CHKERRQ(ierr);
74503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
75503cfb0cSBarry Smith   ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr);
76c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
779b54502bSHong Zhang   PetscFunctionReturn(0);
789b54502bSHong Zhang }
799b54502bSHong Zhang 
809b54502bSHong Zhang #undef __FUNCT__
819b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
829b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
839b54502bSHong Zhang {
849b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
859b54502bSHong Zhang   PetscErrorCode ierr;
869b54502bSHong Zhang 
879b54502bSHong Zhang   PetscFunctionBegin;
88075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
899b54502bSHong Zhang   PetscFunctionReturn(0);
909b54502bSHong Zhang }
919b54502bSHong Zhang 
929b54502bSHong Zhang #undef __FUNCT__
939b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
949b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
959b54502bSHong Zhang {
969b54502bSHong Zhang   PetscErrorCode ierr;
979b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
989b54502bSHong Zhang 
999b54502bSHong Zhang   PetscFunctionBegin;
100075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1019b54502bSHong Zhang   PetscFunctionReturn(0);
1029b54502bSHong Zhang }
1039b54502bSHong Zhang 
1049b54502bSHong Zhang #undef __FUNCT__
1059b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
1069b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
1079b54502bSHong Zhang {
1089b54502bSHong Zhang   PetscErrorCode ierr;
1099b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1109b54502bSHong Zhang 
1119b54502bSHong Zhang   PetscFunctionBegin;
112075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
1139b54502bSHong Zhang   PetscFunctionReturn(0);
1149b54502bSHong Zhang }
1159b54502bSHong Zhang 
1169b54502bSHong Zhang #undef __FUNCT__
1179b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
1184416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ICC(PetscOptionItems *PetscOptionsObject,PC pc)
1199b54502bSHong Zhang {
1209b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
121ace3abfcSBarry Smith   PetscBool      flg;
1229b54502bSHong Zhang   PetscErrorCode ierr;
1232fa5cd67SKarl Rupp   /* PetscReal      dt[3];*/
1249b54502bSHong Zhang 
1259b54502bSHong Zhang   PetscFunctionBegin;
126e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"ICC Options");CHKERRQ(ierr);
127e55864a3SBarry Smith   ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr);
1289b54502bSHong Zhang 
1298ff23777SHong Zhang   ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
1305a586d82SBarry Smith   /*dt[0] = ((PC_Factor*)icc)->info.dt;
1314c9036c7SBarry Smith   dt[1] = ((PC_Factor*)icc)->info.dtcol;
1324c9036c7SBarry Smith   dt[2] = ((PC_Factor*)icc)->info.dtcount;
13378fc6b22SHong Zhang   PetscInt       dtmax = 3;
134b7c853c4SBarry Smith   ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr);
1354c9036c7SBarry Smith   if (flg) {
136b7c853c4SBarry Smith     ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr);
1374c9036c7SBarry Smith   }
13878fc6b22SHong Zhang   */
1399b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1409b54502bSHong Zhang   PetscFunctionReturn(0);
1419b54502bSHong Zhang }
1429b54502bSHong Zhang 
1439b54502bSHong Zhang #undef __FUNCT__
1449b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1459b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1469b54502bSHong Zhang {
1479b54502bSHong Zhang   PetscErrorCode ierr;
1489b54502bSHong Zhang 
1499b54502bSHong Zhang   PetscFunctionBegin;
150914a5d51SHong Zhang   ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr);
1519b54502bSHong Zhang   PetscFunctionReturn(0);
1529b54502bSHong Zhang }
1539b54502bSHong Zhang 
1547087cfbeSBarry Smith extern PetscErrorCode  PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt);
1554c9036c7SBarry Smith 
156c60c7ad4SBarry Smith #undef __FUNCT__
157c60c7ad4SBarry Smith #define __FUNCT__ "PCFactorGetUseInPlace_ICC"
158c60c7ad4SBarry Smith PetscErrorCode  PCFactorGetUseInPlace_ICC(PC pc,PetscBool *flg)
159c60c7ad4SBarry Smith {
160c60c7ad4SBarry Smith   PetscFunctionBegin;
161c60c7ad4SBarry Smith   *flg = PETSC_FALSE;
162c60c7ad4SBarry Smith   PetscFunctionReturn(0);
163c60c7ad4SBarry Smith }
164c60c7ad4SBarry Smith 
1659b54502bSHong Zhang /*MC
1669b54502bSHong Zhang      PCICC - Incomplete Cholesky factorization preconditioners.
1679b54502bSHong Zhang 
1689b54502bSHong Zhang    Options Database Keys:
1692401956bSBarry Smith +  -pc_factor_levels <k> - number of levels of fill for ICC(k)
1702401956bSBarry Smith .  -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for
1719b54502bSHong Zhang                       its factorization (overwrites original matrix)
17255ba2a51SBarry Smith .  -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1
173145b9266SHong Zhang -  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
1749b54502bSHong Zhang 
1759b54502bSHong Zhang    Level: beginner
1769b54502bSHong Zhang 
1779b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1789b54502bSHong Zhang 
179d9ba8547SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel.
1809b54502bSHong Zhang 
1819b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1829b54502bSHong Zhang 
1839b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1849b54502bSHong Zhang 
18514d2772eSHong Zhang           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE);
1869b54502bSHong Zhang           to turn off the shift.
1879b54502bSHong Zhang 
188c582cd25SBarry Smith    References:
18996a0c994SBarry Smith .  1. - TONY F. CHAN AND HENK A. VAN DER VORST, Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS,
19096a0c994SBarry Smith       Chapter in Parallel Numerical Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
19196a0c994SBarry Smith       Science and Engineering, Kluwer.
192c582cd25SBarry Smith 
1939b54502bSHong Zhang 
1949b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
19581f96dccSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(),
196e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1978ff23777SHong Zhang            PCFactorSetLevels()
1989b54502bSHong Zhang 
1999b54502bSHong Zhang M*/
2009b54502bSHong Zhang 
2019b54502bSHong Zhang #undef __FUNCT__
2029b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
2038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ICC(PC pc)
2049b54502bSHong Zhang {
2059b54502bSHong Zhang   PetscErrorCode ierr;
2069b54502bSHong Zhang   PC_ICC         *icc;
2079b54502bSHong Zhang 
2089b54502bSHong Zhang   PetscFunctionBegin;
209b00a9115SJed Brown   ierr = PetscNewLog(pc,&icc);CHKERRQ(ierr);
2109b54502bSHong Zhang 
211075768bcSBarry Smith   ((PC_Factor*)icc)->fact = 0;
2122fa5cd67SKarl Rupp 
21319fd82e9SBarry Smith   ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
214075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
2152fa5cd67SKarl Rupp 
216879e8a4dSBarry Smith   ((PC_Factor*)icc)->factortype  = MAT_FACTOR_ICC;
21775567043SBarry Smith   ((PC_Factor*)icc)->info.levels = 0.;
218075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill   = 1.0;
2199b54502bSHong Zhang   icc->implctx                   = 0;
2209b54502bSHong Zhang 
221075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol       = PETSC_DEFAULT;
222f4db908eSBarry Smith   ((PC_Factor*)icc)->info.shifttype   = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE;
2230ed735ceSHong Zhang   ((PC_Factor*)icc)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON;
2240ed735ceSHong Zhang   ((PC_Factor*)icc)->info.zeropivot   = 100.0*PETSC_MACHINE_EPSILON;
2259b54502bSHong Zhang 
226915743fcSHong Zhang   pc->data                     = (void*)icc;
2279b54502bSHong Zhang   pc->ops->apply               = PCApply_ICC;
228d674540eSJed Brown   pc->ops->applytranspose      = PCApply_ICC;
2294cccfbddSHong Zhang   pc->ops->setup               = PCSetUp_ICC;
230574deadeSBarry Smith   pc->ops->reset               = PCReset_ICC;
2319b54502bSHong Zhang   pc->ops->destroy             = PCDestroy_ICC;
2329b54502bSHong Zhang   pc->ops->setfromoptions      = PCSetFromOptions_ICC;
2339b54502bSHong Zhang   pc->ops->view                = PCView_ICC;
23485317021SBarry Smith   pc->ops->getfactoredmatrix   = PCFactorGetMatrix_Factor;
2359b54502bSHong Zhang   pc->ops->applysymmetricleft  = PCApplySymmetricLeft_ICC;
2369b54502bSHong Zhang   pc->ops->applysymmetricright = PCApplySymmetricRight_ICC;
2379b54502bSHong Zhang 
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr);
239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
240bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr);
2442591b318SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr);
245bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr);
246bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
247bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr);
248bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr);
249c60c7ad4SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_ICC);CHKERRQ(ierr);
2509b54502bSHong Zhang   PetscFunctionReturn(0);
2519b54502bSHong Zhang }
2529b54502bSHong Zhang 
2539b54502bSHong Zhang 
254