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