1 2 #include "../src/ksp/pc/impls/mg/mgimpl.h" /*I "petscksp.h" I*/ 3 /*I "petscpcmg.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "PCMGDefaultResidual" 7 /*@C 8 PCMGDefaultResidual - Default routine to calculate the residual. 9 10 Collective on Mat and Vec 11 12 Input Parameters: 13 + mat - the matrix 14 . b - the right-hand-side 15 - x - the approximate solution 16 17 Output Parameter: 18 . r - location to store the residual 19 20 Level: advanced 21 22 .keywords: MG, default, multigrid, residual 23 24 .seealso: PCMGSetResidual() 25 @*/ 26 PetscErrorCode PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 27 { 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 ierr = MatMult(mat,x,r);CHKERRQ(ierr); 32 ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 /* ---------------------------------------------------------------------------*/ 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "PCMGGetCoarseSolve" 40 /*@ 41 PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 42 43 Not Collective 44 45 Input Parameter: 46 . pc - the multigrid context 47 48 Output Parameter: 49 . ksp - the coarse grid solver context 50 51 Level: advanced 52 53 .keywords: MG, multigrid, get, coarse grid 54 @*/ 55 PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 56 { 57 PC_MG *mg = (PC_MG*)pc->data; 58 PC_MG_Levels **mglevels = mg->levels; 59 60 PetscFunctionBegin; 61 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 62 *ksp = mglevels[0]->smoothd; 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCMGSetResidual" 68 /*@C 69 PCMGSetResidual - Sets the function to be used to calculate the residual 70 on the lth level. 71 72 Logically Collective on PC and Mat 73 74 Input Parameters: 75 + pc - the multigrid context 76 . l - the level (0 is coarsest) to supply 77 . residual - function used to form residual (usually PCMGDefaultResidual) 78 - mat - matrix associated with residual 79 80 Level: advanced 81 82 .keywords: MG, set, multigrid, residual, level 83 84 .seealso: PCMGDefaultResidual() 85 @*/ 86 PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 87 { 88 PC_MG *mg = (PC_MG*)pc->data; 89 PC_MG_Levels **mglevels = mg->levels; 90 91 PetscFunctionBegin; 92 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 93 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 94 mglevels[l]->residual = residual; 95 mglevels[l]->A = mat; 96 PetscFunctionReturn(0); 97 } 98 99 #undef __FUNCT__ 100 #define __FUNCT__ "PCMGSetInterpolation" 101 /*@ 102 PCMGSetInterpolation - Sets the function to be used to calculate the 103 interpolation from l-1 to the lth level 104 105 Logically Collective on PC and Mat 106 107 Input Parameters: 108 + pc - the multigrid context 109 . mat - the interpolation operator 110 - l - the level (0 is coarsest) to supply [do not supply 0] 111 112 Level: advanced 113 114 Notes: 115 Usually this is the same matrix used also to set the restriction 116 for the same level. 117 118 One can pass in the interpolation matrix or its transpose; PETSc figures 119 out from the matrix size which one it is. 120 121 .keywords: multigrid, set, interpolate, level 122 123 .seealso: PCMGSetRestriction() 124 @*/ 125 PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 126 { 127 PC_MG *mg = (PC_MG*)pc->data; 128 PC_MG_Levels **mglevels = mg->levels; 129 PetscErrorCode ierr; 130 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 133 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 134 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 135 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 136 if (mglevels[l]->interpolate) {ierr = MatDestroy(mglevels[l]->interpolate);CHKERRQ(ierr);} 137 mglevels[l]->interpolate = mat; 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "PCMGSetRestriction" 143 /*@ 144 PCMGSetRestriction - Sets the function to be used to restrict vector 145 from level l to l-1. 146 147 Logically Collective on PC and Mat 148 149 Input Parameters: 150 + pc - the multigrid context 151 . mat - the restriction matrix 152 - l - the level (0 is coarsest) to supply [Do not supply 0] 153 154 Level: advanced 155 156 Notes: 157 Usually this is the same matrix used also to set the interpolation 158 for the same level. 159 160 One can pass in the interpolation matrix or its transpose; PETSc figures 161 out from the matrix size which one it is. 162 163 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 164 is used. 165 166 .keywords: MG, set, multigrid, restriction, level 167 168 .seealso: PCMGSetInterpolation() 169 @*/ 170 PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 171 { 172 PetscErrorCode ierr; 173 PC_MG *mg = (PC_MG*)pc->data; 174 PC_MG_Levels **mglevels = mg->levels; 175 176 PetscFunctionBegin; 177 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 178 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 179 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 180 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 181 if (mglevels[l]->restrct) {ierr = MatDestroy(mglevels[l]->restrct);CHKERRQ(ierr);} 182 mglevels[l]->restrct = mat; 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "PCMGGetSmoother" 188 /*@ 189 PCMGGetSmoother - Gets the KSP context to be used as smoother for 190 both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 191 PCMGGetSmootherDown() to use different functions for pre- and 192 post-smoothing. 193 194 Not Collective, KSP returned is parallel if PC is 195 196 Input Parameters: 197 + pc - the multigrid context 198 - l - the level (0 is coarsest) to supply 199 200 Ouput Parameters: 201 . ksp - the smoother 202 203 Level: advanced 204 205 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 206 207 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 208 @*/ 209 PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 210 { 211 PC_MG *mg = (PC_MG*)pc->data; 212 PC_MG_Levels **mglevels = mg->levels; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 216 *ksp = mglevels[l]->smoothd; 217 PetscFunctionReturn(0); 218 } 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "PCMGGetSmootherUp" 222 /*@ 223 PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 224 coarse grid correction (post-smoother). 225 226 Not Collective, KSP returned is parallel if PC is 227 228 Input Parameters: 229 + pc - the multigrid context 230 - l - the level (0 is coarsest) to supply 231 232 Ouput Parameters: 233 . ksp - the smoother 234 235 Level: advanced 236 237 .keywords: MG, multigrid, get, smoother, up, post-smoother, level 238 239 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 240 @*/ 241 PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 242 { 243 PC_MG *mg = (PC_MG*)pc->data; 244 PC_MG_Levels **mglevels = mg->levels; 245 PetscErrorCode ierr; 246 const char *prefix; 247 MPI_Comm comm; 248 249 PetscFunctionBegin; 250 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 251 /* 252 This is called only if user wants a different pre-smoother from post. 253 Thus we check if a different one has already been allocated, 254 if not we allocate it. 255 */ 256 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 257 if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 258 ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 259 ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 260 ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 261 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 262 ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 263 ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 264 ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr); 265 } 266 if (ksp) *ksp = mglevels[l]->smoothu; 267 PetscFunctionReturn(0); 268 } 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "PCMGGetSmootherDown" 272 /*@ 273 PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 274 coarse grid correction (pre-smoother). 275 276 Not Collective, KSP returned is parallel if PC is 277 278 Input Parameters: 279 + pc - the multigrid context 280 - l - the level (0 is coarsest) to supply 281 282 Ouput Parameters: 283 . ksp - the smoother 284 285 Level: advanced 286 287 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 288 289 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 290 @*/ 291 PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 292 { 293 PetscErrorCode ierr; 294 PC_MG *mg = (PC_MG*)pc->data; 295 PC_MG_Levels **mglevels = mg->levels; 296 297 PetscFunctionBegin; 298 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 299 /* make sure smoother up and down are different */ 300 if (l) { 301 ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 302 } 303 *ksp = mglevels[l]->smoothd; 304 PetscFunctionReturn(0); 305 } 306 307 #undef __FUNCT__ 308 #define __FUNCT__ "PCMGSetCyclesOnLevel" 309 /*@ 310 PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 311 312 Logically Collective on PC 313 314 Input Parameters: 315 + pc - the multigrid context 316 . l - the level (0 is coarsest) this is to be used for 317 - n - the number of cycles 318 319 Level: advanced 320 321 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 322 323 .seealso: PCMGSetCycles() 324 @*/ 325 PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 326 { 327 PC_MG *mg = (PC_MG*)pc->data; 328 PC_MG_Levels **mglevels = mg->levels; 329 330 PetscFunctionBegin; 331 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 332 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 333 PetscValidLogicalCollectiveInt(pc,l,2); 334 PetscValidLogicalCollectiveInt(pc,c,3); 335 mglevels[l]->cycles = c; 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "PCMGSetRhs" 341 /*@ 342 PCMGSetRhs - Sets the vector space to be used to store the right-hand side 343 on a particular level. 344 345 Logically Collective on PC and Vec 346 347 Input Parameters: 348 + pc - the multigrid context 349 . l - the level (0 is coarsest) this is to be used for 350 - c - the space 351 352 Level: advanced 353 354 Notes: If this is not provided PETSc will automatically generate one. 355 356 You do not need to keep a reference to this vector if you do 357 not need it PCDestroy() will properly free it. 358 359 .keywords: MG, multigrid, set, right-hand-side, rhs, level 360 361 .seealso: PCMGSetX(), PCMGSetR() 362 @*/ 363 PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 364 { 365 PetscErrorCode ierr; 366 PC_MG *mg = (PC_MG*)pc->data; 367 PC_MG_Levels **mglevels = mg->levels; 368 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 371 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 372 if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 373 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 374 if (mglevels[l]->b) {ierr = VecDestroy(mglevels[l]->b);CHKERRQ(ierr);} 375 mglevels[l]->b = c; 376 PetscFunctionReturn(0); 377 } 378 379 #undef __FUNCT__ 380 #define __FUNCT__ "PCMGSetX" 381 /*@ 382 PCMGSetX - Sets the vector space to be used to store the solution on a 383 particular level. 384 385 Logically Collective on PC and Vec 386 387 Input Parameters: 388 + pc - the multigrid context 389 . l - the level (0 is coarsest) this is to be used for 390 - c - the space 391 392 Level: advanced 393 394 Notes: If this is not provided PETSc will automatically generate one. 395 396 You do not need to keep a reference to this vector if you do 397 not need it PCDestroy() will properly free it. 398 399 .keywords: MG, multigrid, set, solution, level 400 401 .seealso: PCMGSetRhs(), PCMGSetR() 402 @*/ 403 PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 404 { 405 PetscErrorCode ierr; 406 PC_MG *mg = (PC_MG*)pc->data; 407 PC_MG_Levels **mglevels = mg->levels; 408 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 411 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 412 if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 413 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 414 if (mglevels[l]->x) {ierr = VecDestroy(mglevels[l]->x);CHKERRQ(ierr);} 415 mglevels[l]->x = c; 416 PetscFunctionReturn(0); 417 } 418 419 #undef __FUNCT__ 420 #define __FUNCT__ "PCMGSetR" 421 /*@ 422 PCMGSetR - Sets the vector space to be used to store the residual on a 423 particular level. 424 425 Logically Collective on PC and Vec 426 427 Input Parameters: 428 + pc - the multigrid context 429 . l - the level (0 is coarsest) this is to be used for 430 - c - the space 431 432 Level: advanced 433 434 Notes: If this is not provided PETSc will automatically generate one. 435 436 You do not need to keep a reference to this vector if you do 437 not need it PCDestroy() will properly free it. 438 439 .keywords: MG, multigrid, set, residual, level 440 @*/ 441 PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 442 { 443 PetscErrorCode ierr; 444 PC_MG *mg = (PC_MG*)pc->data; 445 PC_MG_Levels **mglevels = mg->levels; 446 447 PetscFunctionBegin; 448 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 449 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 450 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 451 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 452 if (mglevels[l]->r) {ierr = VecDestroy(mglevels[l]->r);CHKERRQ(ierr);} 453 mglevels[l]->r = c; 454 PetscFunctionReturn(0); 455 } 456