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