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