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