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