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