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 PetscErrorCode PCFactorReorderForNonzeroDiagonal_LU(PC pc,PetscReal z) 12 { 13 PC_LU *lu = (PC_LU*)pc->data; 14 15 PetscFunctionBegin; 16 lu->nonzerosalongdiagonal = PETSC_TRUE; 17 if (z == PETSC_DECIDE) lu->nonzerosalongdiagonaltol = 1.e-10; 18 else lu->nonzerosalongdiagonaltol = z; 19 PetscFunctionReturn(0); 20 } 21 22 static PetscErrorCode PCSetFromOptions_LU(PetscOptionItems *PetscOptionsObject,PC pc) 23 { 24 PC_LU *lu = (PC_LU*)pc->data; 25 PetscErrorCode ierr; 26 PetscBool flg = PETSC_FALSE; 27 PetscReal tol; 28 29 PetscFunctionBegin; 30 ierr = PetscOptionsHead(PetscOptionsObject,"LU options");CHKERRQ(ierr); 31 ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 32 33 ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 34 if (flg) { 35 tol = PETSC_DECIDE; 36 ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",lu->nonzerosalongdiagonaltol,&tol,NULL);CHKERRQ(ierr); 37 ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 38 } 39 ierr = PetscOptionsTail();CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 static PetscErrorCode PCSetUp_LU(PC pc) 44 { 45 PetscErrorCode ierr; 46 PC_LU *dir = (PC_LU*)pc->data; 47 MatSolverType stype; 48 MatFactorError err; 49 50 PetscFunctionBegin; 51 pc->failedreason = PC_NOERROR; 52 if (dir->hdr.reusefill && pc->setupcalled) ((PC_Factor*)dir)->info.fill = dir->hdr.actualfill; 53 54 ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 55 if (dir->hdr.inplace) { 56 MatFactorType ftype; 57 58 ierr = MatGetFactorType(pc->pmat, &ftype);CHKERRQ(ierr); 59 if (ftype == MAT_FACTOR_NONE) { 60 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 61 ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 62 /* This should only get the ordering if needed, but since MatGetFactor() is not called we can't know if it is needed */ 63 ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 64 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 65 if (dir->row) { 66 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 67 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 68 } 69 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 70 ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 71 if (err) { /* Factor() fails */ 72 pc->failedreason = (PCFailedReason)err; 73 PetscFunctionReturn(0); 74 } 75 } 76 ((PC_Factor*)dir)->fact = pc->pmat; 77 } else { 78 MatInfo info; 79 80 if (!pc->setupcalled) { 81 PetscBool canuseordering; 82 if (!((PC_Factor*)dir)->fact) { 83 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 84 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 85 } 86 ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr); 87 if (canuseordering) { 88 ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 89 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 90 if (dir->nonzerosalongdiagonal) { 91 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 92 } 93 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 94 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 95 } 96 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 97 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 98 dir->hdr.actualfill = info.fill_ratio_needed; 99 } else if (pc->flag != SAME_NONZERO_PATTERN) { 100 PetscBool canuseordering; 101 if (!dir->hdr.reuseordering) { 102 ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 103 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 104 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 105 ierr = MatFactorGetCanUseOrdering(((PC_Factor*)dir)->fact,&canuseordering);CHKERRQ(ierr); 106 if (canuseordering) { 107 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 108 ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 109 ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 110 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 111 if (dir->nonzerosalongdiagonal) { 112 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 113 } 114 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 115 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 116 } 117 } 118 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 119 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 120 dir->hdr.actualfill = info.fill_ratio_needed; 121 } else { 122 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 123 if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 124 ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 125 pc->failedreason = PC_NOERROR; 126 } 127 } 128 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 129 if (err) { /* FactorSymbolic() fails */ 130 pc->failedreason = (PCFailedReason)err; 131 PetscFunctionReturn(0); 132 } 133 134 ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 135 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 136 if (err) { /* FactorNumeric() fails */ 137 pc->failedreason = (PCFailedReason)err; 138 } 139 140 } 141 142 ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 143 if (!stype) { 144 MatSolverType solverpackage; 145 ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 146 ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 147 } 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode PCReset_LU(PC pc) 152 { 153 PC_LU *dir = (PC_LU*)pc->data; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 158 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 159 ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 160 PetscFunctionReturn(0); 161 } 162 163 static PetscErrorCode PCDestroy_LU(PC pc) 164 { 165 PC_LU *dir = (PC_LU*)pc->data; 166 PetscErrorCode ierr; 167 168 PetscFunctionBegin; 169 ierr = PCReset_LU(pc);CHKERRQ(ierr); 170 ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 171 ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 172 ierr = PetscFree(pc->data);CHKERRQ(ierr); 173 PetscFunctionReturn(0); 174 } 175 176 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 177 { 178 PC_LU *dir = (PC_LU*)pc->data; 179 PetscErrorCode ierr; 180 181 PetscFunctionBegin; 182 if (dir->hdr.inplace) { 183 ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 184 } else { 185 ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 186 } 187 PetscFunctionReturn(0); 188 } 189 190 static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y) 191 { 192 PC_LU *dir = (PC_LU*)pc->data; 193 PetscErrorCode ierr; 194 195 PetscFunctionBegin; 196 if (dir->hdr.inplace) { 197 ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr); 198 } else { 199 ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr); 200 } 201 PetscFunctionReturn(0); 202 } 203 204 static PetscErrorCode PCApplyTranspose_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->hdr.inplace) { 211 ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 212 } else { 213 ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 214 } 215 PetscFunctionReturn(0); 216 } 217 218 /* -----------------------------------------------------------------------------------*/ 219 220 /*MC 221 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 222 223 Options Database Keys: 224 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 225 . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 226 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 227 . -pc_factor_fill <fill> - Sets fill amount 228 . -pc_factor_in_place - Activates in-place factorization 229 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 230 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 231 stability of factorization. 232 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 233 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 234 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 235 diagonal. 236 237 Notes: 238 Not all options work for all matrix formats 239 Run with -help to see additional options for particular matrix formats or factorization 240 algorithms 241 242 Level: beginner 243 244 Notes: 245 Usually this will compute an "exact" solution in one iteration and does 246 not need a Krylov method (i.e. you can use -ksp_type preonly, or 247 KSPSetType(ksp,KSPPREONLY) for the Krylov method 248 249 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 250 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 251 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 252 PCFactorSetPivotInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 253 PCFactorReorderForNonzeroDiagonal() 254 M*/ 255 256 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 257 { 258 PetscErrorCode ierr; 259 PC_LU *dir; 260 261 PetscFunctionBegin; 262 ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 263 pc->data = (void*)dir; 264 ierr = PCFactorInitialize(pc,MAT_FACTOR_LU);CHKERRQ(ierr); 265 dir->nonzerosalongdiagonal = PETSC_FALSE; 266 267 ((PC_Factor*)dir)->info.fill = 5.0; 268 ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 269 ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 270 dir->col = NULL; 271 dir->row = NULL; 272 273 pc->ops->reset = PCReset_LU; 274 pc->ops->destroy = PCDestroy_LU; 275 pc->ops->apply = PCApply_LU; 276 pc->ops->matapply = PCMatApply_LU; 277 pc->ops->applytranspose = PCApplyTranspose_LU; 278 pc->ops->setup = PCSetUp_LU; 279 pc->ops->setfromoptions = PCSetFromOptions_LU; 280 pc->ops->view = PCView_Factor; 281 pc->ops->applyrichardson = NULL; 282 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 283 PetscFunctionReturn(0); 284 } 285