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