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 = PetscObjectTypeCompare((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 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 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 if ((!pc->setupcalled) || (pc->setupcalled && !pc->reuse)) { 160 ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 161 } 162 } 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "PCReset_LU" 168 static PetscErrorCode PCReset_LU(PC pc) 169 { 170 PC_LU *dir = (PC_LU*)pc->data; 171 PetscErrorCode ierr; 172 173 PetscFunctionBegin; 174 if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 175 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 176 ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "PCDestroy_LU" 182 static PetscErrorCode PCDestroy_LU(PC pc) 183 { 184 PC_LU *dir = (PC_LU*)pc->data; 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 ierr = PCReset_LU(pc);CHKERRQ(ierr); 189 ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 190 ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 191 ierr = PetscFree(pc->data);CHKERRQ(ierr); 192 PetscFunctionReturn(0); 193 } 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "PCApply_LU" 197 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 198 { 199 PC_LU *dir = (PC_LU*)pc->data; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 204 else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 205 PetscFunctionReturn(0); 206 } 207 208 #undef __FUNCT__ 209 #define __FUNCT__ "PCApplyTranspose_LU" 210 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 211 { 212 PC_LU *dir = (PC_LU*)pc->data; 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 217 else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 218 PetscFunctionReturn(0); 219 } 220 221 /* -----------------------------------------------------------------------------------*/ 222 223 EXTERN_C_BEGIN 224 #undef __FUNCT__ 225 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 226 PetscErrorCode PCFactorSetUseInPlace_LU(PC pc) 227 { 228 PC_LU *dir = (PC_LU*)pc->data; 229 230 PetscFunctionBegin; 231 dir->inplace = PETSC_TRUE; 232 PetscFunctionReturn(0); 233 } 234 EXTERN_C_END 235 236 /* ------------------------------------------------------------------------ */ 237 238 /*MC 239 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 240 241 Options Database Keys: 242 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 243 . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like superlu 244 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 245 . -pc_factor_fill <fill> - Sets fill amount 246 . -pc_factor_in_place - Activates in-place factorization 247 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 248 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 249 stability of factorization. 250 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 251 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 252 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 253 diagonal. 254 255 Notes: Not all options work for all matrix formats 256 Run with -help to see additional options for particular matrix formats or factorization 257 algorithms 258 259 Level: beginner 260 261 Concepts: LU factorization, direct solver 262 263 Notes: Usually this will compute an "exact" solution in one iteration and does 264 not need a Krylov method (i.e. you can use -ksp_type preonly, or 265 KSPSetType(ksp,KSPPREONLY) for the Krylov method 266 267 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 268 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 269 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 270 PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 271 PCFactorReorderForNonzeroDiagonal() 272 M*/ 273 274 EXTERN_C_BEGIN 275 #undef __FUNCT__ 276 #define __FUNCT__ "PCCreate_LU" 277 PetscErrorCode PCCreate_LU(PC pc) 278 { 279 PetscErrorCode ierr; 280 PetscMPIInt size; 281 PC_LU *dir; 282 283 PetscFunctionBegin; 284 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 285 286 ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 287 ((PC_Factor*)dir)->fact = PETSC_NULL; 288 ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 289 dir->inplace = PETSC_FALSE; 290 dir->nonzerosalongdiagonal = PETSC_FALSE; 291 292 ((PC_Factor*)dir)->info.fill = 5.0; 293 ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 294 ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 295 ((PC_Factor*)dir)->info.shiftamount = 0.0; 296 ((PC_Factor*)dir)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 297 ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 298 dir->col = 0; 299 dir->row = 0; 300 301 ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); /* default SolverPackage */ 302 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 303 if (size == 1) { 304 ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 305 } else { 306 ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 307 } 308 dir->reusefill = PETSC_FALSE; 309 dir->reuseordering = PETSC_FALSE; 310 pc->data = (void*)dir; 311 312 pc->ops->reset = PCReset_LU; 313 pc->ops->destroy = PCDestroy_LU; 314 pc->ops->apply = PCApply_LU; 315 pc->ops->applytranspose = PCApplyTranspose_LU; 316 pc->ops->setup = PCSetUp_LU; 317 pc->ops->setfromoptions = PCSetFromOptions_LU; 318 pc->ops->view = PCView_LU; 319 pc->ops->applyrichardson = 0; 320 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 321 322 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C","PCFactorSetUpMatSolverPackage_Factor", 323 PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 324 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 325 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 326 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 327 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 328 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 329 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 330 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftType_C","PCFactorSetShiftType_Factor", 331 PCFactorSetShiftType_Factor);CHKERRQ(ierr); 332 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftAmount_C","PCFactorSetShiftAmount_Factor", 333 PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 334 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 335 PCFactorSetFill_Factor);CHKERRQ(ierr); 336 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 337 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 338 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 339 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 340 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 341 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 342 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 343 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 344 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 345 PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 346 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 347 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 348 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 349 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 350 PetscFunctionReturn(0); 351 } 352 EXTERN_C_END 353