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