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