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