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