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