1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 3 4 const char *const SNESFASTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "SNESFASType", "SNES_FAS", NULL}; 5 6 static PetscErrorCode SNESReset_FAS(SNES snes) 7 { 8 SNES_FAS *fas = (SNES_FAS *)snes->data; 9 10 PetscFunctionBegin; 11 PetscCall(SNESDestroy(&fas->smoothu)); 12 PetscCall(SNESDestroy(&fas->smoothd)); 13 PetscCall(MatDestroy(&fas->inject)); 14 PetscCall(MatDestroy(&fas->interpolate)); 15 PetscCall(MatDestroy(&fas->restrct)); 16 PetscCall(VecDestroy(&fas->rscale)); 17 PetscCall(VecDestroy(&fas->Xg)); 18 PetscCall(VecDestroy(&fas->Fg)); 19 if (fas->next) PetscCall(SNESReset(fas->next)); 20 PetscFunctionReturn(PETSC_SUCCESS); 21 } 22 23 static PetscErrorCode SNESDestroy_FAS(SNES snes) 24 { 25 SNES_FAS *fas = (SNES_FAS *)snes->data; 26 27 PetscFunctionBegin; 28 /* recursively resets and then destroys */ 29 PetscCall(SNESReset_FAS(snes)); 30 PetscCall(SNESDestroy(&fas->next)); 31 PetscCall(PetscFree(fas)); 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 35 static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth) 36 { 37 SNESLineSearch linesearch; 38 SNESLineSearch slinesearch; 39 void *lsprectx, *lspostctx; 40 PetscErrorCode (*precheck)(SNESLineSearch, Vec, Vec, PetscBool *, void *); 41 PetscErrorCode (*postcheck)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 42 43 PetscFunctionBegin; 44 if (!snes->linesearch) PetscFunctionReturn(PETSC_SUCCESS); 45 PetscCall(SNESGetLineSearch(snes, &linesearch)); 46 PetscCall(SNESGetLineSearch(smooth, &slinesearch)); 47 PetscCall(SNESLineSearchGetPreCheck(linesearch, &precheck, &lsprectx)); 48 PetscCall(SNESLineSearchGetPostCheck(linesearch, &postcheck, &lspostctx)); 49 PetscCall(SNESLineSearchSetPreCheck(slinesearch, precheck, lsprectx)); 50 PetscCall(SNESLineSearchSetPostCheck(slinesearch, postcheck, lspostctx)); 51 PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch)); 52 PetscFunctionReturn(PETSC_SUCCESS); 53 } 54 55 static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth) 56 { 57 SNES_FAS *fas = (SNES_FAS *)snes->data; 58 59 PetscFunctionBegin; 60 PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth)); 61 PetscCall(SNESSetFromOptions(smooth)); 62 PetscCall(SNESFASSetUpLineSearch_Private(snes, smooth)); 63 64 PetscCall(PetscObjectReference((PetscObject)snes->vec_sol)); 65 PetscCall(PetscObjectReference((PetscObject)snes->vec_sol_update)); 66 PetscCall(PetscObjectReference((PetscObject)snes->vec_func)); 67 smooth->vec_sol = snes->vec_sol; 68 smooth->vec_sol_update = snes->vec_sol_update; 69 smooth->vec_func = snes->vec_func; 70 71 if (fas->eventsmoothsetup) PetscCall(PetscLogEventBegin(fas->eventsmoothsetup, smooth, 0, 0, 0)); 72 PetscCall(SNESSetUp(smooth)); 73 if (fas->eventsmoothsetup) PetscCall(PetscLogEventEnd(fas->eventsmoothsetup, smooth, 0, 0, 0)); 74 PetscFunctionReturn(PETSC_SUCCESS); 75 } 76 77 static PetscErrorCode SNESSetUp_FAS(SNES snes) 78 { 79 SNES_FAS *fas = (SNES_FAS *)snes->data; 80 PetscInt dm_levels; 81 SNES next; 82 PetscBool isFine, hasCreateRestriction, hasCreateInjection; 83 84 PetscFunctionBegin; 85 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 86 if (fas->usedmfornumberoflevels && isFine) { 87 PetscCall(DMGetRefineLevel(snes->dm, &dm_levels)); 88 dm_levels++; 89 if (dm_levels > fas->levels) { 90 /* reset the number of levels */ 91 PetscCall(SNESFASSetLevels(snes, dm_levels, NULL)); 92 PetscCall(SNESSetFromOptions(snes)); 93 } 94 } 95 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 96 if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 97 98 PetscCall(SNESSetWorkVecs(snes, 2)); /* work vectors used for intergrid transfers */ 99 100 /* set up the smoothers if they haven't already been set up */ 101 if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 102 103 if (snes->dm) { 104 /* set the smoother DMs properly */ 105 if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm)); 106 PetscCall(SNESSetDM(fas->smoothd, snes->dm)); 107 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 108 if (next) { 109 /* for now -- assume the DM and the evaluation functions have been set externally */ 110 if (!next->dm) { 111 PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm)); 112 PetscCall(SNESSetDM(next, next->dm)); 113 } 114 /* set the interpolation and restriction from the DM */ 115 if (!fas->interpolate) { 116 PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale)); 117 if (!fas->restrct) { 118 PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction)); 119 /* DM can create restrictions, use that */ 120 if (hasCreateRestriction) { 121 PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct)); 122 } else { 123 PetscCall(PetscObjectReference((PetscObject)fas->interpolate)); 124 fas->restrct = fas->interpolate; 125 } 126 } 127 } 128 /* set the injection from the DM */ 129 if (!fas->inject) { 130 PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection)); 131 if (hasCreateInjection) PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject)); 132 } 133 } 134 } 135 136 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 137 if (fas->galerkin) { 138 if (next) PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next)); 139 if (fas->smoothd && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes)); 140 if (fas->smoothu && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes)); 141 } 142 143 /* sets the down (pre) smoother's default norm and sets it from options */ 144 if (fas->smoothd) { 145 if (fas->level == 0) { 146 PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_ALWAYS)); 147 } else { 148 PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY)); 149 } 150 PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd)); 151 } 152 153 /* sets the up (post) smoother's default norm and sets it from options */ 154 if (fas->smoothu) { 155 if (fas->level != fas->levels - 1) { 156 PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE)); 157 } else { 158 PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY)); 159 } 160 PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu)); 161 } 162 163 if (next) { 164 /* gotta set up the solution vector for this to work */ 165 if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_sol)); 166 if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_rhs)); 167 PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next)); 168 PetscCall(SNESFASSetUpLineSearch_Private(snes, next)); 169 PetscCall(SNESSetUp(next)); 170 } 171 172 /* setup FAS work vectors */ 173 if (fas->galerkin) { 174 PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg)); 175 PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg)); 176 } 177 PetscFunctionReturn(PETSC_SUCCESS); 178 } 179 180 static PetscErrorCode SNESSetFromOptions_FAS(SNES snes, PetscOptionItems PetscOptionsObject) 181 { 182 SNES_FAS *fas = (SNES_FAS *)snes->data; 183 PetscInt levels = 1; 184 PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE, continuationflg = PETSC_FALSE; 185 SNESFASType fastype; 186 const char *optionsprefix; 187 SNESLineSearch linesearch; 188 PetscInt m, n_up, n_down; 189 SNES next; 190 PetscBool isFine; 191 192 PetscFunctionBegin; 193 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 194 PetscOptionsHeadBegin(PetscOptionsObject, "SNESFAS Options-----------------------------------"); 195 196 /* number of levels -- only process most options on the finest level */ 197 if (isFine) { 198 PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg)); 199 if (!flg && snes->dm) { 200 PetscCall(DMGetRefineLevel(snes->dm, &levels)); 201 levels++; 202 fas->usedmfornumberoflevels = PETSC_TRUE; 203 } 204 PetscCall(SNESFASSetLevels(snes, levels, NULL)); 205 fastype = fas->fastype; 206 PetscCall(PetscOptionsEnum("-snes_fas_type", "FAS correction type", "SNESFASSetType", SNESFASTypes, (PetscEnum)fastype, (PetscEnum *)&fastype, &flg)); 207 if (flg) PetscCall(SNESFASSetType(snes, fastype)); 208 209 PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 210 PetscCall(PetscOptionsInt("-snes_fas_cycles", "Number of cycles", "SNESFASSetCycles", fas->n_cycles, &m, &flg)); 211 if (flg) PetscCall(SNESFASSetCycles(snes, m)); 212 PetscCall(PetscOptionsBool("-snes_fas_continuation", "Corrected grid-sequence continuation", "SNESFASSetContinuation", fas->continuation, &continuationflg, &flg)); 213 if (flg) PetscCall(SNESFASSetContinuation(snes, continuationflg)); 214 215 PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin", "SNESFASSetGalerkin", fas->galerkin, &galerkinflg, &flg)); 216 if (flg) PetscCall(SNESFASSetGalerkin(snes, galerkinflg)); 217 218 if (fas->fastype == SNES_FAS_FULL) { 219 PetscCall(PetscOptionsBool("-snes_fas_full_downsweep", "Smooth on the initial down sweep for full FAS cycles", "SNESFASFullSetDownSweep", fas->full_downsweep, &fas->full_downsweep, &flg)); 220 if (flg) PetscCall(SNESFASFullSetDownSweep(snes, fas->full_downsweep)); 221 PetscCall(PetscOptionsBool("-snes_fas_full_total", "Use total restriction and interpolaton on the indial down and up sweeps for the full FAS cycle", "SNESFASFullSetUseTotal", fas->full_total, &fas->full_total, &flg)); 222 if (flg) PetscCall(SNESFASFullSetTotal(snes, fas->full_total)); 223 } 224 225 PetscCall(PetscOptionsInt("-snes_fas_smoothup", "Number of post-smoothing steps", "SNESFASSetNumberSmoothUp", fas->max_up_it, &n_up, &upflg)); 226 227 PetscCall(PetscOptionsInt("-snes_fas_smoothdown", "Number of pre-smoothing steps", "SNESFASSetNumberSmoothDown", fas->max_down_it, &n_down, &downflg)); 228 229 { 230 PetscViewer viewer; 231 PetscViewerFormat format; 232 PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_fas_monitor", &viewer, &format, &monflg)); 233 if (monflg) { 234 PetscViewerAndFormat *vf; 235 PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 236 PetscCall(PetscViewerDestroy(&viewer)); 237 PetscCall(SNESFASSetMonitor(snes, vf, PETSC_TRUE)); 238 } 239 } 240 flg = PETSC_FALSE; 241 monflg = PETSC_TRUE; 242 PetscCall(PetscOptionsBool("-snes_fas_log", "Log times for each FAS level", "SNESFASSetLog", monflg, &monflg, &flg)); 243 if (flg) PetscCall(SNESFASSetLog(snes, monflg)); 244 } 245 246 PetscOptionsHeadEnd(); 247 248 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 249 if (upflg) PetscCall(SNESFASSetNumberSmoothUp(snes, n_up)); 250 if (downflg) PetscCall(SNESFASSetNumberSmoothDown(snes, n_down)); 251 252 /* set up the default line search for coarse grid corrections */ 253 if (fas->fastype == SNES_FAS_ADDITIVE) { 254 if (!snes->linesearch) { 255 PetscCall(SNESGetLineSearch(snes, &linesearch)); 256 PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT)); 257 } 258 } 259 260 /* recursive option setting for the smoothers */ 261 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 262 if (next) PetscCall(SNESSetFromOptions(next)); 263 PetscFunctionReturn(PETSC_SUCCESS); 264 } 265 266 #include <petscdraw.h> 267 static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 268 { 269 SNES_FAS *fas = (SNES_FAS *)snes->data; 270 PetscBool isFine, isascii, isdraw; 271 PetscInt i; 272 SNES smoothu, smoothd, levelsnes; 273 274 PetscFunctionBegin; 275 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 276 if (isFine) { 277 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 278 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 279 if (isascii) { 280 PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles)); 281 if (fas->galerkin) { 282 PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid function evaluation\n")); 283 } else { 284 PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid function evaluation\n")); 285 } 286 for (i = 0; i < fas->levels; i++) { 287 PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 288 PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu)); 289 PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd)); 290 if (!i) { 291 PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i)); 292 } else { 293 PetscCall(PetscViewerASCIIPrintf(viewer, " Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 294 } 295 PetscCall(PetscViewerASCIIPushTab(viewer)); 296 if (smoothd) { 297 PetscCall(SNESView(smoothd, viewer)); 298 } else { 299 PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n")); 300 } 301 PetscCall(PetscViewerASCIIPopTab(viewer)); 302 if (i && (smoothd == smoothu)) { 303 PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) same as down solver (pre-smoother)\n")); 304 } else if (i) { 305 PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 306 PetscCall(PetscViewerASCIIPushTab(viewer)); 307 if (smoothu) { 308 PetscCall(SNESView(smoothu, viewer)); 309 } else { 310 PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n")); 311 } 312 PetscCall(PetscViewerASCIIPopTab(viewer)); 313 } 314 } 315 } else if (isdraw) { 316 PetscDraw draw; 317 PetscReal x, w, y, bottom, th, wth; 318 SNES_FAS *curfas = fas; 319 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 320 PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 321 PetscCall(PetscDrawStringGetSize(draw, &wth, &th)); 322 bottom = y - th; 323 while (curfas) { 324 if (!curfas->smoothu) { 325 PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 326 if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer)); 327 PetscCall(PetscDrawPopCurrentPoint(draw)); 328 } else { 329 w = 0.5 * PetscMin(1.0 - x, x); 330 PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom)); 331 if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer)); 332 PetscCall(PetscDrawPopCurrentPoint(draw)); 333 PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom)); 334 if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu, viewer)); 335 PetscCall(PetscDrawPopCurrentPoint(draw)); 336 } 337 /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 338 bottom -= 5 * th; 339 if (curfas->next) curfas = (SNES_FAS *)curfas->next->data; 340 else curfas = NULL; 341 } 342 } 343 } 344 PetscFunctionReturn(PETSC_SUCCESS); 345 } 346 347 /* 348 Defines the action of the downsmoother 349 */ 350 static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 351 { 352 SNESConvergedReason reason; 353 Vec FPC; 354 SNES smoothd; 355 PetscBool flg; 356 SNES_FAS *fas = (SNES_FAS *)snes->data; 357 358 PetscFunctionBegin; 359 PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd)); 360 PetscCall(SNESSetInitialFunction(smoothd, F)); 361 if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothd, B, X, 0)); 362 PetscCall(SNESSolve(smoothd, B, X)); 363 if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothd, B, X, 0)); 364 /* check convergence reason for the smoother */ 365 PetscCall(SNESGetConvergedReason(smoothd, &reason)); 366 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 367 snes->reason = SNES_DIVERGED_INNER; 368 PetscFunctionReturn(PETSC_SUCCESS); 369 } 370 371 PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL)); 372 PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg)); 373 if (!flg) PetscCall(SNESComputeFunction(smoothd, X, FPC)); 374 PetscCall(VecCopy(FPC, F)); 375 if (fnorm) { 376 PetscCall(VecNorm(F, NORM_2, fnorm)); 377 SNESCheckFunctionDomainError(snes, *fnorm); 378 } 379 PetscFunctionReturn(PETSC_SUCCESS); 380 } 381 382 /* 383 Defines the action of the upsmoother 384 */ 385 static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 386 { 387 SNESConvergedReason reason; 388 Vec FPC; 389 SNES smoothu; 390 PetscBool flg; 391 SNES_FAS *fas = (SNES_FAS *)snes->data; 392 393 PetscFunctionBegin; 394 PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu)); 395 if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothu, 0, 0, 0)); 396 PetscCall(SNESSolve(smoothu, B, X)); 397 if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothu, 0, 0, 0)); 398 /* check convergence reason for the smoother */ 399 PetscCall(SNESGetConvergedReason(smoothu, &reason)); 400 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 401 snes->reason = SNES_DIVERGED_INNER; 402 PetscFunctionReturn(PETSC_SUCCESS); 403 } 404 PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL)); 405 PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg)); 406 if (!flg) PetscCall(SNESComputeFunction(smoothu, X, FPC)); 407 PetscCall(VecCopy(FPC, F)); 408 if (fnorm) { 409 PetscCall(VecNorm(F, NORM_2, fnorm)); 410 SNESCheckFunctionDomainError(snes, *fnorm); 411 } 412 PetscFunctionReturn(PETSC_SUCCESS); 413 } 414 415 /*@ 416 SNESFASCreateCoarseVec - create a `Vec` corresponding to a state vector on one level coarser than the current level 417 418 Collective 419 420 Input Parameter: 421 . snes - `SNESFAS` object 422 423 Output Parameter: 424 . Xcoarse - vector on level one coarser than the current level 425 426 Level: developer 427 428 .seealso: [](ch_snes), `SNESFASSetRestriction()`, `SNESFASRestrict()`, `SNESFAS` 429 @*/ 430 PetscErrorCode SNESFASCreateCoarseVec(SNES snes, Vec *Xcoarse) 431 { 432 SNES_FAS *fas; 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 436 PetscAssertPointer(Xcoarse, 2); 437 fas = (SNES_FAS *)snes->data; 438 if (fas->rscale) { 439 PetscCall(VecDuplicate(fas->rscale, Xcoarse)); 440 } else if (fas->interpolate) { 441 PetscCall(MatCreateVecs(fas->interpolate, Xcoarse, NULL)); 442 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set rscale or interpolation"); 443 PetscFunctionReturn(PETSC_SUCCESS); 444 } 445 446 /*@ 447 SNESFASRestrict - restrict a `Vec` to the next coarser level 448 449 Collective 450 451 Input Parameters: 452 + fine - `SNES` from which to restrict 453 - Xfine - vector to restrict 454 455 Output Parameter: 456 . Xcoarse - result of restriction 457 458 Level: developer 459 460 .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`, `SNESFASCreateCoarseVec()` 461 @*/ 462 PetscErrorCode SNESFASRestrict(SNES fine, Vec Xfine, Vec Xcoarse) 463 { 464 SNES_FAS *fas; 465 466 PetscFunctionBegin; 467 PetscValidHeaderSpecificType(fine, SNES_CLASSID, 1, SNESFAS); 468 PetscValidHeaderSpecific(Xfine, VEC_CLASSID, 2); 469 PetscValidHeaderSpecific(Xcoarse, VEC_CLASSID, 3); 470 fas = (SNES_FAS *)fine->data; 471 if (fas->inject) { 472 PetscCall(MatRestrict(fas->inject, Xfine, Xcoarse)); 473 } else { 474 PetscCall(MatRestrict(fas->restrct, Xfine, Xcoarse)); 475 PetscCall(VecPointwiseMult(Xcoarse, fas->rscale, Xcoarse)); 476 } 477 PetscFunctionReturn(PETSC_SUCCESS); 478 } 479 480 /* 481 Performs a variant of FAS using the interpolated total coarse solution 482 483 fine problem: F(x) = b 484 coarse problem: F^c(x^c) = Rb, Initial guess Rx 485 interpolated solution: x^f = I x^c (total solution interpolation 486 */ 487 static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new) 488 { 489 Vec X_c, B_c; 490 SNESConvergedReason reason; 491 SNES next; 492 Mat restrct, interpolate; 493 SNES_FAS *fasc; 494 495 PetscFunctionBegin; 496 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 497 if (next) { 498 fasc = (SNES_FAS *)next->data; 499 500 PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 501 PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 502 503 X_c = next->vec_sol; 504 505 if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 506 /* restrict the total solution: Rb */ 507 PetscCall(SNESFASRestrict(snes, X, X_c)); 508 B_c = next->vec_rhs; 509 if (snes->vec_rhs) { 510 /* restrict the total rhs defect: Rb */ 511 PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c)); 512 } else { 513 PetscCall(VecSet(B_c, 0.)); 514 } 515 if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 516 517 PetscCall(SNESSolve(next, B_c, X_c)); 518 PetscCall(SNESGetConvergedReason(next, &reason)); 519 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 520 snes->reason = SNES_DIVERGED_INNER; 521 PetscFunctionReturn(PETSC_SUCCESS); 522 } 523 /* x^f <- Ix^c*/ 524 DM dmc, dmf; 525 526 PetscCall(SNESGetDM(next, &dmc)); 527 PetscCall(SNESGetDM(snes, &dmf)); 528 if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 529 PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new)); 530 if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 531 PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse solution")); 532 PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 533 PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution")); 534 PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 535 } 536 PetscFunctionReturn(PETSC_SUCCESS); 537 } 538 539 /* 540 Performs the FAS coarse correction as: 541 542 fine problem: F(x) = b 543 coarse problem: F^c(x^c) = b^c 544 545 b^c = F^c(Rx) - R(F(x) - b) 546 */ 547 static PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 548 { 549 Vec X_c, Xo_c, F_c, B_c; 550 SNESConvergedReason reason; 551 SNES next; 552 Mat restrct, interpolate; 553 SNES_FAS *fasc; 554 555 PetscFunctionBegin; 556 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 557 if (next) { 558 fasc = (SNES_FAS *)next->data; 559 560 PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 561 PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 562 563 X_c = next->vec_sol; 564 Xo_c = next->work[0]; 565 F_c = next->vec_func; 566 B_c = next->vec_rhs; 567 568 if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 569 PetscCall(SNESFASRestrict(snes, X, Xo_c)); 570 /* restrict the defect: R(F(x) - b) */ 571 PetscCall(MatRestrict(restrct, F, B_c)); 572 if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 573 574 if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0)); 575 /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 576 PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 577 if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0)); 578 579 /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 580 PetscCall(VecCopy(B_c, X_c)); 581 PetscCall(VecCopy(F_c, B_c)); 582 PetscCall(VecCopy(X_c, F_c)); 583 /* set initial guess of the coarse problem to the projected fine solution */ 584 PetscCall(VecCopy(Xo_c, X_c)); 585 586 /* recurse to the next level */ 587 PetscCall(SNESSetInitialFunction(next, F_c)); 588 PetscCall(SNESSolve(next, B_c, X_c)); 589 PetscCall(SNESGetConvergedReason(next, &reason)); 590 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 591 snes->reason = SNES_DIVERGED_INNER; 592 PetscFunctionReturn(PETSC_SUCCESS); 593 } 594 /* correct as x <- x + I(x^c - Rx)*/ 595 PetscCall(VecAXPY(X_c, -1.0, Xo_c)); 596 597 if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 598 PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new)); 599 if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 600 PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse correction")); 601 PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 602 PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution")); 603 PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 604 } 605 PetscFunctionReturn(PETSC_SUCCESS); 606 } 607 608 /* 609 The additive cycle is: 610 611 xhat = x 612 xhat = dS(x, b) 613 x = coarsecorrection(xhat, b_d) 614 x = x + nu*(xhat - x); 615 (optional) x = uS(x, b) 616 617 With the coarse RHS (defect correction) as below. 618 */ 619 static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 620 { 621 Vec F, B, Xhat; 622 Vec X_c, Xo_c, F_c, B_c; 623 SNESConvergedReason reason; 624 PetscReal xnorm, fnorm, ynorm; 625 SNES next; 626 Mat restrct, interpolate; 627 SNES_FAS *fas = (SNES_FAS *)snes->data, *fasc; 628 629 PetscFunctionBegin; 630 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 631 F = snes->vec_func; 632 B = snes->vec_rhs; 633 Xhat = snes->work[1]; 634 PetscCall(VecCopy(X, Xhat)); 635 /* recurse first */ 636 if (next) { 637 fasc = (SNES_FAS *)next->data; 638 PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 639 PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 640 if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0)); 641 PetscCall(SNESComputeFunction(snes, Xhat, F)); 642 if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0)); 643 PetscCall(VecNorm(F, NORM_2, &fnorm)); 644 SNESCheckFunctionDomainError(snes, fnorm); 645 X_c = next->vec_sol; 646 Xo_c = next->work[0]; 647 F_c = next->vec_func; 648 B_c = next->vec_rhs; 649 650 PetscCall(SNESFASRestrict(snes, Xhat, Xo_c)); 651 /* restrict the defect */ 652 PetscCall(MatRestrict(restrct, F, B_c)); 653 654 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 655 if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0)); 656 PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 657 if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0)); 658 PetscCall(VecCopy(B_c, X_c)); 659 PetscCall(VecCopy(F_c, B_c)); 660 PetscCall(VecCopy(X_c, F_c)); 661 /* set initial guess of the coarse problem to the projected fine solution */ 662 PetscCall(VecCopy(Xo_c, X_c)); 663 664 /* recurse */ 665 PetscCall(SNESSetInitialFunction(next, F_c)); 666 PetscCall(SNESSolve(next, B_c, X_c)); 667 668 /* smooth on this level */ 669 PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm)); 670 671 PetscCall(SNESGetConvergedReason(next, &reason)); 672 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 673 snes->reason = SNES_DIVERGED_INNER; 674 PetscFunctionReturn(PETSC_SUCCESS); 675 } 676 677 /* correct as x <- x + I(x^c - Rx)*/ 678 PetscCall(VecAYPX(X_c, -1.0, Xo_c)); 679 PetscCall(MatInterpolate(interpolate, X_c, Xhat)); 680 681 /* additive correction of the coarse direction*/ 682 PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat)); 683 PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm)); 684 SNESCheckLineSearchFailure(snes); 685 } else { 686 PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 687 } 688 PetscFunctionReturn(PETSC_SUCCESS); 689 } 690 691 /* 692 Defines the FAS cycle as: 693 694 fine problem: F(x) = b 695 coarse problem: F^c(x) = b^c 696 697 b^c = F^c(Rx) - R(F(x) - b) 698 699 correction: 700 701 x = x + I(x^c - Rx) 702 */ 703 static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 704 { 705 Vec F, B; 706 SNES next; 707 708 PetscFunctionBegin; 709 F = snes->vec_func; 710 B = snes->vec_rhs; 711 /* pre-smooth -- just update using the pre-smoother */ 712 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 713 PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 714 if (next) { 715 PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 716 PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 717 } 718 PetscFunctionReturn(PETSC_SUCCESS); 719 } 720 721 static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 722 { 723 SNES next; 724 SNES_FAS *fas = (SNES_FAS *)snes->data; 725 PetscBool isFine; 726 727 PetscFunctionBegin; 728 /* pre-smooth -- just update using the pre-smoother */ 729 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 730 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 731 fas->full_stage = 0; 732 if (next) PetscCall(SNESFASCycleSetupPhase_Full(next)); 733 PetscFunctionReturn(PETSC_SUCCESS); 734 } 735 736 static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 737 { 738 Vec F, B; 739 SNES_FAS *fas = (SNES_FAS *)snes->data; 740 PetscBool isFine; 741 SNES next; 742 743 PetscFunctionBegin; 744 F = snes->vec_func; 745 B = snes->vec_rhs; 746 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 747 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 748 749 if (isFine) PetscCall(SNESFASCycleSetupPhase_Full(snes)); 750 751 if (fas->full_stage == 0) { 752 /* downsweep */ 753 if (next) { 754 if (fas->level != 1) next->max_its += 1; 755 if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 756 fas->full_downsweep = PETSC_TRUE; 757 if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes, X, X)); 758 else PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 759 fas->full_total = PETSC_FALSE; 760 PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 761 if (fas->level != 1) next->max_its -= 1; 762 } else { 763 /* The smoother on the coarse level is the coarse solver */ 764 PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 765 } 766 fas->full_stage = 1; 767 } else if (fas->full_stage == 1) { 768 if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 769 if (next) { 770 PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 771 PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 772 } 773 } 774 /* final v-cycle */ 775 if (isFine) { 776 if (next) { 777 PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 778 PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 779 } 780 } 781 PetscFunctionReturn(PETSC_SUCCESS); 782 } 783 784 static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 785 { 786 Vec F, B; 787 SNES next; 788 789 PetscFunctionBegin; 790 F = snes->vec_func; 791 B = snes->vec_rhs; 792 PetscCall(SNESFASCycleGetCorrection(snes, &next)); 793 if (next) { 794 PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 795 PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 796 } else { 797 PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 798 } 799 PetscFunctionReturn(PETSC_SUCCESS); 800 } 801 802 PetscBool SNEScite = PETSC_FALSE; 803 const char SNESCitation[] = "@Article{bruneknepleysmithtu15," 804 " title = {Composing Scalable Nonlinear Algebraic Solvers}," 805 " author = {Peter R. Brune and Matthew G. Knepley and Barry F. Smith and Xuemin Tu}," 806 " journal = {SIAM Review}," 807 " volume = {57}," 808 " number = {4}," 809 " pages = {535--565}," 810 " doi = {10.1137/130936725}," 811 " year = {2015}\n"; 812 813 static PetscErrorCode SNESSolve_FAS(SNES snes) 814 { 815 PetscInt i; 816 Vec X, F; 817 PetscReal fnorm; 818 SNES_FAS *fas = (SNES_FAS *)snes->data, *ffas; 819 DM dm; 820 PetscBool isFine; 821 822 PetscFunctionBegin; 823 PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 824 825 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 826 snes->reason = SNES_CONVERGED_ITERATING; 827 X = snes->vec_sol; 828 F = snes->vec_func; 829 830 PetscCall(SNESFASCycleIsFine(snes, &isFine)); 831 /* norm setup */ 832 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 833 snes->iter = 0; 834 snes->norm = 0; 835 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 836 if (!snes->vec_func_init_set) { 837 if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0)); 838 PetscCall(SNESComputeFunction(snes, X, F)); 839 if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0)); 840 } else snes->vec_func_init_set = PETSC_FALSE; 841 842 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 843 SNESCheckFunctionDomainError(snes, fnorm); 844 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 845 snes->norm = fnorm; 846 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 847 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 848 849 /* test convergence */ 850 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 851 PetscCall(SNESMonitor(snes, snes->iter, fnorm)); 852 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 853 854 if (isFine) { 855 /* propagate scale-dependent data up the hierarchy */ 856 PetscCall(SNESGetDM(snes, &dm)); 857 for (ffas = fas; ffas->next; ffas = (SNES_FAS *)ffas->next->data) { 858 DM dmcoarse; 859 PetscCall(SNESGetDM(ffas->next, &dmcoarse)); 860 PetscCall(DMRestrict(dm, ffas->restrct, ffas->rscale, ffas->inject, dmcoarse)); 861 dm = dmcoarse; 862 } 863 } 864 865 for (i = 0; i < snes->max_its; i++) { 866 /* Call general purpose update function */ 867 PetscTryTypeMethod(snes, update, snes->iter); 868 869 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 870 PetscCall(SNESFASCycle_Multiplicative(snes, X)); 871 } else if (fas->fastype == SNES_FAS_ADDITIVE) { 872 PetscCall(SNESFASCycle_Additive(snes, X)); 873 } else if (fas->fastype == SNES_FAS_FULL) { 874 PetscCall(SNESFASCycle_Full(snes, X)); 875 } else if (fas->fastype == SNES_FAS_KASKADE) { 876 PetscCall(SNESFASCycle_Kaskade(snes, X)); 877 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported FAS type"); 878 879 /* check for FAS cycle divergence */ 880 if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 881 882 /* Monitor convergence */ 883 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 884 snes->iter = i + 1; 885 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 886 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 887 PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, snes->norm)); 888 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 889 if (snes->reason) break; 890 } 891 PetscFunctionReturn(PETSC_SUCCESS); 892 } 893 894 /*MC 895 SNESFAS - An implementation of the Full Approximation Scheme nonlinear multigrid solver, FAS, or nonlinear multigrid {cite}`bruneknepleysmithtu15` for 896 solving nonlinear systems of equations with `SNES`. 897 898 The nonlinear problem is solved by correction using coarse versions 899 of the nonlinear problem. This problem is perturbed so that a projected 900 solution of the fine problem elicits no correction from the coarse problem. 901 902 Options Database Keys and Prefixes: 903 + -snes_fas_levels <l> - The number of levels 904 . -snes_fas_cycles <c> - The number of cycles -- 1 for V, 2 for W 905 . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 906 . -snes_fas_galerkin <false,true> - Form coarse problems by projection back upon the fine problem 907 . -snes_fas_smoothup <u> - The number of iterations of the post-smoother 908 . -snes_fas_smoothdown <d> - The number of iterations of the pre-smoother 909 . -snes_fas_monitor - Monitor progress of all of the levels 910 . -snes_fas_full_downsweep <false,true> - call the downsmooth on the initial downsweep of full FAS 911 . -fas_levels_snes_ - prefix for `SNES` options for all smoothers 912 . -fas_levels_cycle_snes_ - prefix for `SNES` options for all cycles 913 . -fas_levels_i_snes_ - prefix `SNES` options for the smoothers on level i 914 . -fas_levels_i_cycle_snes_ - prefix for `SNES` options for the cycle on level i 915 - -fas_coarse_snes_ - prefix for `SNES` options for the coarsest smoother 916 917 Level: beginner 918 919 Note: 920 The organization of the `SNESFAS` solver is slightly different from the organization of `PCMG` 921 As each level has smoother `SNES` instances(down and potentially up) and a cycle `SNES` instance. 922 The cycle `SNES` instance may be used for monitoring convergence on a particular level. 923 924 .seealso: [](ch_snes), `PCMG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`, 925 `SNESFASFullGetTotal()`, `SNESFASSetType()`, `SNESFASGetType()`, `SNESFASSetLevels()`, `SNESFASGetLevels()`, `SNESFASGetCycleSNES()`, 926 `SNESFASSetNumberSmoothUp()`, `SNESFASSetNumberSmoothDown()`, `SNESFASSetContinuation()`, `SNESFASSetCycles()`, `SNESFASSetMonitor()`, 927 `SNESFASSetLog()`, `SNESFASCycleSetCycles()`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`, 928 `SNESFASCycleGetCorrection()`, `SNESFASCycleGetInterpolation()`, `SNESFASCycleGetRestriction()`, `SNESFASCycleGetInjection()`, 929 `SNESFASCycleGetRScale()`, `SNESFASCycleIsFine()`, `SNESFASSetInterpolation()`, `SNESFASGetInterpolation()`, `SNESFASSetRestriction()`, 930 `SNESFASGetRestriction()`, `SNESFASSetInjection()`, `SNESFASGetInjection()`, `SNESFASSetRScale()`,`SNESFASGetSmoother()`, 931 `SNESFASGetSmootherDown()`, `SNESFASGetSmootherUp()`, `SNESFASGetCoarseSolve()`, `SNESFASFullSetDownSweep()`, `SNESFASFullSetTotal()`, `SNESFASFullGetTotal()` 932 M*/ 933 934 PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 935 { 936 SNES_FAS *fas; 937 938 PetscFunctionBegin; 939 snes->ops->destroy = SNESDestroy_FAS; 940 snes->ops->setup = SNESSetUp_FAS; 941 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 942 snes->ops->view = SNESView_FAS; 943 snes->ops->solve = SNESSolve_FAS; 944 snes->ops->reset = SNESReset_FAS; 945 946 snes->usesksp = PETSC_FALSE; 947 snes->usesnpc = PETSC_FALSE; 948 949 PetscCall(SNESParametersInitialize(snes)); 950 PetscObjectParameterSetDefault(snes, max_funcs, 30000); 951 PetscObjectParameterSetDefault(snes, max_its, 10000); 952 953 snes->alwayscomputesfinalresidual = PETSC_TRUE; 954 955 PetscCall(PetscNew(&fas)); 956 957 snes->data = (void *)fas; 958 fas->level = 0; 959 fas->levels = 1; 960 fas->n_cycles = 1; 961 fas->max_up_it = 1; 962 fas->max_down_it = 1; 963 fas->smoothu = NULL; 964 fas->smoothd = NULL; 965 fas->next = NULL; 966 fas->previous = NULL; 967 fas->fine = snes; 968 fas->interpolate = NULL; 969 fas->restrct = NULL; 970 fas->inject = NULL; 971 fas->usedmfornumberoflevels = PETSC_FALSE; 972 fas->fastype = SNES_FAS_MULTIPLICATIVE; 973 fas->full_downsweep = PETSC_FALSE; 974 fas->full_total = PETSC_FALSE; 975 976 fas->eventsmoothsetup = 0; 977 fas->eventsmoothsolve = 0; 978 fas->eventresidual = 0; 979 fas->eventinterprestrict = 0; 980 PetscFunctionReturn(PETSC_SUCCESS); 981 } 982