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