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