1 #define PETSCKSP_DLL 2 3 /* 4 Defines a direct factorization preconditioner for any Mat implementation 5 Note: this need not be consided a preconditioner since it supplies 6 a direct solver. 7 */ 8 9 #include "../src/ksp/pc/impls/factor/lu/lu.h" /*I "petscpc.h" I*/ 10 11 EXTERN_C_BEGIN 12 #undef __FUNCT__ 13 #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_LU" 14 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 15 { 16 PC_LU *lu = (PC_LU*)pc->data; 17 18 PetscFunctionBegin; 19 lu->nonzerosalongdiagonal = PETSC_TRUE; 20 if (z == PETSC_DECIDE) { 21 lu->nonzerosalongdiagonaltol = 1.e-10; 22 } else { 23 lu->nonzerosalongdiagonaltol = z; 24 } 25 PetscFunctionReturn(0); 26 } 27 EXTERN_C_END 28 29 EXTERN_C_BEGIN 30 #undef __FUNCT__ 31 #define __FUNCT__ "PCFactorSetReuseOrdering_LU" 32 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_LU(PC pc,PetscTruth flag) 33 { 34 PC_LU *lu = (PC_LU*)pc->data; 35 36 PetscFunctionBegin; 37 lu->reuseordering = flag; 38 PetscFunctionReturn(0); 39 } 40 EXTERN_C_END 41 42 EXTERN_C_BEGIN 43 #undef __FUNCT__ 44 #define __FUNCT__ "PCFactorSetReuseFill_LU" 45 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_LU(PC pc,PetscTruth flag) 46 { 47 PC_LU *lu = (PC_LU*)pc->data; 48 49 PetscFunctionBegin; 50 lu->reusefill = flag; 51 PetscFunctionReturn(0); 52 } 53 EXTERN_C_END 54 55 #undef __FUNCT__ 56 #define __FUNCT__ "PCSetFromOptions_LU" 57 static PetscErrorCode PCSetFromOptions_LU(PC pc) 58 { 59 PC_LU *lu = (PC_LU*)pc->data; 60 PetscErrorCode ierr; 61 PetscTruth flg = PETSC_FALSE; 62 PetscReal tol; 63 64 PetscFunctionBegin; 65 ierr = PetscOptionsHead("LU options");CHKERRQ(ierr); 66 ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 67 68 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 69 if (flg) { 70 tol = PETSC_DECIDE; 71 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 72 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 73 } 74 ierr = PetscOptionsTail();CHKERRQ(ierr); 75 PetscFunctionReturn(0); 76 } 77 78 #undef __FUNCT__ 79 #define __FUNCT__ "PCView_LU" 80 static PetscErrorCode PCView_LU(PC pc,PetscViewer viewer) 81 { 82 PC_LU *lu = (PC_LU*)pc->data; 83 PetscErrorCode ierr; 84 PetscTruth iascii,isstring; 85 86 PetscFunctionBegin; 87 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 88 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 89 if (iascii) { 90 91 if (lu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," LU: in-place factorization\n");CHKERRQ(ierr);} 92 else {ierr = PetscViewerASCIIPrintf(viewer," LU: out-of-place factorization\n");CHKERRQ(ierr);} 93 ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)lu)->ordering);CHKERRQ(ierr); 94 ierr = PetscViewerASCIIPrintf(viewer," LU: tolerance for zero pivot %G\n",((PC_Factor*)lu)->info.zeropivot);CHKERRQ(ierr); 95 if (((PC_Factor*)lu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," LU: using Manteuffel shift\n");CHKERRQ(ierr);} 96 if (((PC_Factor*)lu)->fact) { 97 ierr = PetscViewerASCIIPrintf(viewer," LU: factor fill ratio needed %G\n",lu->actualfill);CHKERRQ(ierr); 98 ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 99 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 100 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 101 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 102 ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 103 ierr = MatView(((PC_Factor*)lu)->fact,viewer);CHKERRQ(ierr); 104 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 105 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 106 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 107 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 108 } 109 if (lu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 110 if (lu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 111 } else if (isstring) { 112 ierr = PetscViewerStringSPrintf(viewer," order=%s",((PC_Factor*)lu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 113 } else { 114 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCLU",((PetscObject)viewer)->type_name); 115 } 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "PCSetUp_LU" 121 static PetscErrorCode PCSetUp_LU(PC pc) 122 { 123 PetscErrorCode ierr; 124 PC_LU *dir = (PC_LU*)pc->data; 125 126 PetscFunctionBegin; 127 if (dir->reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->actualfill; 128 129 if (dir->inplace) { 130 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 131 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 132 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 133 if (dir->row) { 134 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 135 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 136 } 137 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 138 ((PC_Factor*)dir)->fact = pc->pmat; 139 } else { 140 MatInfo info; 141 if (!pc->setupcalled) { 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(pc,dir->row);CHKERRQ(ierr); 148 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 149 } 150 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 151 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 152 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 153 dir->actualfill = info.fill_ratio_needed; 154 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 155 } else if (pc->flag != SAME_NONZERO_PATTERN) { 156 if (!dir->reuseordering) { 157 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 158 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 159 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 160 if (dir->nonzerosalongdiagonal) { 161 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 162 } 163 if (dir->row) { 164 ierr = PetscLogObjectParent(pc,dir->row);CHKERRQ(ierr); 165 ierr = PetscLogObjectParent(pc,dir->col);CHKERRQ(ierr); 166 } 167 } 168 ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 169 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 170 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 171 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 172 dir->actualfill = info.fill_ratio_needed; 173 ierr = PetscLogObjectParent(pc,((PC_Factor*)dir)->fact);CHKERRQ(ierr); 174 } 175 ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 176 } 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 if (!dir->inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 189 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(dir->row);CHKERRQ(ierr);} 190 if (dir->col) {ierr = ISDestroy(dir->col);CHKERRQ(ierr);} 191 ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 192 ierr = PetscStrfree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 193 ierr = PetscFree(dir);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "PCApply_LU" 199 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 200 { 201 PC_LU *dir = (PC_LU*)pc->data; 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 if (dir->inplace) {ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr);} 206 else {ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "PCApplyTranspose_LU" 212 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 213 { 214 PC_LU *dir = (PC_LU*)pc->data; 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 if (dir->inplace) {ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr);} 219 else {ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr);} 220 PetscFunctionReturn(0); 221 } 222 223 /* -----------------------------------------------------------------------------------*/ 224 225 EXTERN_C_BEGIN 226 #undef __FUNCT__ 227 #define __FUNCT__ "PCFactorSetUseInPlace_LU" 228 PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_LU(PC pc) 229 { 230 PC_LU *dir = (PC_LU*)pc->data; 231 232 PetscFunctionBegin; 233 dir->inplace = PETSC_TRUE; 234 PetscFunctionReturn(0); 235 } 236 EXTERN_C_END 237 238 /* ------------------------------------------------------------------------ */ 239 240 /*MC 241 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 242 243 Options Database Keys: 244 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 245 . -pc_factor_mat_solver_package - Actives PCFactorSetMatSolverPackage() to choose the direct solver, like spooles 246 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 247 . -pc_factor_fill <fill> - Sets fill amount 248 . -pc_factor_in_place - Activates in-place factorization 249 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 250 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 251 stability of factorization. 252 . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 253 . -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 254 is optional with PETSC_TRUE being the default 255 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 256 diagonal. 257 258 Notes: Not all options work for all matrix formats 259 Run with -help to see additional options for particular matrix formats or factorization 260 algorithms 261 262 Level: beginner 263 264 Concepts: LU factorization, direct solver 265 266 Notes: Usually this will compute an "exact" solution in one iteration and does 267 not need a Krylov method (i.e. you can use -ksp_type preonly, or 268 KSPSetType(ksp,KSPPREONLY) for the Krylov method 269 270 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 271 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 272 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 273 PCFactorSetPivotingInBlocks(),PCFactorSetShiftNonzero(),PCFactorSetShiftPd(),PCFactorSetShiftInBlocks() 274 PCFactorReorderForNonzeroDiagonal() 275 M*/ 276 277 EXTERN_C_BEGIN 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCCreate_LU" 280 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_LU(PC pc) 281 { 282 PetscErrorCode ierr; 283 PetscMPIInt size; 284 PC_LU *dir; 285 286 PetscFunctionBegin; 287 ierr = PetscNewLog(pc,PC_LU,&dir);CHKERRQ(ierr); 288 289 ierr = MatFactorInfoInitialize(&((PC_Factor*)dir)->info);CHKERRQ(ierr); 290 ((PC_Factor*)dir)->fact = 0; 291 dir->inplace = PETSC_FALSE; 292 dir->nonzerosalongdiagonal = PETSC_FALSE; 293 294 ((PC_Factor*)dir)->info.fill = 5.0; 295 ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 296 ((PC_Factor*)dir)->info.shiftnz = 0.0; 297 ((PC_Factor*)dir)->info.zeropivot = 1.e-12; 298 ((PC_Factor*)dir)->info.pivotinblocks = 1.0; 299 ((PC_Factor*)dir)->info.shiftpd = 0.0; /* false */ 300 dir->col = 0; 301 dir->row = 0; 302 303 ierr = PetscStrallocpy(MAT_SOLVER_PETSC,&((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 304 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 305 if (size == 1) { 306 ierr = PetscStrallocpy(MATORDERING_ND,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 307 } else { 308 ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 309 } 310 dir->reusefill = PETSC_FALSE; 311 dir->reuseordering = PETSC_FALSE; 312 pc->data = (void*)dir; 313 314 pc->ops->destroy = PCDestroy_LU; 315 pc->ops->apply = PCApply_LU; 316 pc->ops->applytranspose = PCApplyTranspose_LU; 317 pc->ops->setup = PCSetUp_LU; 318 pc->ops->setfromoptions = PCSetFromOptions_LU; 319 pc->ops->view = PCView_LU; 320 pc->ops->applyrichardson = 0; 321 pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 322 323 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorGetMatSolverPackage_C","PCFactorGetMatSolverPackage_Factor", 324 PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 325 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatSolverPackage_C","PCFactorSetMatSolverPackage_Factor", 326 PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 327 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_Factor", 328 PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 329 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_Factor", 330 PCFactorSetShiftNonzero_Factor);CHKERRQ(ierr); 331 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_Factor", 332 PCFactorSetShiftPd_Factor);CHKERRQ(ierr); 333 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_Factor", 334 PCFactorSetShiftInBlocks_Factor);CHKERRQ(ierr); 335 336 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_Factor", 337 PCFactorSetFill_Factor);CHKERRQ(ierr); 338 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_LU", 339 PCFactorSetUseInPlace_LU);CHKERRQ(ierr); 340 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_Factor", 341 PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 342 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_LU", 343 PCFactorSetReuseOrdering_LU);CHKERRQ(ierr); 344 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_LU", 345 PCFactorSetReuseFill_LU);CHKERRQ(ierr); 346 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetColumnPivot_C","PCFactorSetColumnPivot_Factor", 347 PCFactorSetColumnPivot_Factor);CHKERRQ(ierr); 348 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_Factor", 349 PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 350 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_LU", 351 PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 EXTERN_C_END 355