1 2 /* 3 Defines a direct factorization preconditioner for any Mat implementation 4 Note: this need not be consided a preconditioner since it supplies 5 a direct solver. 6 */ 7 8 #include "../src/ksp/pc/impls/factor/lu/lu.h" /*I "petscpc.h" I*/ 9 10 EXTERN_C_BEGIN 11 #undef __FUNCT__ 12 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 13 PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 14 { 15 PC_LU *lu = (PC_LU*)pc->data; 16 17 PetscFunctionBegin; 18 lu->nonzerosalongdiagonal = PETSC_TRUE; 19 if (z == PETSC_DECIDE) { 20 lu->nonzerosalongdiagonaltol = 1.e-10; 21 } else { 22 lu->nonzerosalongdiagonaltol = z; 23 } 24 PetscFunctionReturn(0); 25 } 26 EXTERN_C_END 27 28 EXTERN_C_BEGIN 29 #undef __FUNCT__ 30 #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 31 PetscErrorCode PCFactorSetReuseOrdering_LU(PC pc,PetscBool flag) 32 { 33 PC_LU *lu = (PC_LU*)pc->data; 34 35 PetscFunctionBegin; 36 lu->reuseordering = flag; 37 PetscFunctionReturn(0); 38 } 39 EXTERN_C_END 40 41 EXTERN_C_BEGIN 42 #undef __FUNCT__ 43 #define __FUNCT__ "PCFactorSetReuseFill_LU" 44 PetscErrorCode PCFactorSetReuseFill_LU(PC pc,PetscBool flag) 45 { 46 PC_LU *lu = (PC_LU*)pc->data; 47 48 PetscFunctionBegin; 49 lu->reusefill = flag; 50 PetscFunctionReturn(0); 51 } 52 EXTERN_C_END 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "PCSetFromOptions_LU" 56 static PetscErrorCode PCSetFromOptions_LU(PC pc) 57 { 58 PC_LU *lu = (PC_LU*)pc->data; 59 PetscErrorCode ierr; 60 PetscBool flg = PETSC_FALSE; 61 PetscReal tol; 62 63 PetscFunctionBegin; 64 ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 65 ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 66 67 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 68 if (flg) { 69 tol = PETSC_DECIDE; 70 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 71 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 72 } 73 ierr = PetscOptionsTail();CHKERRQ(ierr); 74 PetscFunctionReturn(0); 75 } 76 77 #undef __FUNCT__ 78 #define __FUNCT__ "PCView_LU" 79 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 80 { 81 PC_LU *lu = (PC_LU*)pc->data; 82 PetscErrorCode ierr; 83 PetscBool iascii; 84 85 PetscFunctionBegin; 86 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 87 if (iascii) { 88 if (lu->inplace) { 89 ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr); 90 } else { 91 ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr); 92 } 93 94 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 95 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 96 } 97 ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "PCSetUp_LU" 103 static PetscErrorCode PCSetUp_LU(PC pc) 104 { 105 PetscErrorCode ierr; 106 PC_LU *dir = (PC_LU*)pc->data; 107 108 PetscFunctionBegin; 109 if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 110 111 if (dir->inplace) { 112 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 113 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 114 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 115 if (dir->row) { 116 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 117 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 118 } 119 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 120 ((PC_Factor*)dir)->fact = pc->pmat; 121 } else { 122 MatInfo info; 123 if (!pc->setupcalled) { 124 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 125 if (dir->nonzerosalongdiagonal) { 126 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 127 } 128 if (dir->row) { 129 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 130 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 131 } 132 if (!((PC_Factor*)dir)->fact){ 133 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 134 } 135 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 136 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 137 dir->actualfill = info.fill_ratio_needed; 138 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 139 } else if (pc->flag != SAME_NONZERO_PATTERN) { 140 if (!dir->reuseordering) { 141 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 142 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 143 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 144 if (dir->nonzerosalongdiagonal) { 145 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 146 } 147 if (dir->row) { 148 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 149 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 150 } 151 } 152 ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 153 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 154 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 155 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 156 dir->actualfill = info.fill_ratio_needed; 157 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 158 } 159 ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 160 } 161 PetscFunctionReturn(0); 162 } 163 164 #undef __FUNCT__ 165 #define __FUNCT__ "PCReset_LU" 166 static PetscErrorCode PCReset_LU(PC pc) 167 { 168 PC_LU *dir = (PC_LU*)pc->data; 169 PetscErrorCode ierr; 170 171 PetscFunctionBegin; 172 if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 173 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 174 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "PCDestroy_LU" 180 static PetscErrorCode PCDestroy_LU(PC pc) 181 { 182 PC_LU *dir = (PC_LU*)pc->data; 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 ierr = PCReset_LU(pc);CHKERRQ(ierr); 187 ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 188 ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 189 ierr = PetscFree(pc->data);CHKERRQ(ierr); 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "PCApply_LU" 195 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 196 { 197 PC_LU *dir = (PC_LU*)pc->data; 198 PetscErrorCode ierr; 199 200 PetscFunctionBegin; 201 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 202 else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "PCApplyTranspose_LU" 208 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 209 { 210 PC_LU *dir = (PC_LU*)pc->data; 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 215 else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 216 PetscFunctionReturn(0); 217 } 218 219 /* -----------------------------------------------------------------------------------*/ 220 221 EXTERN_C_BEGIN 222 #undef __FUNCT__ 223 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 224 PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 225 { 226 PC_LU *dir = (PC_LU*)pc->data; 227 228 PetscFunctionBegin; 229 dir->inplace = PETSC_TRUE; 230 PetscFunctionReturn(0); 231 } 232 EXTERN_C_END 233 234 /* ------------------------------------------------------------------------ */ 235 236 /*MC 237 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 238 239 Options Database Keys: 240 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 241 . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 242 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 243 . -pc_factor_fill <fill> - Sets fill amount 244 . -pc_factor_in_place - Activates in-place factorization 245 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 246 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 247 stability of factorization. 248 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 249 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 250 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 251 diagonal. 252 253 Notes: Not all options work for all matrix formats 254 Run with -help to see additional options for particular matrix formats or factorization 255 algorithms 256 257 Level: beginner 258 259 Concepts: LU factorization, direct solver 260 261 Notes: Usually this will compute an "exact" solution in one iteration and does 262 not need a Krylov method (i.e. you can use -ksp_type preonly, or 263 KSPSetType(ksp,KSPPREONLY) for the Krylov method 264 265 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 266 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 267 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 268 PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 269 PCFactorReorderForNonzeroDiagonal() 270 M*/ 271 272 EXTERN_C_BEGIN 273 #undef __FUNCT__ 274 #define __FUNCT__ "PCCreate_LU" 275 PetscErrorCode PCCreate_LU(PC pc) 276 { 277 PetscErrorCode ierr; 278 PetscMPIInt size; 279 PC_LU *dir; 280 281 PetscFunctionBegin; 282 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 283 284 ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 285 ((PC_Factor*)dir)->fact = PETSC_NULL; 286 ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 287 dir->inplace = PETSC_FALSE; 288 dir->nonzerosalongdiagonal = PETSC_FALSE; 289 290 ((PC_Factor*)dir)->info.fill = 5.0; 291 ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 292 ((PC_Factor*)dir)->info.shifttype = (PetscReal) MAT_SHIFT_NONE; 293 ((PC_Factor*)dir)->info.shiftamount = 0.0; 294 ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 295 ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 296 dir->col = 0; 297 dir->row = 0; 298 299 ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 300 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 301 if (size == 1) { 302 ierr = PetscStrallocpy(MATORDERINGND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 303 } else { 304 ierr = PetscStrallocpy(MATORDERINGNATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 305 } 306 dir->reusefill = PETSC_FALSE; 307 dir->reuseordering = PETSC_FALSE; 308 pc->data = (void*)dir; 309 310 pc->ops->reset = PCReset_LU; 311 pc->ops->destroy = PCDestroy_LU; 312 pc->ops->apply = PCApply_LU; 313 pc->ops->applytranspose = PCApplyTranspose_LU; 314 pc->ops->setup = PCSetUp_LU; 315 pc->ops->setfromoptions = PCSetFromOptions_LU; 316 pc->ops->view = PCView_LU; 317 pc->ops->applyrichardson = 0; 318 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 319 320 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 321 PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 322 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 323 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 324 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 325 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 326 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 327 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 328 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 329 PCFactorSetShiftType_Factor);CHKERRQ(ierr); 330 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 331 PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 332 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 333 PCFactorSetFill_Factor);CHKERRQ(ierr); 334 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 335 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 336 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 337 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 338 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 339 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 340 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 341 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 342 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 343 PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 344 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 345 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 346 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 347 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 348 PetscFunctionReturn(0); 349 } 350 EXTERN_C_END 351