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