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