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