1 #define PETSCKSP_DLL 2 3 #include "../src/ksp/pc/impls/factor/icc/icc.h" /*I "petscpc.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "PCSetup_ICC" 7 static PetscErrorCode PCSetup_ICC(PC pc) 8 { 9 PC_ICC *icc = (PC_ICC*)pc->data; 10 IS perm,cperm; 11 PetscErrorCode ierr; 12 MatInfo info; 13 14 PetscFunctionBegin; 15 ierr = MatGetOrdering(pc->pmat, ((PC_Factor*)icc)->ordering,&perm,&cperm);CHKERRQ(ierr); 16 17 if (!pc->setupcalled) { 18 if (!((PC_Factor*)icc)->fact){ 19 ierr = MatGetFactor(pc->pmat,((PC_Factor*)icc)->solvertype,MAT_FACTOR_ICC,& ((PC_Factor*)icc)->fact);CHKERRQ(ierr); 20 } 21 ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 22 } else if (pc->flag != SAME_NONZERO_PATTERN) { 23 ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr); 24 ierr = MatGetFactor(pc->pmat,MATSOLVERPETSC,MAT_FACTOR_ICC,&((PC_Factor*)icc)->fact);CHKERRQ(ierr); 25 ierr = MatICCFactorSymbolic(((PC_Factor*)icc)->fact,pc->pmat,perm,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 26 } 27 ierr = MatGetInfo(((PC_Factor*)icc)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 28 icc->actualfill = info.fill_ratio_needed; 29 30 ierr = ISDestroy(cperm);CHKERRQ(ierr); 31 ierr = ISDestroy(perm);CHKERRQ(ierr); 32 ierr = MatCholeskyFactorNumeric(((PC_Factor*)icc)->fact,pc->pmat,&((PC_Factor*)icc)->info);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "PCDestroy_ICC" 38 static PetscErrorCode PCDestroy_ICC(PC pc) 39 { 40 PC_ICC *icc = (PC_ICC*)pc->data; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 if (((PC_Factor*)icc)->fact) {ierr = MatDestroy(((PC_Factor*)icc)->fact);CHKERRQ(ierr);} 45 ierr = PetscFree(((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 46 ierr = PetscFree(((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 47 ierr = PetscFree(icc);CHKERRQ(ierr); 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "PCApply_ICC" 53 static PetscErrorCode PCApply_ICC(PC pc,Vec x,Vec y) 54 { 55 PC_ICC *icc = (PC_ICC*)pc->data; 56 PetscErrorCode ierr; 57 58 PetscFunctionBegin; 59 ierr = MatSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "PCApplySymmetricLeft_ICC" 65 static PetscErrorCode PCApplySymmetricLeft_ICC(PC pc,Vec x,Vec y) 66 { 67 PetscErrorCode ierr; 68 PC_ICC *icc = (PC_ICC*)pc->data; 69 70 PetscFunctionBegin; 71 ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNCT__ 76 #define __FUNCT__ "PCApplySymmetricRight_ICC" 77 static PetscErrorCode PCApplySymmetricRight_ICC(PC pc,Vec x,Vec y) 78 { 79 PetscErrorCode ierr; 80 PC_ICC *icc = (PC_ICC*)pc->data; 81 82 PetscFunctionBegin; 83 ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "PCSetFromOptions_ICC" 89 static PetscErrorCode PCSetFromOptions_ICC(PC pc) 90 { 91 PC_ICC *icc = (PC_ICC*)pc->data; 92 PetscTruth flg; 93 /* PetscReal dt[3];*/ 94 PetscErrorCode ierr; 95 96 PetscFunctionBegin; 97 ierr = PetscOptionsHead("ICC Options");CHKERRQ(ierr); 98 ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 99 100 ierr = PetscOptionsReal("-pc_factor_levels","levels of fill","PCFactorSetLevels",((PC_Factor*)icc)->info.levels,&((PC_Factor*)icc)->info.levels,&flg);CHKERRQ(ierr); 101 /*dt[0] = ((PC_Factor*)icc)->info.dt; 102 dt[1] = ((PC_Factor*)icc)->info.dtcol; 103 dt[2] = ((PC_Factor*)icc)->info.dtcount; 104 PetscInt dtmax = 3; 105 ierr = PetscOptionsRealArray("-pc_factor_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 106 if (flg) { 107 ierr = PCFactorSetDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 108 } 109 */ 110 ierr = PetscOptionsTail();CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 #undef __FUNCT__ 115 #define __FUNCT__ "PCView_ICC" 116 static PetscErrorCode PCView_ICC(PC pc,PetscViewer viewer) 117 { 118 PetscErrorCode ierr; 119 120 PetscFunctionBegin; 121 ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 122 PetscFunctionReturn(0); 123 } 124 125 EXTERN_C_BEGIN 126 extern PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetDropTolerance_ILU(PC,PetscReal,PetscReal,PetscInt); 127 EXTERN_C_END 128 129 /*MC 130 PCICC - Incomplete Cholesky factorization preconditioners. 131 132 Options Database Keys: 133 + -pc_factor_levels <k> - number of levels of fill for ICC(k) 134 . -pc_factor_in_place - only for ICC(0) with natural ordering, reuses the space of the matrix for 135 its factorization (overwrites original matrix) 136 . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 137 - -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 138 139 Level: beginner 140 141 Concepts: incomplete Cholesky factorization 142 143 Notes: Only implemented for some matrix formats. Not implemented in parallel. 144 145 For BAIJ matrices this implements a point block ICC. 146 147 The Manteuffel shift is only implemented for matrices with block size 1 148 149 By default, the Manteuffel is applied (for matrices with block size 1). Call PCFactorSetShiftType(pc,MAT_SHIFT_POSITIVE_DEFINITE); 150 to turn off the shift. 151 152 References: 153 Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 154 http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 155 Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 156 Science and Engineering, Kluwer, pp. 167--202. 157 158 159 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 160 PCFactorSetZeroPivot(), PCFactorSetShiftType(), PCFactorSetShiftAmount(), 161 PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 162 PCFactorSetLevels() 163 164 M*/ 165 166 EXTERN_C_BEGIN 167 #undef __FUNCT__ 168 #define __FUNCT__ "PCCreate_ICC" 169 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ICC(PC pc) 170 { 171 PetscErrorCode ierr; 172 PC_ICC *icc; 173 174 PetscFunctionBegin; 175 ierr = PetscNewLog(pc,PC_ICC,&icc);CHKERRQ(ierr); 176 177 ((PC_Factor*)icc)->fact = 0; 178 ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)icc)->ordering);CHKERRQ(ierr); 179 ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)icc)->solvertype);CHKERRQ(ierr); 180 ierr = MatFactorInfoInitialize(&((PC_Factor*)icc)->info);CHKERRQ(ierr); 181 ((PC_Factor*)icc)->factortype = MAT_FACTOR_ICC; 182 ((PC_Factor*)icc)->info.levels = 0.; 183 ((PC_Factor*)icc)->info.fill = 1.0; 184 icc->implctx = 0; 185 186 ((PC_Factor*)icc)->info.dtcol = PETSC_DEFAULT; 187 ((PC_Factor*)icc)->info.shifttype = (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE; 188 ((PC_Factor*)icc)->info.shiftamount = 1.e-12; 189 ((PC_Factor*)icc)->info.zeropivot = 1.e-12; 190 191 pc->data = (void*)icc; 192 pc->ops->apply = PCApply_ICC; 193 pc->ops->setup = PCSetup_ICC; 194 pc->ops->destroy = PCDestroy_ICC; 195 pc->ops->setfromoptions = PCSetFromOptions_ICC; 196 pc->ops->view = PCView_ICC; 197 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 198 pc->ops->applysymmetricleft = PCApplySymmetricLeft_ICC; 199 pc->ops->applysymmetricright = PCApplySymmetricRight_ICC; 200 201 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 202 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 203 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 204 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 205 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 206 PCFactorSetShiftType_Factor);CHKERRQ(ierr); 207 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 208 PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 209 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_Factor", 210 PCFactorSetLevels_Factor);CHKERRQ(ierr); 211 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 212 PCFactorSetFill_Factor);CHKERRQ(ierr); 213 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 214 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 215 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 216 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 217 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetDropTolerance_C","PCFactorSetDropTolerance_ILU", 218 PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 219 PetscFunctionReturn(0); 220 } 221 EXTERN_C_END 222 223 224