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 = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 64 if (dir->row) { 65 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 66 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 67 } 68 ierr = MatLUFactor(pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 69 ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 70 if (err) { /* Factor() fails */ 71 pc->failedreason = (PCFailedReason)err; 72 PetscFunctionReturn(0); 73 } 74 } 75 ((PC_Factor*)dir)->fact = pc->pmat; 76 } else { 77 MatInfo info; 78 79 if (!pc->setupcalled) { 80 PetscBool useordering; 81 if (!((PC_Factor*)dir)->fact) { 82 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 83 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 84 } 85 ierr = MatFactorGetUseOrdering(((PC_Factor*)dir)->fact,&useordering);CHKERRQ(ierr); 86 if (useordering) { 87 ierr = MatGetOrdering(pc->pmat,((PC_Factor*)dir)->ordering,&dir->row,&dir->col);CHKERRQ(ierr); 88 if (dir->nonzerosalongdiagonal) { 89 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 90 } 91 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 92 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 93 } 94 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 95 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 96 dir->hdr.actualfill = info.fill_ratio_needed; 97 } else if (pc->flag != SAME_NONZERO_PATTERN) { 98 PetscBool useordering; 99 if (!dir->hdr.reuseordering) { 100 ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 101 ierr = MatGetFactor(pc->pmat,((PC_Factor*)dir)->solvertype,MAT_FACTOR_LU,&((PC_Factor*)dir)->fact);CHKERRQ(ierr); 102 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)dir)->fact);CHKERRQ(ierr); 103 ierr = MatFactorGetUseOrdering(((PC_Factor*)dir)->fact,&useordering);CHKERRQ(ierr); 104 if (useordering) { 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->nonzerosalongdiagonal) { 109 ierr = MatReorderForNonzeroDiagonal(pc->pmat,dir->nonzerosalongdiagonaltol,dir->row,dir->col);CHKERRQ(ierr); 110 } 111 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->row);CHKERRQ(ierr); 112 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)dir->col);CHKERRQ(ierr); 113 } 114 } 115 ierr = MatLUFactorSymbolic(((PC_Factor*)dir)->fact,pc->pmat,dir->row,dir->col,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 116 ierr = MatGetInfo(((PC_Factor*)dir)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 117 dir->hdr.actualfill = info.fill_ratio_needed; 118 } else { 119 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 120 if (err == MAT_FACTOR_NUMERIC_ZEROPIVOT) { 121 ierr = MatFactorClearError(((PC_Factor*)dir)->fact);CHKERRQ(ierr); 122 pc->failedreason = PC_NOERROR; 123 } 124 } 125 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 126 if (err) { /* FactorSymbolic() fails */ 127 pc->failedreason = (PCFailedReason)err; 128 PetscFunctionReturn(0); 129 } 130 131 ierr = MatLUFactorNumeric(((PC_Factor*)dir)->fact,pc->pmat,&((PC_Factor*)dir)->info);CHKERRQ(ierr); 132 ierr = MatFactorGetError(((PC_Factor*)dir)->fact,&err);CHKERRQ(ierr); 133 if (err) { /* FactorNumeric() fails */ 134 pc->failedreason = (PCFailedReason)err; 135 } 136 137 } 138 139 ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 140 if (!stype) { 141 MatSolverType solverpackage; 142 ierr = MatFactorGetSolverType(((PC_Factor*)dir)->fact,&solverpackage);CHKERRQ(ierr); 143 ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 144 } 145 PetscFunctionReturn(0); 146 } 147 148 static PetscErrorCode PCReset_LU(PC pc) 149 { 150 PC_LU *dir = (PC_LU*)pc->data; 151 PetscErrorCode ierr; 152 153 PetscFunctionBegin; 154 if (!dir->hdr.inplace && ((PC_Factor*)dir)->fact) {ierr = MatDestroy(&((PC_Factor*)dir)->fact);CHKERRQ(ierr);} 155 if (dir->row && dir->col && dir->row != dir->col) {ierr = ISDestroy(&dir->row);CHKERRQ(ierr);} 156 ierr = ISDestroy(&dir->col);CHKERRQ(ierr); 157 PetscFunctionReturn(0); 158 } 159 160 static PetscErrorCode PCDestroy_LU(PC pc) 161 { 162 PC_LU *dir = (PC_LU*)pc->data; 163 PetscErrorCode ierr; 164 165 PetscFunctionBegin; 166 ierr = PCReset_LU(pc);CHKERRQ(ierr); 167 ierr = PetscFree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 168 ierr = PetscFree(((PC_Factor*)dir)->solvertype);CHKERRQ(ierr); 169 ierr = PetscFree(pc->data);CHKERRQ(ierr); 170 PetscFunctionReturn(0); 171 } 172 173 static PetscErrorCode PCApply_LU(PC pc,Vec x,Vec y) 174 { 175 PC_LU *dir = (PC_LU*)pc->data; 176 PetscErrorCode ierr; 177 178 PetscFunctionBegin; 179 if (dir->hdr.inplace) { 180 ierr = MatSolve(pc->pmat,x,y);CHKERRQ(ierr); 181 } else { 182 ierr = MatSolve(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 183 } 184 PetscFunctionReturn(0); 185 } 186 187 static PetscErrorCode PCMatApply_LU(PC pc,Mat X,Mat Y) 188 { 189 PC_LU *dir = (PC_LU*)pc->data; 190 PetscErrorCode ierr; 191 192 PetscFunctionBegin; 193 if (dir->hdr.inplace) { 194 ierr = MatMatSolve(pc->pmat,X,Y);CHKERRQ(ierr); 195 } else { 196 ierr = MatMatSolve(((PC_Factor*)dir)->fact,X,Y);CHKERRQ(ierr); 197 } 198 PetscFunctionReturn(0); 199 } 200 201 static PetscErrorCode PCApplyTranspose_LU(PC pc,Vec x,Vec y) 202 { 203 PC_LU *dir = (PC_LU*)pc->data; 204 PetscErrorCode ierr; 205 206 PetscFunctionBegin; 207 if (dir->hdr.inplace) { 208 ierr = MatSolveTranspose(pc->pmat,x,y);CHKERRQ(ierr); 209 } else { 210 ierr = MatSolveTranspose(((PC_Factor*)dir)->fact,x,y);CHKERRQ(ierr); 211 } 212 PetscFunctionReturn(0); 213 } 214 215 /* -----------------------------------------------------------------------------------*/ 216 217 /*MC 218 PCLU - Uses a direct solver, based on LU factorization, as a preconditioner 219 220 Options Database Keys: 221 + -pc_factor_reuse_ordering - Activate PCFactorSetReuseOrdering() 222 . -pc_factor_mat_solver_type - Actives PCFactorSetMatSolverType() to choose the direct solver, like superlu 223 . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 224 . -pc_factor_fill <fill> - Sets fill amount 225 . -pc_factor_in_place - Activates in-place factorization 226 . -pc_factor_mat_ordering_type <nd,rcm,...> - Sets ordering routine 227 . -pc_factor_pivot_in_blocks <true,false> - allow pivoting within the small blocks during factorization (may increase 228 stability of factorization. 229 . -pc_factor_shift_type <shifttype> - Sets shift type or PETSC_DECIDE for the default; use '-help' for a list of available types 230 . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 231 - -pc_factor_nonzeros_along_diagonal - permutes the rows and columns to try to put nonzero value along the 232 diagonal. 233 234 Notes: 235 Not all options work for all matrix formats 236 Run with -help to see additional options for particular matrix formats or factorization 237 algorithms 238 239 Level: beginner 240 241 Notes: 242 Usually this will compute an "exact" solution in one iteration and does 243 not need a Krylov method (i.e. you can use -ksp_type preonly, or 244 KSPSetType(ksp,KSPPREONLY) for the Krylov method 245 246 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 247 PCILU, PCCHOLESKY, PCICC, PCFactorSetReuseOrdering(), PCFactorSetReuseFill(), PCFactorGetMatrix(), 248 PCFactorSetFill(), PCFactorSetUseInPlace(), PCFactorSetMatOrderingType(), PCFactorSetColumnPivot(), 249 PCFactorSetPivotingInBlocks(),PCFactorSetShiftType(),PCFactorSetShiftAmount() 250 PCFactorReorderForNonzeroDiagonal() 251 M*/ 252 253 PETSC_EXTERN PetscErrorCode PCCreate_LU(PC pc) 254 { 255 PetscErrorCode ierr; 256 PetscMPIInt size; 257 PC_LU *dir; 258 259 PetscFunctionBegin; 260 ierr = PetscNewLog(pc,&dir);CHKERRQ(ierr); 261 pc->data = (void*)dir; 262 ierr = PCFactorInitialize(pc);CHKERRQ(ierr); 263 ((PC_Factor*)dir)->factortype = MAT_FACTOR_LU; 264 dir->nonzerosalongdiagonal = PETSC_FALSE; 265 266 ((PC_Factor*)dir)->info.fill = 5.0; 267 ((PC_Factor*)dir)->info.dtcol = 1.e-6; /* default to pivoting; this is only thing PETSc LU supports */ 268 ((PC_Factor*)dir)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 269 dir->col = NULL; 270 dir->row = NULL; 271 272 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 273 if (size == 1) { 274 ierr = PetscStrallocpy(MATORDERINGND,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 275 } else { 276 ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 277 } 278 279 pc->ops->reset = PCReset_LU; 280 pc->ops->destroy = PCDestroy_LU; 281 pc->ops->apply = PCApply_LU; 282 pc->ops->matapply = PCMatApply_LU; 283 pc->ops->applytranspose = PCApplyTranspose_LU; 284 pc->ops->setup = PCSetUp_LU; 285 pc->ops->setfromoptions = PCSetFromOptions_LU; 286 pc->ops->view = PCView_Factor; 287 pc->ops->applyrichardson = NULL; 288 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_LU);CHKERRQ(ierr); 289 PetscFunctionReturn(0); 290 } 291