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 PCMGSetOperators - Sets operator and preconditioning matrix for lth level 132 133 Logically Collective on PC and Mat 134 135 Input Parameters: 136 + pc - the multigrid context 137 . Amat - the operator 138 . pmat - the preconditioning operator 139 - l - the level (0 is the coarsest) to supply 140 141 Level: advanced 142 143 .keywords: multigrid, set, interpolate, level 144 145 .seealso: PCMGSetRestriction(), PCMGSetInterpolation() 146 @*/ 147 PetscErrorCode PCMGSetOperators(PC pc,PetscInt l,Mat Amat,Mat Pmat) 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 PetscValidHeaderSpecific(Amat,MAT_CLASSID,3); 156 PetscValidHeaderSpecific(Pmat,MAT_CLASSID,4); 157 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 158 ierr = KSPSetOperators(mglevels[l]->smoothd,Amat,Pmat);CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 /*@ 163 PCMGGetInterpolation - Gets the function to be used to calculate the 164 interpolation from l-1 to the lth level 165 166 Logically Collective on PC 167 168 Input Parameters: 169 + pc - the multigrid context 170 - l - the level (0 is coarsest) to supply [Do not supply 0] 171 172 Output Parameter: 173 . mat - the interpolation matrix, can be NULL 174 175 Level: advanced 176 177 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale() 178 @*/ 179 PetscErrorCode PCMGGetInterpolation(PC pc,PetscInt l,Mat *mat) 180 { 181 PC_MG *mg = (PC_MG*)pc->data; 182 PC_MG_Levels **mglevels = mg->levels; 183 PetscErrorCode ierr; 184 185 PetscFunctionBegin; 186 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 187 if (mat) PetscValidPointer(mat,3); 188 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 189 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); 190 if (!mglevels[l]->interpolate) { 191 if (!mglevels[l]->restrct) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetInterpolation() or PCMGSetRestriction()"); 192 ierr = PCMGSetInterpolation(pc,l,mglevels[l]->restrct);CHKERRQ(ierr); 193 } 194 if (mat) *mat = mglevels[l]->interpolate; 195 PetscFunctionReturn(0); 196 } 197 198 /*@ 199 PCMGGetInterpolations - Gets interpolation matrices for all levels (except level 0) 200 201 Logically Collective on PC 202 203 Input Parameters: 204 + pc - the multigrid context 205 206 Output Parameter: 207 - num_levels - the number of levels 208 . interpolations - the interpolation matrices (size of num_levels-1) 209 210 Notes: 211 No new matrices are created, and the interpolations are the references to the original ones 212 213 Level: advanced 214 215 .keywords: MG, get, multigrid, interpolation, level 216 217 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale(), PCMGGetInterpolation() 218 @*/ 219 PetscErrorCode PCMGGetInterpolations(PC pc,PetscInt *num_levels,Mat *interpolations[]) 220 { 221 PC_MG *mg = (PC_MG*)pc->data; 222 PC_MG_Levels **mglevels = mg->levels; 223 Mat *mat; 224 PetscInt l; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 229 PetscValidPointer(num_levels,2); 230 PetscValidPointer(interpolations,3); 231 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 232 ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 233 for (l=1; l< mg->nlevels; l++) { 234 mat[l-1] = mglevels[l]->interpolate; 235 ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 236 } 237 *num_levels = mg->nlevels; 238 *interpolations = mat; 239 PetscFunctionReturn(0); 240 } 241 242 /*@ 243 PCMGGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 244 245 Logically Collective on PC 246 247 Input Parameters: 248 + pc - the multigrid context 249 250 Output Parameter: 251 - num_levels - the number of levels 252 . coarseOperators - the coarse operator matrices (size of num_levels-1) 253 254 Notes: 255 No new matrices are created, and the coarse operator matrices are the references to the original ones 256 257 Level: advanced 258 259 .keywords: MG, get, multigrid, interpolation, level 260 261 .seealso: PCMGGetRestriction(), PCMGSetInterpolation(), PCMGGetRScale(), PCMGGetInterpolation(), PCMGGetInterpolations() 262 @*/ 263 PetscErrorCode PCMGGetCoarseOperators(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 264 { 265 PC_MG *mg = (PC_MG*)pc->data; 266 PC_MG_Levels **mglevels = mg->levels; 267 PetscInt l; 268 Mat *mat; 269 PetscErrorCode ierr; 270 271 PetscFunctionBegin; 272 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 273 PetscValidPointer(num_levels,2); 274 PetscValidPointer(coarseOperators,3); 275 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 276 ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 277 for (l=0; l<mg->nlevels-1; l++) { 278 ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 279 ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 280 } 281 *num_levels = mg->nlevels; 282 *coarseOperators = mat; 283 PetscFunctionReturn(0); 284 } 285 286 /*@ 287 PCMGSetRestriction - Sets the function to be used to restrict dual vectors 288 from level l to l-1. 289 290 Logically Collective on PC and Mat 291 292 Input Parameters: 293 + pc - the multigrid context 294 . l - the level (0 is coarsest) to supply [Do not supply 0] 295 - mat - the restriction matrix 296 297 Level: advanced 298 299 Notes: 300 Usually this is the same matrix used also to set the interpolation 301 for the same level. 302 303 One can pass in the interpolation matrix or its transpose; PETSc figures 304 out from the matrix size which one it is. 305 306 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 307 is used. 308 309 .seealso: PCMGSetInterpolation(), PCMGetSetInjection() 310 @*/ 311 PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 312 { 313 PetscErrorCode ierr; 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 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 321 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 322 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 323 ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 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 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 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 { 348 PC_MG *mg = (PC_MG*)pc->data; 349 PC_MG_Levels **mglevels = mg->levels; 350 PetscErrorCode ierr; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 354 if (mat) PetscValidPointer(mat,3); 355 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 356 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); 357 if (!mglevels[l]->restrct) { 358 if (!mglevels[l]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call PCMGSetRestriction() or PCMGSetInterpolation()"); 359 ierr = PCMGSetRestriction(pc,l,mglevels[l]->interpolate);CHKERRQ(ierr); 360 } 361 if (mat) *mat = mglevels[l]->restrct; 362 PetscFunctionReturn(0); 363 } 364 365 /*@ 366 PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 367 368 Logically Collective on PC and Vec 369 370 Input Parameters: 371 + pc - the multigrid context 372 - l - the level (0 is coarsest) to supply [Do not supply 0] 373 . rscale - the scaling 374 375 Level: advanced 376 377 Notes: 378 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. 379 380 .seealso: PCMGSetInterpolation(), PCMGSetRestriction(), PCMGGetRScale(), PCMGSetInjection() 381 @*/ 382 PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 383 { 384 PetscErrorCode ierr; 385 PC_MG *mg = (PC_MG*)pc->data; 386 PC_MG_Levels **mglevels = mg->levels; 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 391 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); 392 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 393 ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 394 395 mglevels[l]->rscale = rscale; 396 PetscFunctionReturn(0); 397 } 398 399 /*@ 400 PCMGGetRScale - Gets the pointwise scaling for the restriction operator from level l to l-1. 401 402 Collective on PC 403 404 Input Parameters: 405 + pc - the multigrid context 406 . rscale - the scaling 407 - l - the level (0 is coarsest) to supply [Do not supply 0] 408 409 Level: advanced 410 411 Notes: 412 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. 413 414 .seealso: PCMGSetInterpolation(), PCMGGetRestriction(), PCMGGetInjection() 415 @*/ 416 PetscErrorCode PCMGGetRScale(PC pc,PetscInt l,Vec *rscale) 417 { 418 PetscErrorCode ierr; 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 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 425 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); 426 if (!mglevels[l]->rscale) { 427 Mat R; 428 Vec X,Y,coarse,fine; 429 PetscInt M,N; 430 ierr = PCMGGetRestriction(pc,l,&R);CHKERRQ(ierr); 431 ierr = MatCreateVecs(R,&X,&Y);CHKERRQ(ierr); 432 ierr = MatGetSize(R,&M,&N);CHKERRQ(ierr); 433 if (M < N) { 434 fine = X; 435 coarse = Y; 436 } else if (N < M) { 437 fine = Y; coarse = X; 438 } else SETERRQ(PetscObjectComm((PetscObject)R),PETSC_ERR_SUP,"Restriction matrix is square, cannot determine which Vec is coarser"); 439 ierr = VecSet(fine,1.);CHKERRQ(ierr); 440 ierr = MatRestrict(R,fine,coarse);CHKERRQ(ierr); 441 ierr = VecDestroy(&fine);CHKERRQ(ierr); 442 ierr = VecReciprocal(coarse);CHKERRQ(ierr); 443 mglevels[l]->rscale = coarse; 444 } 445 *rscale = mglevels[l]->rscale; 446 PetscFunctionReturn(0); 447 } 448 449 /*@ 450 PCMGSetInjection - Sets the function to be used to inject primal vectors 451 from level l to l-1. 452 453 Logically Collective on PC and Mat 454 455 Input Parameters: 456 + pc - the multigrid context 457 . l - the level (0 is coarsest) to supply [Do not supply 0] 458 - mat - the injection matrix 459 460 Level: advanced 461 462 .seealso: PCMGSetRestriction() 463 @*/ 464 PetscErrorCode PCMGSetInjection(PC pc,PetscInt l,Mat mat) 465 { 466 PetscErrorCode ierr; 467 PC_MG *mg = (PC_MG*)pc->data; 468 PC_MG_Levels **mglevels = mg->levels; 469 470 PetscFunctionBegin; 471 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 472 PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 473 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 474 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 475 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 476 ierr = MatDestroy(&mglevels[l]->inject);CHKERRQ(ierr); 477 478 mglevels[l]->inject = mat; 479 PetscFunctionReturn(0); 480 } 481 482 /*@ 483 PCMGGetInjection - Gets the function to be used to inject primal vectors 484 from level l to l-1. 485 486 Logically Collective on PC and Mat 487 488 Input Parameters: 489 + pc - the multigrid context 490 - l - the level (0 is coarsest) to supply [Do not supply 0] 491 492 Output Parameter: 493 . mat - the restriction matrix (may be NULL if no injection is available). 494 495 Level: advanced 496 497 .seealso: PCMGSetInjection(), PCMGetGetRestriction() 498 @*/ 499 PetscErrorCode PCMGGetInjection(PC pc,PetscInt l,Mat *mat) 500 { 501 PC_MG *mg = (PC_MG*)pc->data; 502 PC_MG_Levels **mglevels = mg->levels; 503 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 506 if (mat) PetscValidPointer(mat,3); 507 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 508 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); 509 if (mat) *mat = mglevels[l]->inject; 510 PetscFunctionReturn(0); 511 } 512 513 /*@ 514 PCMGGetSmoother - Gets the KSP context to be used as smoother for 515 both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 516 PCMGGetSmootherDown() to use different functions for pre- and 517 post-smoothing. 518 519 Not Collective, KSP returned is parallel if PC is 520 521 Input Parameters: 522 + pc - the multigrid context 523 - l - the level (0 is coarsest) to supply 524 525 Ouput Parameters: 526 . ksp - the smoother 527 528 Notes: 529 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. 530 You can also modify smoother options by calling the various KSPSetXXX() options on this ksp. In addition you can call KSPGetPC(ksp,&pc) 531 and modify PC options for the smoother; for example PCSetType(pc,PCSOR); to use SOR smoothing. 532 533 Level: advanced 534 535 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown(), PCMGGetCoarseSolve() 536 @*/ 537 PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 538 { 539 PC_MG *mg = (PC_MG*)pc->data; 540 PC_MG_Levels **mglevels = mg->levels; 541 542 PetscFunctionBegin; 543 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 544 *ksp = mglevels[l]->smoothd; 545 PetscFunctionReturn(0); 546 } 547 548 /*@ 549 PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 550 coarse grid correction (post-smoother). 551 552 Not Collective, KSP returned is parallel if PC is 553 554 Input Parameters: 555 + pc - the multigrid context 556 - l - the level (0 is coarsest) to supply 557 558 Ouput Parameters: 559 . ksp - the smoother 560 561 Level: advanced 562 563 Notes: 564 calling this will result in a different pre and post smoother so you may need to 565 set options on the pre smoother also 566 567 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 568 @*/ 569 PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 570 { 571 PC_MG *mg = (PC_MG*)pc->data; 572 PC_MG_Levels **mglevels = mg->levels; 573 PetscErrorCode ierr; 574 const char *prefix; 575 MPI_Comm comm; 576 577 PetscFunctionBegin; 578 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 579 /* 580 This is called only if user wants a different pre-smoother from post. 581 Thus we check if a different one has already been allocated, 582 if not we allocate it. 583 */ 584 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 585 if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 586 KSPType ksptype; 587 PCType pctype; 588 PC ipc; 589 PetscReal rtol,abstol,dtol; 590 PetscInt maxits; 591 KSPNormType normtype; 592 ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 593 ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 594 ierr = KSPGetTolerances(mglevels[l]->smoothd,&rtol,&abstol,&dtol,&maxits);CHKERRQ(ierr); 595 ierr = KSPGetType(mglevels[l]->smoothd,&ksptype);CHKERRQ(ierr); 596 ierr = KSPGetNormType(mglevels[l]->smoothd,&normtype);CHKERRQ(ierr); 597 ierr = KSPGetPC(mglevels[l]->smoothd,&ipc);CHKERRQ(ierr); 598 ierr = PCGetType(ipc,&pctype);CHKERRQ(ierr); 599 600 ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 601 ierr = KSPSetErrorIfNotConverged(mglevels[l]->smoothu,pc->erroriffailure);CHKERRQ(ierr); 602 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 603 ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 604 ierr = KSPSetTolerances(mglevels[l]->smoothu,rtol,abstol,dtol,maxits);CHKERRQ(ierr); 605 ierr = KSPSetType(mglevels[l]->smoothu,ksptype);CHKERRQ(ierr); 606 ierr = KSPSetNormType(mglevels[l]->smoothu,normtype);CHKERRQ(ierr); 607 ierr = KSPSetConvergenceTest(mglevels[l]->smoothu,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 608 ierr = KSPGetPC(mglevels[l]->smoothu,&ipc);CHKERRQ(ierr); 609 ierr = PCSetType(ipc,pctype);CHKERRQ(ierr); 610 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[l]->smoothu);CHKERRQ(ierr); 611 ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[l]->smoothu, PetscMGLevelId, mglevels[l]->level);CHKERRQ(ierr); 612 } 613 if (ksp) *ksp = mglevels[l]->smoothu; 614 PetscFunctionReturn(0); 615 } 616 617 /*@ 618 PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 619 coarse grid correction (pre-smoother). 620 621 Not Collective, KSP returned is parallel if PC is 622 623 Input Parameters: 624 + pc - the multigrid context 625 - l - the level (0 is coarsest) to supply 626 627 Ouput Parameters: 628 . ksp - the smoother 629 630 Level: advanced 631 632 Notes: 633 calling this will result in a different pre and post smoother so you may need to 634 set options on the post smoother also 635 636 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 637 @*/ 638 PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 639 { 640 PetscErrorCode ierr; 641 PC_MG *mg = (PC_MG*)pc->data; 642 PC_MG_Levels **mglevels = mg->levels; 643 644 PetscFunctionBegin; 645 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 646 /* make sure smoother up and down are different */ 647 if (l) { 648 ierr = PCMGGetSmootherUp(pc,l,NULL);CHKERRQ(ierr); 649 } 650 *ksp = mglevels[l]->smoothd; 651 PetscFunctionReturn(0); 652 } 653 654 /*@ 655 PCMGSetCycleTypeOnLevel - Sets the type of cycle (aka cycle index) to run on the specified level. 656 657 Logically Collective on PC 658 659 Input Parameters: 660 + pc - the multigrid context 661 . l - the level (0 is coarsest) 662 - c - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 663 664 Level: advanced 665 666 .seealso: PCMGSetCycleType() 667 @*/ 668 PetscErrorCode PCMGSetCycleTypeOnLevel(PC pc,PetscInt l,PCMGCycleType c) 669 { 670 PC_MG *mg = (PC_MG*)pc->data; 671 PC_MG_Levels **mglevels = mg->levels; 672 673 PetscFunctionBegin; 674 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 675 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 676 PetscValidLogicalCollectiveInt(pc,l,2); 677 PetscValidLogicalCollectiveEnum(pc,c,3); 678 mglevels[l]->cycles = c; 679 PetscFunctionReturn(0); 680 } 681 682 /*@ 683 PCMGSetRhs - Sets the vector space to be used to store the right-hand side 684 on a particular level. 685 686 Logically Collective on PC and Vec 687 688 Input Parameters: 689 + pc - the multigrid context 690 . l - the level (0 is coarsest) this is to be used for 691 - c - the space 692 693 Level: advanced 694 695 Notes: 696 If this is not provided PETSc will automatically generate one. 697 698 You do not need to keep a reference to this vector if you do 699 not need it PCDestroy() will properly free it. 700 701 .seealso: PCMGSetX(), PCMGSetR() 702 @*/ 703 PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 704 { 705 PetscErrorCode ierr; 706 PC_MG *mg = (PC_MG*)pc->data; 707 PC_MG_Levels **mglevels = mg->levels; 708 709 PetscFunctionBegin; 710 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 711 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 712 if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 713 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 714 ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 715 716 mglevels[l]->b = c; 717 PetscFunctionReturn(0); 718 } 719 720 /*@ 721 PCMGSetX - Sets the vector space to be used to store the solution on a 722 particular level. 723 724 Logically Collective on PC and Vec 725 726 Input Parameters: 727 + pc - the multigrid context 728 . l - the level (0 is coarsest) this is to be used for (do not supply the finest level) 729 - c - the space 730 731 Level: advanced 732 733 Notes: 734 If this is not provided PETSc will automatically generate one. 735 736 You do not need to keep a reference to this vector if you do 737 not need it PCDestroy() will properly free it. 738 739 .seealso: PCMGSetRhs(), PCMGSetR() 740 @*/ 741 PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 742 { 743 PetscErrorCode ierr; 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 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 750 if (l == mglevels[0]->levels-1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 751 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 752 ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 753 754 mglevels[l]->x = c; 755 PetscFunctionReturn(0); 756 } 757 758 /*@ 759 PCMGSetR - Sets the vector space to be used to store the residual on a 760 particular level. 761 762 Logically Collective on PC and Vec 763 764 Input Parameters: 765 + pc - the multigrid context 766 . l - the level (0 is coarsest) this is to be used for 767 - c - the space 768 769 Level: advanced 770 771 Notes: 772 If this is not provided PETSc will automatically generate one. 773 774 You do not need to keep a reference to this vector if you do 775 not need it PCDestroy() will properly free it. 776 777 @*/ 778 PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 779 { 780 PetscErrorCode ierr; 781 PC_MG *mg = (PC_MG*)pc->data; 782 PC_MG_Levels **mglevels = mg->levels; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 786 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 787 if (!l) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 788 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 789 ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 790 791 mglevels[l]->r = c; 792 PetscFunctionReturn(0); 793 } 794