1 2 #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 3 4 /* ---------------------------------------------------------------------------*/ 5 /*@C 6 PCMGResidualDefault - Default routine to calculate the residual. 7 8 Collective on Mat and Vec 9 10 Input Parameters: 11 + mat - the matrix 12 . b - the right-hand-side 13 - x - the approximate solution 14 15 Output Parameter: 16 . r - location to store the residual 17 18 Level: developer 19 20 .keywords: MG, default, multigrid, residual 21 22 .seealso: PCMGSetResidual() 23 @*/ 24 PetscErrorCode PCMGResidualDefault(Mat mat,Vec b,Vec x,Vec r) 25 { 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 ierr = MatResidual(mat,b,x,r);CHKERRQ(ierr); 30 PetscFunctionReturn(0); 31 } 32 33 /*@ 34 PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 35 36 Not Collective 37 38 Input Parameter: 39 . pc - the multigrid context 40 41 Output Parameter: 42 . ksp - the coarse grid solver context 43 44 Level: advanced 45 46 .keywords: MG, multigrid, get, coarse grid 47 48 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetSmoother() 49 @*/ 50 PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 51 { 52 PC_MG *mg = (PC_MG*)pc->data; 53 PC_MG_Levels **mglevels = mg->levels; 54 55 PetscFunctionBegin; 56 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 57 *ksp = mglevels[0]->smoothd; 58 PetscFunctionReturn(0); 59 } 60 61 /*@C 62 PCMGSetResidual - Sets the function to be used to calculate the residual 63 on the lth level. 64 65 Logically Collective on PC and Mat 66 67 Input Parameters: 68 + pc - the multigrid context 69 . l - the level (0 is coarsest) to supply 70 . residual - function used to form residual, if none is provided the previously provide one is used, if no 71 previous one were provided then a default is used 72 - mat - matrix associated with residual 73 74 Level: advanced 75 76 .keywords: MG, set, multigrid, residual, level 77 78 .seealso: PCMGResidualDefault() 79 @*/ 80 PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 81 { 82 PC_MG *mg = (PC_MG*)pc->data; 83 PC_MG_Levels **mglevels = mg->levels; 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 89 if (residual) mglevels[l]->residual = residual; 90 if (!mglevels[l]->residual) mglevels[l]->residual = PCMGResidualDefault; 91 if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 92 ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 93 mglevels[l]->A = mat; 94 PetscFunctionReturn(0); 95 } 96 97 /*@ 98 PCMGSetInterpolation - Sets the function to be used to calculate the 99 interpolation from l-1 to the lth level 100 101 Logically Collective on PC and Mat 102 103 Input Parameters: 104 + pc - the multigrid context 105 . mat - the interpolation operator 106 - l - the level (0 is coarsest) to supply [do not supply 0] 107 108 Level: advanced 109 110 Notes: 111 Usually this is the same matrix used also to set the restriction 112 for the same level. 113 114 One can pass in the interpolation matrix or its transpose; PETSc figures 115 out from the matrix size which one it is. 116 117 .keywords: multigrid, set, interpolate, level 118 119 .seealso: PCMGSetRestriction() 120 @*/ 121 PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 122 { 123 PC_MG *mg = (PC_MG*)pc->data; 124 PC_MG_Levels **mglevels = mg->levels; 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 129 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 130 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 131 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 132 ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 133 134 mglevels[l]->interpolate = mat; 135 PetscFunctionReturn(0); 136 } 137 138 /*@ 139 PCMGGetInterpolation - Gets the function to be used to calculate the 140 interpolation from l-1 to the lth level 141 142 Logically Collective on PC 143 144 Input Parameters: 145 + pc - the multigrid context 146 - l - the level (0 is coarsest) to supply [Do not supply 0] 147 148 Output Parameter: 149 . mat - the interpolation matrix, can be NULL 150 151 Level: advanced 152 153 .keywords: MG, get, multigrid, interpolation, level 154 155 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 156 @*/ 157 PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 158 { 159 PC_MG *mg = (PC_MG*)pc->data; 160 PC_MG_Levels **mglevels = mg->levels; 161 PetscErrorCode ierr; 162 163 PetscFunctionBegin; 164 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 165 if (mat) PetscValidPointer(mat,3); 166 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 167 if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 168 if (!mglevels[l]->interpolate) { 169 if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 170 ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 171 } 172 if (mat) *mat = mglevels[l]->interpolate; 173 PetscFunctionReturn(0); 174 } 175 176 /*@ 177 PCMGSetRestriction - Sets the function to be used to restrict dual vectors 178 from level l to l-1. 179 180 Logically Collective on PC and Mat 181 182 Input Parameters: 183 + pc - the multigrid context 184 . l - the level (0 is coarsest) to supply [Do not supply 0] 185 - mat - the restriction matrix 186 187 Level: advanced 188 189 Notes: 190 Usually this is the same matrix used also to set the interpolation 191 for the same level. 192 193 One can pass in the interpolation matrix or its transpose; PETSc figures 194 out from the matrix size which one it is. 195 196 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 197 is used. 198 199 .keywords: MG, set, multigrid, restriction, level 200 201 .seealso: PCMGSetInterpolation(), PCMGetSetInjection() 202 @*/ 203 PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 204 { 205 PetscErrorCode ierr; 206 PC_MG *mg = (PC_MG*)pc->data; 207 PC_MG_Levels **mglevels = mg->levels; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 211 PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 212 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 213 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 214 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 215 ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 216 217 mglevels[l]->restrct = mat; 218 PetscFunctionReturn(0); 219 } 220 221 /*@ 222 PCMGGetRestriction - Gets the function to be used to restrict dual vectors 223 from level l to l-1. 224 225 Logically Collective on PC and Mat 226 227 Input Parameters: 228 + pc - the multigrid context 229 - l - the level (0 is coarsest) to supply [Do not supply 0] 230 231 Output Parameter: 232 . mat - the restriction matrix 233 234 Level: advanced 235 236 .keywords: MG, get, multigrid, restriction, level 237 238 .seealso: PCMGGetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGGetInjection() 239 @*/ 240 PetscErrorCode PCMGGetRestriction(PC pc,PetscInt l,Mat *mat) 241 { 242 PC_MG *mg = (PC_MG*)pc->data; 243 PC_MG_Levels **mglevels = mg->levels; 244 PetscErrorCode ierr; 245 246 PetscFunctionBegin; 247 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 248 if (mat) PetscValidPointer(mat,3); 249 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 250 if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 251 if (!mglevels[l]->restrct) { 252 if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 253 ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 254 } 255 if (mat) *mat = mglevels[l]->restrct; 256 PetscFunctionReturn(0); 257 } 258 259 /*@ 260 PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 261 262 Logically Collective on PC and Vec 263 264 Input Parameters: 265 + pc - the multigrid context 266 - l - the level (0 is coarsest) to supply [Do not supply 0] 267 . rscale - the scaling 268 269 Level: advanced 270 271 Notes: 272 When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGSetInjection() to control moving primal vectors. 273 274 .keywords: MG, set, multigrid, restriction, level 275 276 .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 277 @*/ 278 PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 279 { 280 PetscErrorCode ierr; 281 PC_MG *mg = (PC_MG*)pc->data; 282 PC_MG_Levels **mglevels = mg->levels; 283 284 PetscFunctionBegin; 285 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 286 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 287 if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 288 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 289 ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 290 291 mglevels[l]->rscale = rscale; 292 PetscFunctionReturn(0); 293 } 294 295 /*@ 296 PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 297 298 Collective on PC 299 300 Input Parameters: 301 + pc - the multigrid context 302 . rscale - the scaling 303 - l - the level (0 is coarsest) to supply [Do not supply 0] 304 305 Level: advanced 306 307 Notes: 308 When evaluating a function on a coarse level one does not want to do F(R * x) one does F(rscale * R * x) where rscale is 1 over the row sums of R. It is preferable to use PCMGGetInjection() to control moving primal vectors. 309 310 .keywords: MG, set, multigrid, restriction, level 311 312 .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 313 @*/ 314 PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 315 { 316 PetscErrorCode ierr; 317 PC_MG *mg = (PC_MG*)pc->data; 318 PC_MG_Levels **mglevels = mg->levels; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 322 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 323 if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 324 if (!mglevels[l]->rscale) { 325 Mat R; 326 Vec X,Y,coarse,fine; 327 PetscInt M,N; 328 ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 329 ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 330 ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 331 if (M < N) { 332 fine = X; 333 coarse = Y; 334 } else if (N < M) { 335 fine = Y; coarse = X; 336 } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 337 ierr = VecSet(fine,1.);CHKERRQ(ierr); 338 ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 339 ierr = VecDestroy(&fine);CHKERRQ(ierr); 340 ierr = VecReciprocal(coarse);CHKERRQ(ierr); 341 mglevels[l]->rscale = coarse; 342 } 343 *rscale = mglevels[l]->rscale; 344 PetscFunctionReturn(0); 345 } 346 347 /*@ 348 PCMGSetInjection - Sets the function to be used to inject primal vectors 349 from level l to l-1. 350 351 Logically Collective on PC and Mat 352 353 Input Parameters: 354 + pc - the multigrid context 355 . l - the level (0 is coarsest) to supply [Do not supply 0] 356 - mat - the injection matrix 357 358 Level: advanced 359 360 .keywords: MG, set, multigrid, restriction, injection, level 361 362 .seealso: PCMGSetRestriction() 363 @*/ 364 PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 365 { 366 PetscErrorCode ierr; 367 PC_MG *mg = (PC_MG*)pc->data; 368 PC_MG_Levels **mglevels = mg->levels; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 372 PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 373 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 374 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 375 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 376 ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 377 378 mglevels[l]->inject = mat; 379 PetscFunctionReturn(0); 380 } 381 382 /*@ 383 PCMGGetInjection - Gets the function to be used to inject primal vectors 384 from level l to l-1. 385 386 Logically Collective on PC and Mat 387 388 Input Parameters: 389 + pc - the multigrid context 390 - l - the level (0 is coarsest) to supply [Do not supply 0] 391 392 Output Parameter: 393 . mat - the restriction matrix (may be NULL if no injection is available). 394 395 Level: advanced 396 397 .keywords: MG, get, multigrid, restriction, injection, level 398 399 .seealso: PCMGSetInjection(), PCMGetGetRestriction() 400 @*/ 401 PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 402 { 403 PC_MG *mg = (PC_MG*)pc->data; 404 PC_MG_Levels **mglevels = mg->levels; 405 406 PetscFunctionBegin; 407 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 408 if (mat) PetscValidPointer(mat,3); 409 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 410 if (l <= 0 || mg->nlevels <= l) SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Level %D must be in range {1,...,%D}",l,mg->nlevels-1); 411 if (mat) *mat = mglevels[l]->inject; 412 PetscFunctionReturn(0); 413 } 414 415 /*@ 416 PCMGGetSmoother - Gets the KSP context to be used as smoother for 417 both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 418 PCMGGetSmootherDown() to use different functions for pre- and 419 post-smoothing. 420 421 Not Collective, KSP returned is parallel if PC is 422 423 Input Parameters: 424 + pc - the multigrid context 425 - l - the level (0 is coarsest) to supply 426 427 Ouput Parameters: 428 . ksp - the smoother 429 430 Notes: 431 Once you have called this routine, you can call KSPSetOperators(ksp,...) on the resulting ksp to provide the operators for the smoother for this level. 432 You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 433 and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 434 435 Level: advanced 436 437 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 438 439 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 440 @*/ 441 PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 442 { 443 PC_MG *mg = (PC_MG*)pc->data; 444 PC_MG_Levels **mglevels = mg->levels; 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 448 *ksp = mglevels[l]->smoothd; 449 PetscFunctionReturn(0); 450 } 451 452 /*@ 453 PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 454 coarse grid correction (post-smoother). 455 456 Not Collective, KSP returned is parallel if PC is 457 458 Input Parameters: 459 + pc - the multigrid context 460 - l - the level (0 is coarsest) to supply 461 462 Ouput Parameters: 463 . ksp - the smoother 464 465 Level: advanced 466 467 Notes: 468 calling this will result in a different pre and post smoother so you may need to 469 set options on the pre smoother also 470 471 .keywords: MG, multigrid, get, smoother, up, post-smoother, level 472 473 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 474 @*/ 475 PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 476 { 477 PC_MG *mg = (PC_MG*)pc->data; 478 PC_MG_Levels **mglevels = mg->levels; 479 PetscErrorCode ierr; 480 const char *prefix; 481 MPI_Comm comm; 482 483 PetscFunctionBegin; 484 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 485 /* 486 This is called only if user wants a different pre-smoother from post. 487 Thus we check if a different one has already been allocated, 488 if not we allocate it. 489 */ 490 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 491 if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 492 KSPType ksptype; 493 PCType pctype; 494 PC ipc; 495 PetscReal rtol,abstol,dtol; 496 PetscInt maxits; 497 KSPNormType normtype; 498 ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 499 ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 500 ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 501 ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 502 ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 503 ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 504 ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 505 506 ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 507 ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 508 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 509 ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 510 ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 511 ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 512 ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 513 ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 514 ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 515 ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 516 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 517 ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 518 } 519 if (ksp) *ksp = mglevels[l]->smoothu; 520 PetscFunctionReturn(0); 521 } 522 523 /*@ 524 PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 525 coarse grid correction (pre-smoother). 526 527 Not Collective, KSP returned is parallel if PC is 528 529 Input Parameters: 530 + pc - the multigrid context 531 - l - the level (0 is coarsest) to supply 532 533 Ouput Parameters: 534 . ksp - the smoother 535 536 Level: advanced 537 538 Notes: 539 calling this will result in a different pre and post smoother so you may need to 540 set options on the post smoother also 541 542 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 543 544 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 545 @*/ 546 PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 547 { 548 PetscErrorCode ierr; 549 PC_MG *mg = (PC_MG*)pc->data; 550 PC_MG_Levels **mglevels = mg->levels; 551 552 PetscFunctionBegin; 553 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 554 /* make sure smoother up and down are different */ 555 if (l) { 556 ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 557 } 558 *ksp = mglevels[l]->smoothd; 559 PetscFunctionReturn(0); 560 } 561 562 /*@ 563 PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 564 565 Logically Collective on PC 566 567 Input Parameters: 568 + pc - the multigrid context 569 . l - the level (0 is coarsest) 570 - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 571 572 Level: advanced 573 574 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 575 576 .seealso: PCMGSetCycleType() 577 @*/ 578 PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 579 { 580 PC_MG *mg = (PC_MG*)pc->data; 581 PC_MG_Levels **mglevels = mg->levels; 582 583 PetscFunctionBegin; 584 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 585 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 586 PetscValidLogicalCollectiveInt(pc,l,2); 587 PetscValidLogicalCollectiveEnum(pc,c,3); 588 mglevels[l]->cycles = c; 589 PetscFunctionReturn(0); 590 } 591 592 /*@ 593 PCMGSetRhs - Sets the vector space to be used to store the right-hand side 594 on a particular level. 595 596 Logically Collective on PC and Vec 597 598 Input Parameters: 599 + pc - the multigrid context 600 . l - the level (0 is coarsest) this is to be used for 601 - c - the space 602 603 Level: advanced 604 605 Notes: 606 If this is not provided PETSc will automatically generate one. 607 608 You do not need to keep a reference to this vector if you do 609 not need it PCDestroy() will properly free it. 610 611 .keywords: MG, multigrid, set, right-hand-side, rhs, level 612 613 .seealso: PCMGSetX(), PCMGSetR() 614 @*/ 615 PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 616 { 617 PetscErrorCode ierr; 618 PC_MG *mg = (PC_MG*)pc->data; 619 PC_MG_Levels **mglevels = mg->levels; 620 621 PetscFunctionBegin; 622 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 623 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 624 if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 625 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 626 ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 627 628 mglevels[l]->b = c; 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 PCMGSetX - Sets the vector space to be used to store the solution on a 634 particular level. 635 636 Logically Collective on PC and Vec 637 638 Input Parameters: 639 + pc - the multigrid context 640 . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 641 - c - the space 642 643 Level: advanced 644 645 Notes: 646 If this is not provided PETSc will automatically generate one. 647 648 You do not need to keep a reference to this vector if you do 649 not need it PCDestroy() will properly free it. 650 651 .keywords: MG, multigrid, set, solution, level 652 653 .seealso: PCMGSetRhs(), PCMGSetR() 654 @*/ 655 PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 656 { 657 PetscErrorCode ierr; 658 PC_MG *mg = (PC_MG*)pc->data; 659 PC_MG_Levels **mglevels = mg->levels; 660 661 PetscFunctionBegin; 662 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 663 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 664 if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 665 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 666 ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 667 668 mglevels[l]->x = c; 669 PetscFunctionReturn(0); 670 } 671 672 /*@ 673 PCMGSetR - Sets the vector space to be used to store the residual on a 674 particular level. 675 676 Logically Collective on PC and Vec 677 678 Input Parameters: 679 + pc - the multigrid context 680 . l - the level (0 is coarsest) this is to be used for 681 - c - the space 682 683 Level: advanced 684 685 Notes: 686 If this is not provided PETSc will automatically generate one. 687 688 You do not need to keep a reference to this vector if you do 689 not need it PCDestroy() will properly free it. 690 691 .keywords: MG, multigrid, set, residual, level 692 @*/ 693 PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 694 { 695 PetscErrorCode ierr; 696 PC_MG *mg = (PC_MG*)pc->data; 697 PC_MG_Levels **mglevels = mg->levels; 698 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 701 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 702 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 703 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 704 ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 705 706 mglevels[l]->r = c; 707 PetscFunctionReturn(0); 708 } 709