1 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 2 3 4 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES*); 5 6 /* -------------- functions called on the fine level -------------- */ 7 8 /*@ 9 SNESFASSetType - Sets the update and correction type used for FAS. 10 11 Logically Collective 12 13 Input Parameters: 14 + snes - FAS context 15 - fastype - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE 16 17 Level: intermediate 18 19 .seealso: PCMGSetType() 20 @*/ 21 PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 22 { 23 SNES_FAS *fas = (SNES_FAS*)snes->data; 24 PetscErrorCode ierr; 25 26 PetscFunctionBegin; 27 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 28 PetscValidLogicalCollectiveEnum(snes,fastype,2); 29 fas->fastype = fastype; 30 if (fas->next) { 31 ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr); 32 } 33 PetscFunctionReturn(0); 34 } 35 36 37 /*@ 38 SNESFASGetType - Sets the update and correction type used for FAS. 39 40 Logically Collective 41 42 Input Parameters: 43 . snes - FAS context 44 45 Output Parameters: 46 . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 47 48 Level: intermediate 49 50 .seealso: PCMGSetType() 51 @*/ 52 PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype) 53 { 54 SNES_FAS *fas = (SNES_FAS*)snes->data; 55 56 PetscFunctionBegin; 57 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 58 PetscValidPointer(fastype, 2); 59 *fastype = fas->fastype; 60 PetscFunctionReturn(0); 61 } 62 63 /*@C 64 SNESFASSetLevels - Sets the number of levels to use with FAS. 65 Must be called before any other FAS routine. 66 67 Input Parameters: 68 + snes - the snes context 69 . levels - the number of levels 70 - comms - optional communicators for each level; this is to allow solving the coarser 71 problems on smaller sets of processors. Use NULL_OBJECT for default in 72 Fortran. 73 74 Level: intermediate 75 76 Notes: 77 If the number of levels is one then the multigrid uses the -fas_levels prefix 78 for setting the level options rather than the -fas_coarse prefix. 79 80 .keywords: FAS, MG, set, levels, multigrid 81 82 .seealso: SNESFASGetLevels() 83 @*/ 84 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) 85 { 86 PetscErrorCode ierr; 87 PetscInt i; 88 const char *optionsprefix; 89 char tprefix[128]; 90 SNES_FAS *fas = (SNES_FAS*)snes->data; 91 SNES prevsnes; 92 MPI_Comm comm; 93 94 PetscFunctionBegin; 95 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 96 if (levels == fas->levels) { 97 if (!comms) PetscFunctionReturn(0); 98 } 99 /* user has changed the number of levels; reset */ 100 ierr = SNESReset(snes);CHKERRQ(ierr); 101 /* destroy any coarser levels if necessary */ 102 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 103 fas->next = NULL; 104 fas->previous = NULL; 105 prevsnes = snes; 106 /* setup the finest level */ 107 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 108 ierr = PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);CHKERRQ(ierr); 109 for (i = levels - 1; i >= 0; i--) { 110 if (comms) comm = comms[i]; 111 fas->level = i; 112 fas->levels = levels; 113 fas->fine = snes; 114 fas->next = NULL; 115 if (i > 0) { 116 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 117 ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 118 sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level); 119 ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 120 ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 121 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 122 ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 123 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 124 ierr = PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);CHKERRQ(ierr); 125 126 ((SNES_FAS*)fas->next->data)->previous = prevsnes; 127 128 prevsnes = fas->next; 129 fas = (SNES_FAS*)prevsnes->data; 130 } 131 } 132 PetscFunctionReturn(0); 133 } 134 135 136 /*@ 137 SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 138 139 Input Parameter: 140 . snes - the nonlinear solver context 141 142 Output parameter: 143 . levels - the number of levels 144 145 Level: advanced 146 147 .keywords: MG, get, levels, multigrid 148 149 .seealso: SNESFASSetLevels(), PCMGGetLevels() 150 @*/ 151 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 152 { 153 SNES_FAS * fas = (SNES_FAS*)snes->data; 154 155 PetscFunctionBegin; 156 *levels = fas->levels; 157 PetscFunctionReturn(0); 158 } 159 160 161 /*@ 162 SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 163 level of the FAS hierarchy. 164 165 Input Parameters: 166 + snes - the multigrid context 167 level - the level to get 168 - lsnes - whether to use the nonlinear smoother or not 169 170 Level: advanced 171 172 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 173 174 .seealso: SNESFASSetLevels(), SNESFASGetLevels() 175 @*/ 176 PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) 177 { 178 SNES_FAS *fas = (SNES_FAS*)snes->data; 179 PetscInt i; 180 181 PetscFunctionBegin; 182 if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 183 if (fas->level != fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level); 184 185 *lsnes = snes; 186 for (i = fas->level; i > level; i--) { 187 *lsnes = fas->next; 188 fas = (SNES_FAS*)(*lsnes)->data; 189 } 190 if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 191 PetscFunctionReturn(0); 192 } 193 194 /*@ 195 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 196 use on all levels. 197 198 Logically Collective on SNES 199 200 Input Parameters: 201 + snes - the multigrid context 202 - n - the number of smoothing steps 203 204 Options Database Key: 205 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 206 207 Level: advanced 208 209 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 210 211 .seealso: SNESFASSetNumberSmoothDown() 212 @*/ 213 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 214 { 215 SNES_FAS *fas = (SNES_FAS*)snes->data; 216 PetscErrorCode ierr; 217 218 PetscFunctionBegin; 219 fas->max_up_it = n; 220 if (!fas->smoothu && fas->level != 0) { 221 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 222 } 223 if (fas->smoothu) { 224 ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 225 } 226 if (fas->next) { 227 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 228 } 229 PetscFunctionReturn(0); 230 } 231 232 /*@ 233 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 234 use on all levels. 235 236 Logically Collective on SNES 237 238 Input Parameters: 239 + snes - the multigrid context 240 - n - the number of smoothing steps 241 242 Options Database Key: 243 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 244 245 Level: advanced 246 247 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 248 249 .seealso: SNESFASSetNumberSmoothUp() 250 @*/ 251 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 252 { 253 SNES_FAS *fas = (SNES_FAS*)snes->data; 254 PetscErrorCode ierr = 0; 255 256 PetscFunctionBegin; 257 if (!fas->smoothd) { 258 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 259 } 260 ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 261 262 fas->max_down_it = n; 263 if (fas->next) { 264 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 265 } 266 PetscFunctionReturn(0); 267 } 268 269 270 /*@ 271 SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep 272 273 Logically Collective on SNES 274 275 Input Parameters: 276 + snes - the multigrid context 277 - n - the number of smoothing steps 278 279 Options Database Key: 280 . -snes_fas_continuation - sets continuation to true 281 282 Level: advanced 283 284 Notes: This sets the prefix on the upsweep smoothers to -fas_continuation 285 286 .keywords: FAS, MG, smoother, continuation 287 288 .seealso: SNESFAS 289 @*/ 290 PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation) 291 { 292 const char *optionsprefix; 293 char tprefix[128]; 294 SNES_FAS *fas = (SNES_FAS*)snes->data; 295 PetscErrorCode ierr = 0; 296 297 PetscFunctionBegin; 298 ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 299 if (!fas->smoothu) { 300 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 301 } 302 sprintf(tprefix,"fas_levels_continuation_"); 303 ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr); 304 ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr); 305 ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr); 306 ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr); 307 fas->continuation = continuation; 308 if (fas->next) { 309 ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr); 310 } 311 PetscFunctionReturn(0); 312 } 313 314 315 /*@ 316 SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 317 complicated cycling. 318 319 Logically Collective on SNES 320 321 Input Parameters: 322 + snes - the multigrid context 323 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 324 325 Options Database Key: 326 . -snes_fas_cycles 1 or 2 327 328 Level: advanced 329 330 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 331 332 .seealso: SNESFASSetCyclesOnLevel() 333 @*/ 334 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 335 { 336 SNES_FAS *fas = (SNES_FAS*)snes->data; 337 PetscErrorCode ierr; 338 PetscBool isFine; 339 340 PetscFunctionBegin; 341 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 342 343 fas->n_cycles = cycles; 344 if (!isFine) { 345 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 346 } 347 if (fas->next) { 348 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 349 } 350 PetscFunctionReturn(0); 351 } 352 353 354 /*@ 355 SNESFASSetMonitor - Sets the method-specific cycle monitoring 356 357 Logically Collective on SNES 358 359 Input Parameters: 360 + snes - the FAS context 361 . vf - viewer and format structure (may be NULL if flg is FALSE) 362 - flg - monitor or not 363 364 Level: advanced 365 366 .keywords: FAS, monitor 367 368 .seealso: SNESFASSetCyclesOnLevel() 369 @*/ 370 PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 371 { 372 SNES_FAS *fas = (SNES_FAS*)snes->data; 373 PetscErrorCode ierr; 374 PetscBool isFine; 375 PetscInt i, levels = fas->levels; 376 SNES levelsnes; 377 378 PetscFunctionBegin; 379 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 380 if (isFine) { 381 for (i = 0; i < levels; i++) { 382 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 383 fas = (SNES_FAS*)levelsnes->data; 384 if (flg) { 385 /* set the monitors for the upsmoother and downsmoother */ 386 ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 387 /* Only register destroy on finest level */ 388 ierr = SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));CHKERRQ(ierr); 389 } else if (i != fas->levels - 1) { 390 /* unset the monitors on the coarse levels */ 391 ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 392 } 393 } 394 } 395 PetscFunctionReturn(0); 396 } 397 398 /*@ 399 SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels 400 401 Logically Collective on SNES 402 403 Input Parameters: 404 + snes - the FAS context 405 - flg - monitor or not 406 407 Level: advanced 408 409 .keywords: FAS, logging 410 411 .seealso: SNESFASSetMonitor() 412 @*/ 413 PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 414 { 415 SNES_FAS *fas = (SNES_FAS*)snes->data; 416 PetscErrorCode ierr; 417 PetscBool isFine; 418 PetscInt i, levels = fas->levels; 419 SNES levelsnes; 420 char eventname[128]; 421 422 PetscFunctionBegin; 423 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 424 if (isFine) { 425 for (i = 0; i < levels; i++) { 426 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 427 fas = (SNES_FAS*)levelsnes->data; 428 if (flg) { 429 sprintf(eventname,"FASSetup %d",(int)i); 430 ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr); 431 sprintf(eventname,"FASSmooth %d",(int)i); 432 ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr); 433 sprintf(eventname,"FASResid %d",(int)i); 434 ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr); 435 if (i) { 436 sprintf(eventname,"FASInterp %d",(int)i); 437 ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr); 438 } 439 } else { 440 fas->eventsmoothsetup = 0; 441 fas->eventsmoothsolve = 0; 442 fas->eventresidual = 0; 443 fas->eventinterprestrict = 0; 444 } 445 } 446 } 447 PetscFunctionReturn(0); 448 } 449 450 /* 451 Creates the default smoother type. 452 453 This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 454 455 */ 456 PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 457 { 458 SNES_FAS *fas; 459 const char *optionsprefix; 460 char tprefix[128]; 461 PetscErrorCode ierr; 462 SNES nsmooth; 463 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 466 fas = (SNES_FAS*)snes->data; 467 ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 468 /* create the default smoother */ 469 ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr); 470 if (fas->level == 0) { 471 sprintf(tprefix,"fas_coarse_"); 472 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 473 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 474 ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr); 475 ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 476 } else { 477 sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 478 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 479 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 480 ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 481 ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 482 } 483 ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 484 ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr); 485 ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 486 ierr = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr); 487 *smooth = nsmooth; 488 PetscFunctionReturn(0); 489 } 490 491 /* ------------- Functions called on a particular level ----------------- */ 492 493 /*@ 494 SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 495 496 Logically Collective on SNES 497 498 Input Parameters: 499 + snes - the multigrid context 500 . level - the level to set the number of cycles on 501 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 502 503 Level: advanced 504 505 .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 506 507 .seealso: SNESFASSetCycles() 508 @*/ 509 PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 510 { 511 SNES_FAS *fas = (SNES_FAS*)snes->data; 512 PetscErrorCode ierr; 513 514 PetscFunctionBegin; 515 fas->n_cycles = cycles; 516 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 517 PetscFunctionReturn(0); 518 } 519 520 521 /*@ 522 SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 523 524 Logically Collective on SNES 525 526 Input Parameters: 527 . snes - the multigrid context 528 529 Output Parameters: 530 . smooth - the smoother 531 532 Level: advanced 533 534 .keywords: SNES, FAS, get, smoother, multigrid 535 536 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 537 @*/ 538 PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 539 { 540 SNES_FAS *fas; 541 542 PetscFunctionBegin; 543 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 544 fas = (SNES_FAS*)snes->data; 545 *smooth = fas->smoothd; 546 PetscFunctionReturn(0); 547 } 548 /*@ 549 SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 550 551 Logically Collective on SNES 552 553 Input Parameters: 554 . snes - the multigrid context 555 556 Output Parameters: 557 . smoothu - the smoother 558 559 Notes: 560 Returns the downsmoother if no up smoother is available. This enables transparent 561 default behavior in the process of the solve. 562 563 Level: advanced 564 565 .keywords: SNES, FAS, get, smoother, multigrid 566 567 .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 568 @*/ 569 PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 570 { 571 SNES_FAS *fas; 572 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 575 fas = (SNES_FAS*)snes->data; 576 if (!fas->smoothu) *smoothu = fas->smoothd; 577 else *smoothu = fas->smoothu; 578 PetscFunctionReturn(0); 579 } 580 581 /*@ 582 SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 583 584 Logically Collective on SNES 585 586 Input Parameters: 587 . snes - the multigrid context 588 589 Output Parameters: 590 . smoothd - the smoother 591 592 Level: advanced 593 594 .keywords: SNES, FAS, get, smoother, multigrid 595 596 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 597 @*/ 598 PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 599 { 600 SNES_FAS *fas; 601 602 PetscFunctionBegin; 603 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 604 fas = (SNES_FAS*)snes->data; 605 *smoothd = fas->smoothd; 606 PetscFunctionReturn(0); 607 } 608 609 610 /*@ 611 SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 612 613 Logically Collective on SNES 614 615 Input Parameters: 616 . snes - the multigrid context 617 618 Output Parameters: 619 . correction - the coarse correction on this level 620 621 Notes: 622 Returns NULL on the coarsest level. 623 624 Level: advanced 625 626 .keywords: SNES, FAS, get, smoother, multigrid 627 628 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 629 @*/ 630 PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 631 { 632 SNES_FAS *fas; 633 634 PetscFunctionBegin; 635 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 636 fas = (SNES_FAS*)snes->data; 637 *correction = fas->next; 638 PetscFunctionReturn(0); 639 } 640 641 /*@ 642 SNESFASCycleGetInterpolation - Gets the interpolation on this level 643 644 Logically Collective on SNES 645 646 Input Parameters: 647 . snes - the multigrid context 648 649 Output Parameters: 650 . mat - the interpolation operator on this level 651 652 Level: developer 653 654 .keywords: SNES, FAS, get, smoother, multigrid 655 656 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 657 @*/ 658 PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 659 { 660 SNES_FAS *fas; 661 662 PetscFunctionBegin; 663 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 664 fas = (SNES_FAS*)snes->data; 665 *mat = fas->interpolate; 666 PetscFunctionReturn(0); 667 } 668 669 670 /*@ 671 SNESFASCycleGetRestriction - Gets the restriction on this level 672 673 Logically Collective on SNES 674 675 Input Parameters: 676 . snes - the multigrid context 677 678 Output Parameters: 679 . mat - the restriction operator on this level 680 681 Level: developer 682 683 .keywords: SNES, FAS, get, smoother, multigrid 684 685 .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 686 @*/ 687 PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 688 { 689 SNES_FAS *fas; 690 691 PetscFunctionBegin; 692 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 693 fas = (SNES_FAS*)snes->data; 694 *mat = fas->restrct; 695 PetscFunctionReturn(0); 696 } 697 698 699 /*@ 700 SNESFASCycleGetInjection - Gets the injection on this level 701 702 Logically Collective on SNES 703 704 Input Parameters: 705 . snes - the multigrid context 706 707 Output Parameters: 708 . mat - the restriction operator on this level 709 710 Level: developer 711 712 .keywords: SNES, FAS, get, smoother, multigrid 713 714 .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 715 @*/ 716 PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 717 { 718 SNES_FAS *fas; 719 720 PetscFunctionBegin; 721 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 722 fas = (SNES_FAS*)snes->data; 723 *mat = fas->inject; 724 PetscFunctionReturn(0); 725 } 726 727 /*@ 728 SNESFASCycleGetRScale - Gets the injection on this level 729 730 Logically Collective on SNES 731 732 Input Parameters: 733 . snes - the multigrid context 734 735 Output Parameters: 736 . mat - the restriction operator on this level 737 738 Level: developer 739 740 .keywords: SNES, FAS, get, smoother, multigrid 741 742 .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 743 @*/ 744 PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 745 { 746 SNES_FAS *fas; 747 748 PetscFunctionBegin; 749 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 750 fas = (SNES_FAS*)snes->data; 751 *vec = fas->rscale; 752 PetscFunctionReturn(0); 753 } 754 755 /*@ 756 SNESFASCycleIsFine - Determines if a given cycle is the fine level. 757 758 Logically Collective on SNES 759 760 Input Parameters: 761 . snes - the FAS context 762 763 Output Parameters: 764 . flg - indicates if this is the fine level or not 765 766 Level: advanced 767 768 .keywords: SNES, FAS 769 770 .seealso: SNESFASSetLevels() 771 @*/ 772 PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 773 { 774 SNES_FAS *fas; 775 776 PetscFunctionBegin; 777 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 778 fas = (SNES_FAS*)snes->data; 779 if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 780 else *flg = PETSC_FALSE; 781 PetscFunctionReturn(0); 782 } 783 784 /* ---------- functions called on the finest level that return level-specific information ---------- */ 785 786 /*@ 787 SNESFASSetInterpolation - Sets the function to be used to calculate the 788 interpolation from l-1 to the lth level 789 790 Input Parameters: 791 + snes - the multigrid context 792 . mat - the interpolation operator 793 - level - the level (0 is coarsest) to supply [do not supply 0] 794 795 Level: advanced 796 797 Notes: 798 Usually this is the same matrix used also to set the restriction 799 for the same level. 800 801 One can pass in the interpolation matrix or its transpose; PETSc figures 802 out from the matrix size which one it is. 803 804 .keywords: FAS, multigrid, set, interpolate, level 805 806 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 807 @*/ 808 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 809 { 810 SNES_FAS *fas; 811 PetscErrorCode ierr; 812 SNES levelsnes; 813 814 PetscFunctionBegin; 815 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 816 fas = (SNES_FAS*)levelsnes->data; 817 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 818 ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 819 820 fas->interpolate = mat; 821 PetscFunctionReturn(0); 822 } 823 824 /*@ 825 SNESFASGetInterpolation - Gets the matrix used to calculate the 826 interpolation from l-1 to the lth level 827 828 Input Parameters: 829 + snes - the multigrid context 830 - level - the level (0 is coarsest) to supply [do not supply 0] 831 832 Output Parameters: 833 . mat - the interpolation operator 834 835 Level: advanced 836 837 .keywords: FAS, multigrid, get, interpolate, level 838 839 .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 840 @*/ 841 PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 842 { 843 SNES_FAS *fas; 844 PetscErrorCode ierr; 845 SNES levelsnes; 846 847 PetscFunctionBegin; 848 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 849 fas = (SNES_FAS*)levelsnes->data; 850 *mat = fas->interpolate; 851 PetscFunctionReturn(0); 852 } 853 854 /*@ 855 SNESFASSetRestriction - Sets the function to be used to restrict the defect 856 from level l to l-1. 857 858 Input Parameters: 859 + snes - the multigrid context 860 . mat - the restriction matrix 861 - level - the level (0 is coarsest) to supply [Do not supply 0] 862 863 Level: advanced 864 865 Notes: 866 Usually this is the same matrix used also to set the interpolation 867 for the same level. 868 869 One can pass in the interpolation matrix or its transpose; PETSc figures 870 out from the matrix size which one it is. 871 872 If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 873 is used. 874 875 .keywords: FAS, MG, set, multigrid, restriction, level 876 877 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 878 @*/ 879 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 880 { 881 SNES_FAS *fas; 882 PetscErrorCode ierr; 883 SNES levelsnes; 884 885 PetscFunctionBegin; 886 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 887 fas = (SNES_FAS*)levelsnes->data; 888 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 889 ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 890 891 fas->restrct = mat; 892 PetscFunctionReturn(0); 893 } 894 895 /*@ 896 SNESFASGetRestriction - Gets the matrix used to calculate the 897 restriction from l to the l-1th level 898 899 Input Parameters: 900 + snes - the multigrid context 901 - level - the level (0 is coarsest) to supply [do not supply 0] 902 903 Output Parameters: 904 . mat - the interpolation operator 905 906 Level: advanced 907 908 .keywords: FAS, multigrid, get, restrict, level 909 910 .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 911 @*/ 912 PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 913 { 914 SNES_FAS *fas; 915 PetscErrorCode ierr; 916 SNES levelsnes; 917 918 PetscFunctionBegin; 919 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 920 fas = (SNES_FAS*)levelsnes->data; 921 *mat = fas->restrct; 922 PetscFunctionReturn(0); 923 } 924 925 926 /*@ 927 SNESFASSetInjection - Sets the function to be used to inject the solution 928 from level l to l-1. 929 930 Input Parameters: 931 + snes - the multigrid context 932 . mat - the restriction matrix 933 - level - the level (0 is coarsest) to supply [Do not supply 0] 934 935 Level: advanced 936 937 Notes: 938 If you do not set this, the restriction and rscale is used to 939 project the solution instead. 940 941 .keywords: FAS, MG, set, multigrid, restriction, level 942 943 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 944 @*/ 945 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 946 { 947 SNES_FAS *fas; 948 PetscErrorCode ierr; 949 SNES levelsnes; 950 951 PetscFunctionBegin; 952 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 953 fas = (SNES_FAS*)levelsnes->data; 954 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 955 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 956 957 fas->inject = mat; 958 PetscFunctionReturn(0); 959 } 960 961 962 /*@ 963 SNESFASGetInjection - Gets the matrix used to calculate the 964 injection from l-1 to the lth level 965 966 Input Parameters: 967 + snes - the multigrid context 968 - level - the level (0 is coarsest) to supply [do not supply 0] 969 970 Output Parameters: 971 . mat - the injection operator 972 973 Level: advanced 974 975 .keywords: FAS, multigrid, get, restrict, level 976 977 .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 978 @*/ 979 PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 980 { 981 SNES_FAS *fas; 982 PetscErrorCode ierr; 983 SNES levelsnes; 984 985 PetscFunctionBegin; 986 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 987 fas = (SNES_FAS*)levelsnes->data; 988 *mat = fas->inject; 989 PetscFunctionReturn(0); 990 } 991 992 /*@ 993 SNESFASSetRScale - Sets the scaling factor of the restriction 994 operator from level l to l-1. 995 996 Input Parameters: 997 + snes - the multigrid context 998 . rscale - the restriction scaling 999 - level - the level (0 is coarsest) to supply [Do not supply 0] 1000 1001 Level: advanced 1002 1003 Notes: 1004 This is only used in the case that the injection is not set. 1005 1006 .keywords: FAS, MG, set, multigrid, restriction, level 1007 1008 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1009 @*/ 1010 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 1011 { 1012 SNES_FAS *fas; 1013 PetscErrorCode ierr; 1014 SNES levelsnes; 1015 1016 PetscFunctionBegin; 1017 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1018 fas = (SNES_FAS*)levelsnes->data; 1019 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 1020 ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 1021 1022 fas->rscale = rscale; 1023 PetscFunctionReturn(0); 1024 } 1025 1026 /*@ 1027 SNESFASGetSmoother - Gets the default smoother on a level. 1028 1029 Input Parameters: 1030 + snes - the multigrid context 1031 - level - the level (0 is coarsest) to supply 1032 1033 Output Parameters: 1034 smooth - the smoother 1035 1036 Level: advanced 1037 1038 .keywords: FAS, MG, get, multigrid, smoother, level 1039 1040 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1041 @*/ 1042 PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 1043 { 1044 SNES_FAS *fas; 1045 PetscErrorCode ierr; 1046 SNES levelsnes; 1047 1048 PetscFunctionBegin; 1049 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1050 fas = (SNES_FAS*)levelsnes->data; 1051 if (!fas->smoothd) { 1052 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1053 } 1054 *smooth = fas->smoothd; 1055 PetscFunctionReturn(0); 1056 } 1057 1058 /*@ 1059 SNESFASGetSmootherDown - Gets the downsmoother on a level. 1060 1061 Input Parameters: 1062 + snes - the multigrid context 1063 - level - the level (0 is coarsest) to supply 1064 1065 Output Parameters: 1066 smooth - the smoother 1067 1068 Level: advanced 1069 1070 .keywords: FAS, MG, get, multigrid, smoother, level 1071 1072 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1073 @*/ 1074 PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 1075 { 1076 SNES_FAS *fas; 1077 PetscErrorCode ierr; 1078 SNES levelsnes; 1079 1080 PetscFunctionBegin; 1081 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1082 fas = (SNES_FAS*)levelsnes->data; 1083 /* if the user chooses to differentiate smoothers, create them both at this point */ 1084 if (!fas->smoothd) { 1085 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1086 } 1087 if (!fas->smoothu) { 1088 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1089 } 1090 *smooth = fas->smoothd; 1091 PetscFunctionReturn(0); 1092 } 1093 1094 /*@ 1095 SNESFASGetSmootherUp - Gets the upsmoother on a level. 1096 1097 Input Parameters: 1098 + snes - the multigrid context 1099 - level - the level (0 is coarsest) 1100 1101 Output Parameters: 1102 smooth - the smoother 1103 1104 Level: advanced 1105 1106 .keywords: FAS, MG, get, multigrid, smoother, level 1107 1108 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1109 @*/ 1110 PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 1111 { 1112 SNES_FAS *fas; 1113 PetscErrorCode ierr; 1114 SNES levelsnes; 1115 1116 PetscFunctionBegin; 1117 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1118 fas = (SNES_FAS*)levelsnes->data; 1119 /* if the user chooses to differentiate smoothers, create them both at this point */ 1120 if (!fas->smoothd) { 1121 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1122 } 1123 if (!fas->smoothu) { 1124 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1125 } 1126 *smooth = fas->smoothu; 1127 PetscFunctionReturn(0); 1128 } 1129 1130 /*@ 1131 SNESFASGetCoarseSolve - Gets the coarsest solver. 1132 1133 Input Parameters: 1134 . snes - the multigrid context 1135 1136 Output Parameters: 1137 . coarse - the coarse-level solver 1138 1139 Level: advanced 1140 1141 .keywords: FAS, MG, get, multigrid, solver, coarse 1142 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1143 @*/ 1144 PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 1145 { 1146 SNES_FAS *fas; 1147 PetscErrorCode ierr; 1148 SNES levelsnes; 1149 1150 PetscFunctionBegin; 1151 ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1152 fas = (SNES_FAS*)levelsnes->data; 1153 /* if the user chooses to differentiate smoothers, create them both at this point */ 1154 if (!fas->smoothd) { 1155 ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1156 } 1157 *coarse = fas->smoothd; 1158 PetscFunctionReturn(0); 1159 } 1160 1161 /*@ 1162 SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS 1163 1164 Logically Collective on SNES 1165 1166 Input Parameters: 1167 + snes - the multigrid context 1168 - swp - whether to downsweep or not 1169 1170 Options Database Key: 1171 . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1172 1173 Level: advanced 1174 1175 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 1176 1177 .seealso: SNESFASSetNumberSmoothUp() 1178 @*/ 1179 PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp) 1180 { 1181 SNES_FAS *fas = (SNES_FAS*)snes->data; 1182 PetscErrorCode ierr = 0; 1183 1184 PetscFunctionBegin; 1185 fas->full_downsweep = swp; 1186 if (fas->next) { 1187 ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr); 1188 } 1189 PetscFunctionReturn(0); 1190 } 1191