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