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