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