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