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