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