xref: /petsc/src/ksp/pc/impls/factor/icc/icc.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
3*7c4f633dSBarry 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) {
18075768bcSBarry Smith     ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr);
19075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
209b54502bSHong Zhang   } else if (pc->flag != SAME_NONZERO_PATTERN) {
21075768bcSBarry Smith     ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);
22075768bcSBarry Smith     ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr);
23075768bcSBarry Smith     ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
249b54502bSHong Zhang   }
25075768bcSBarry Smith   ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr);
26f3a39becSBarry Smith   icc->actualfill = info.fill_ratio_needed;
27f3a39becSBarry Smith 
289b54502bSHong Zhang   ierr = ISDestroy(cperm);CHKERRQ(ierr);
299b54502bSHong Zhang   ierr = ISDestroy(perm);CHKERRQ(ierr);
30075768bcSBarry Smith   ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr);
319b54502bSHong Zhang   PetscFunctionReturn(0);
329b54502bSHong Zhang }
339b54502bSHong Zhang 
349b54502bSHong Zhang #undef __FUNCT__
359b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ICC"
369b54502bSHong Zhang static PetscErrorCode PCDestroy_ICC(PC pc)
379b54502bSHong Zhang {
389b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
399b54502bSHong Zhang   PetscErrorCode ierr;
409b54502bSHong Zhang 
419b54502bSHong Zhang   PetscFunctionBegin;
42075768bcSBarry Smith   if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);}
43075768bcSBarry Smith   ierr = PetscStrfree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
449b54502bSHong Zhang   ierr = PetscFree(icc);CHKERRQ(ierr);
459b54502bSHong Zhang   PetscFunctionReturn(0);
469b54502bSHong Zhang }
479b54502bSHong Zhang 
489b54502bSHong Zhang #undef __FUNCT__
499b54502bSHong Zhang #define __FUNCT__ "PCApply_ICC"
509b54502bSHong Zhang static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y)
519b54502bSHong Zhang {
529b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
539b54502bSHong Zhang   PetscErrorCode ierr;
549b54502bSHong Zhang 
559b54502bSHong Zhang   PetscFunctionBegin;
56075768bcSBarry Smith   ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
579b54502bSHong Zhang   PetscFunctionReturn(0);
589b54502bSHong Zhang }
599b54502bSHong Zhang 
609b54502bSHong Zhang #undef __FUNCT__
619b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricLeft_ICC"
629b54502bSHong Zhang static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y)
639b54502bSHong Zhang {
649b54502bSHong Zhang   PetscErrorCode ierr;
659b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
669b54502bSHong Zhang 
679b54502bSHong Zhang   PetscFunctionBegin;
68075768bcSBarry Smith   ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
699b54502bSHong Zhang   PetscFunctionReturn(0);
709b54502bSHong Zhang }
719b54502bSHong Zhang 
729b54502bSHong Zhang #undef __FUNCT__
739b54502bSHong Zhang #define __FUNCT__ "PCApplySymmetricRight_ICC"
749b54502bSHong Zhang static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y)
759b54502bSHong Zhang {
769b54502bSHong Zhang   PetscErrorCode ierr;
779b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
789b54502bSHong Zhang 
799b54502bSHong Zhang   PetscFunctionBegin;
80075768bcSBarry Smith   ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr);
819b54502bSHong Zhang   PetscFunctionReturn(0);
829b54502bSHong Zhang }
839b54502bSHong Zhang 
849b54502bSHong Zhang #undef __FUNCT__
859b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ICC"
869b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ICC(PC pc)
879b54502bSHong Zhang {
889b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
899b54502bSHong Zhang   char           tname[256];
909b54502bSHong Zhang   PetscTruth     flg;
919b54502bSHong Zhang   PetscErrorCode ierr;
929b54502bSHong Zhang   PetscFList     ordlist;
939b54502bSHong Zhang 
949b54502bSHong Zhang   PetscFunctionBegin;
959b54502bSHong Zhang   ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr);
969b54502bSHong Zhang   ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr);
97075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr);
98075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)icc)->info.fill,&((PC_Factor*)icc)->info.fill,&flg);CHKERRQ(ierr);
999b54502bSHong Zhang     ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr);
100075768bcSBarry Smith     ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ICC","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)icc)->ordering,tname,256,&flg);CHKERRQ(ierr);
1019b54502bSHong Zhang     if (flg) {
102e5a9bf91SBarry Smith       ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr);
1039b54502bSHong Zhang     }
1049f95998fSHong Zhang     ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr);
1059b54502bSHong Zhang     if (flg) {
106afaefe49SHong Zhang       ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr);
1079b54502bSHong Zhang     }
108075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)icc)->info.shiftnz,&((PC_Factor*)icc)->info.shiftnz,0);CHKERRQ(ierr);
109075768bcSBarry Smith     flg = (((PC_Factor*)icc)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE;
11062bba022SBarry Smith     ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
11162bba022SBarry Smith     ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr);
112075768bcSBarry Smith     ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)icc)->info.zeropivot,&((PC_Factor*)icc)->info.zeropivot,0);CHKERRQ(ierr);
1139b54502bSHong Zhang 
1149b54502bSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
1159b54502bSHong Zhang   PetscFunctionReturn(0);
1169b54502bSHong Zhang }
1179b54502bSHong Zhang 
1189b54502bSHong Zhang #undef __FUNCT__
1199b54502bSHong Zhang #define __FUNCT__ "PCView_ICC"
1209b54502bSHong Zhang static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer)
1219b54502bSHong Zhang {
1229b54502bSHong Zhang   PC_ICC         *icc = (PC_ICC*)pc->data;
1239b54502bSHong Zhang   PetscErrorCode ierr;
1249b54502bSHong Zhang   PetscTruth     isstring,iascii;
1259b54502bSHong Zhang 
1269b54502bSHong Zhang   PetscFunctionBegin;
1279b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
1289b54502bSHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1299b54502bSHong Zhang   if (iascii) {
130075768bcSBarry Smith     if (((PC_Factor*)icc)->info.levels == 1) {
131075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D level of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
1329b54502bSHong Zhang     } else {
133075768bcSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  ICC: %D levels of fill\n",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);
1349b54502bSHong Zhang     }
135075768bcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio allocated %G, ordering used %s\n",((PC_Factor*)icc)->info.fill,((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
136075768bcSBarry Smith     if (((PC_Factor*)icc)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using Manteuffel shift\n");CHKERRQ(ierr);}
137075768bcSBarry Smith     if (((PC_Factor*)icc)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer,"  ICC: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);}
138075768bcSBarry Smith     if (((PC_Factor*)icc)->fact) {
139f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  ICC: factor fill ratio needed %G\n",icc->actualfill);CHKERRQ(ierr);
140f3a39becSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"       Factored matrix follows\n");CHKERRQ(ierr);
141f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
142f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
143f3a39becSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
144f3a39becSBarry Smith       ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
145075768bcSBarry Smith       ierr = MatView(((PC_Factor*)icc)->fact,viewer);CHKERRQ(ierr);
146f3a39becSBarry Smith       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
147f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
148f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
149f3a39becSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
150f3a39becSBarry Smith     }
1519b54502bSHong Zhang   } else if (isstring) {
152075768bcSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," lvls=%D",(PetscInt)((PC_Factor*)icc)->info.levels);CHKERRQ(ierr);CHKERRQ(ierr);
1539b54502bSHong Zhang   } else {
1549b54502bSHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCICC",((PetscObject)viewer)->type_name);
1559b54502bSHong Zhang   }
1569b54502bSHong Zhang   PetscFunctionReturn(0);
1579b54502bSHong Zhang }
1589b54502bSHong Zhang 
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
1672401956bSBarry Smith .  -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix
168f251bdbdSHong Zhang .  -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default
169f251bdbdSHong Zhang -  -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value
170f251bdbdSHong Zhang    is optional with PETSC_TRUE being the default
1719b54502bSHong Zhang 
1729b54502bSHong Zhang    Level: beginner
1739b54502bSHong Zhang 
1749b54502bSHong Zhang   Concepts: incomplete Cholesky factorization
1759b54502bSHong Zhang 
176c582cd25SBarry Smith    Notes: Only implemented for some matrix formats. Not implemented in parallel (for parallel use you
177c582cd25SBarry Smith              must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ICC(0) and this is not recommended
178c582cd25SBarry Smith              unless you really want a parallel ICC).
1799b54502bSHong Zhang 
1809b54502bSHong Zhang           For BAIJ matrices this implements a point block ICC.
1819b54502bSHong Zhang 
1829b54502bSHong Zhang           The Manteuffel shift is only implemented for matrices with block size 1
1839b54502bSHong Zhang 
1842401956bSBarry Smith           By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftPd(pc,PETSC_FALSE);
1859b54502bSHong Zhang           to turn off the shift.
1869b54502bSHong Zhang 
187c582cd25SBarry Smith    References:
188c582cd25SBarry Smith    Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST
189c582cd25SBarry Smith       http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical
190c582cd25SBarry Smith       Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in
191c582cd25SBarry Smith       Science and Engineering, Kluwer, pp. 167--202.
192c582cd25SBarry Smith 
1939b54502bSHong Zhang 
1949b54502bSHong Zhang .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType,
195ee45ca4aSHong Zhang            PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(),
196e5a9bf91SBarry Smith            PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(),
1972401956bSBarry Smith            PCFactorSetLevels(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(),
1989b54502bSHong Zhang 
1999b54502bSHong Zhang M*/
2009b54502bSHong Zhang 
2019b54502bSHong Zhang EXTERN_C_BEGIN
2029b54502bSHong Zhang #undef __FUNCT__
2039b54502bSHong Zhang #define __FUNCT__ "PCCreate_ICC"
204dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc)
2059b54502bSHong Zhang {
2069b54502bSHong Zhang   PetscErrorCode ierr;
2079b54502bSHong Zhang   PC_ICC         *icc;
2089b54502bSHong Zhang 
2099b54502bSHong Zhang   PetscFunctionBegin;
21038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr);
2119b54502bSHong Zhang 
212075768bcSBarry Smith   ((PC_Factor*)icc)->fact	          = 0;
213075768bcSBarry Smith   ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr);
214075768bcSBarry Smith   ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr);
215075768bcSBarry Smith   ((PC_Factor*)icc)->info.levels	  = 0;
216075768bcSBarry Smith   ((PC_Factor*)icc)->info.fill          = 1.0;
2179b54502bSHong Zhang   icc->implctx            = 0;
2189b54502bSHong Zhang 
219075768bcSBarry Smith   ((PC_Factor*)icc)->info.dtcol              = PETSC_DEFAULT;
220075768bcSBarry Smith   ((PC_Factor*)icc)->info.shiftnz            = 0.0;
221075768bcSBarry Smith   ((PC_Factor*)icc)->info.shiftpd            = 1.0; /* true */
222075768bcSBarry Smith   ((PC_Factor*)icc)->info.zeropivot          = 1.e-12;
2239b54502bSHong Zhang   pc->data	               = (void*)icc;
2249b54502bSHong Zhang 
2259b54502bSHong Zhang   pc->ops->apply	       = PCApply_ICC;
2269b54502bSHong Zhang   pc->ops->setup               = PCSetup_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 
2347112b564SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor",
2357112b564SBarry Smith                     PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr);
23685317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor",
23785317021SBarry Smith                     PCFactorSetZeroPivot_Factor);CHKERRQ(ierr);
23885317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor",
23985317021SBarry Smith                     PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr);
24085317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor",
24185317021SBarry Smith                     PCFactorSetShiftPd_Factor);CHKERRQ(ierr);
242afaefe49SHong Zhang 
24385317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor",
24485317021SBarry Smith                     PCFactorSetLevels_Factor);CHKERRQ(ierr);
24585317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor",
24685317021SBarry Smith                     PCFactorSetFill_Factor);CHKERRQ(ierr);
24785317021SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor",
24885317021SBarry Smith                     PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr);
2499b54502bSHong Zhang   PetscFunctionReturn(0);
2509b54502bSHong Zhang }
2519b54502bSHong Zhang EXTERN_C_END
2529b54502bSHong Zhang 
2539b54502bSHong Zhang 
254