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