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