1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3 4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 5 6 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 7 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 8 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 9 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 10 extern PetscErrorCode SNESSolve_FAS(SNES snes); 11 extern PetscErrorCode SNESReset_FAS(SNES snes); 12 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 13 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *); 14 15 /*MC 16 17 SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 18 19 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 20 21 Level: advanced 22 23 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 24 M*/ 25 26 #undef __FUNCT__ 27 #define __FUNCT__ "SNESCreate_FAS" 28 PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes) 29 { 30 SNES_FAS * fas; 31 PetscErrorCode ierr; 32 33 PetscFunctionBegin; 34 snes->ops->destroy = SNESDestroy_FAS; 35 snes->ops->setup = SNESSetUp_FAS; 36 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 37 snes->ops->view = SNESView_FAS; 38 snes->ops->solve = SNESSolve_FAS; 39 snes->ops->reset = SNESReset_FAS; 40 41 snes->usesksp = PETSC_FALSE; 42 snes->usespc = PETSC_FALSE; 43 44 if (!snes->tolerancesset) { 45 snes->max_funcs = 30000; 46 snes->max_its = 10000; 47 } 48 49 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 50 snes->data = (void*) fas; 51 fas->level = 0; 52 fas->levels = 1; 53 fas->n_cycles = 1; 54 fas->max_up_it = 1; 55 fas->max_down_it = 1; 56 fas->smoothu = PETSC_NULL; 57 fas->smoothd = PETSC_NULL; 58 fas->next = PETSC_NULL; 59 fas->previous = PETSC_NULL; 60 fas->fine = snes; 61 fas->interpolate = PETSC_NULL; 62 fas->restrct = PETSC_NULL; 63 fas->inject = PETSC_NULL; 64 fas->monitor = PETSC_NULL; 65 fas->usedmfornumberoflevels = PETSC_FALSE; 66 fas->fastype = SNES_FAS_MULTIPLICATIVE; 67 PetscFunctionReturn(0); 68 } 69 70 #undef __FUNCT__ 71 #define __FUNCT__ "SNESReset_FAS" 72 PetscErrorCode SNESReset_FAS(SNES snes) 73 { 74 PetscErrorCode ierr = 0; 75 SNES_FAS * fas = (SNES_FAS *)snes->data; 76 77 PetscFunctionBegin; 78 if (fas->smoothu != fas->smoothd) { 79 ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 80 ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 81 } else { 82 ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 83 fas->smoothu = PETSC_NULL; 84 } 85 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 86 ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 87 ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 88 ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 89 if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 90 91 PetscFunctionReturn(0); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "SNESDestroy_FAS" 96 PetscErrorCode SNESDestroy_FAS(SNES snes) 97 { 98 SNES_FAS * fas = (SNES_FAS *)snes->data; 99 PetscErrorCode ierr = 0; 100 101 PetscFunctionBegin; 102 /* recursively resets and then destroys */ 103 ierr = SNESReset(snes);CHKERRQ(ierr); 104 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 105 ierr = PetscFree(fas);CHKERRQ(ierr); 106 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "SNESSetUp_FAS" 112 PetscErrorCode SNESSetUp_FAS(SNES snes) 113 { 114 SNES_FAS *fas = (SNES_FAS *) snes->data; 115 PetscErrorCode ierr; 116 VecScatter injscatter; 117 PetscInt dm_levels; 118 Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 119 SNES next; 120 PetscBool isFine; 121 PetscFunctionBegin; 122 123 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 124 if (fas->usedmfornumberoflevels && isFine) { 125 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 126 dm_levels++; 127 if (dm_levels > fas->levels) { 128 /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 129 vec_sol = snes->vec_sol; 130 vec_func = snes->vec_func; 131 vec_sol_update = snes->vec_sol_update; 132 vec_rhs = snes->vec_rhs; 133 snes->vec_sol = PETSC_NULL; 134 snes->vec_func = PETSC_NULL; 135 snes->vec_sol_update = PETSC_NULL; 136 snes->vec_rhs = PETSC_NULL; 137 138 /* reset the number of levels */ 139 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 140 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 141 142 snes->vec_sol = vec_sol; 143 snes->vec_func = vec_func; 144 snes->vec_rhs = vec_rhs; 145 snes->vec_sol_update = vec_sol_update; 146 } 147 } 148 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 149 if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 150 151 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 152 ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 153 } else { 154 ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 155 } 156 157 /* set up the smoothers if they haven't already been set up */ 158 if (!fas->smoothd) { 159 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 160 } 161 162 if (snes->dm) { 163 /* set the smoother DMs properly */ 164 if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 165 ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 166 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 167 if (next) { 168 /* for now -- assume the DM and the evaluation functions have been set externally */ 169 if (!next->dm) { 170 ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr); 171 ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 172 } 173 /* set the interpolation and restriction from the DM */ 174 if (!fas->interpolate) { 175 ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 176 if (!fas->restrct) { 177 ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 178 fas->restrct = fas->interpolate; 179 } 180 } 181 /* set the injection from the DM */ 182 if (!fas->inject) { 183 ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 184 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 185 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 186 } 187 } 188 } 189 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 190 191 if (fas->galerkin) { 192 if (next) 193 ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 194 if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 195 if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 196 } 197 198 /* set the smoothers up here so that precedence is taken for instance-specific options over the whole-solver options */ 199 if(fas->smoothu) ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 200 if(fas->smoothd) ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 201 202 if (next) { 203 /* gotta set up the solution vector for this to work */ 204 if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 205 if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 206 ierr = SNESSetUp(next);CHKERRQ(ierr); 207 } 208 /* setup FAS work vectors */ 209 if (fas->galerkin) { 210 ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 211 ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 212 } 213 PetscFunctionReturn(0); 214 } 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "SNESSetFromOptions_FAS" 218 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 219 { 220 SNES_FAS *fas = (SNES_FAS *) snes->data; 221 PetscInt levels = 1; 222 PetscBool flg, upflg, downflg, monflg, galerkinflg; 223 PetscErrorCode ierr; 224 char monfilename[PETSC_MAX_PATH_LEN]; 225 SNESFASType fastype; 226 const char *optionsprefix; 227 SNESLineSearch linesearch; 228 PetscInt m, n_up, n_down; 229 SNES next; 230 PetscBool isFine; 231 232 PetscFunctionBegin; 233 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 234 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 235 236 /* number of levels -- only process most options on the finest level */ 237 if (isFine) { 238 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 239 if (!flg && snes->dm) { 240 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 241 levels++; 242 fas->usedmfornumberoflevels = PETSC_TRUE; 243 } 244 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 245 fastype = fas->fastype; 246 ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 247 if (flg) { 248 ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 249 } 250 251 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 252 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 253 if (flg) { 254 ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 255 } 256 257 ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 258 if (flg) { 259 ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 260 } 261 262 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 263 264 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 265 266 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 267 if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 268 } 269 270 ierr = PetscOptionsTail();CHKERRQ(ierr); 271 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 272 if (upflg) { 273 ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 274 } 275 if (downflg) { 276 ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 277 } 278 279 /* set up the default line search for coarse grid corrections */ 280 if (fas->fastype == SNES_FAS_ADDITIVE) { 281 if (!snes->linesearch) { 282 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 283 ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr); 284 } 285 } 286 287 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 288 /* recursive option setting for the smoothers */ 289 if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 290 PetscFunctionReturn(0); 291 } 292 293 #undef __FUNCT__ 294 #define __FUNCT__ "SNESView_FAS" 295 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 296 { 297 SNES_FAS *fas = (SNES_FAS *) snes->data; 298 PetscBool isFine, iascii; 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 = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 307 if (iascii) { 308 ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 309 if (fas->galerkin) { 310 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 311 } else { 312 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 313 } 314 for (i=0; i<fas->levels; i++) { 315 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 316 ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 317 ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 318 if (!i) { 319 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 320 } else { 321 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 322 } 323 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 324 ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 325 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 326 if (i && (smoothd == smoothu)) { 327 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 328 } else if (i) { 329 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 330 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 331 ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 332 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 333 } 334 } 335 } else { 336 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 337 } 338 } 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "FASDownSmooth" 344 /* 345 Defines the action of the downsmoother 346 */ 347 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 348 PetscErrorCode ierr = 0; 349 SNESConvergedReason reason; 350 Vec FPC; 351 SNES smoothd; 352 PetscFunctionBegin; 353 ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 354 ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 355 /* check convergence reason for the smoother */ 356 ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 357 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 358 snes->reason = SNES_DIVERGED_INNER; 359 PetscFunctionReturn(0); 360 } 361 ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 362 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "FASUpSmooth" 369 /* 370 Defines the action of the upsmoother 371 */ 372 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 373 PetscErrorCode ierr = 0; 374 SNESConvergedReason reason; 375 Vec FPC; 376 SNES smoothu; 377 PetscFunctionBegin; 378 379 ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 380 ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 381 /* check convergence reason for the smoother */ 382 ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 383 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 384 snes->reason = SNES_DIVERGED_INNER; 385 PetscFunctionReturn(0); 386 } 387 ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 388 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "SNESFASCreateCoarseVec" 394 /*@ 395 SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 396 397 Collective 398 399 Input Arguments: 400 . snes - SNESFAS 401 402 Output Arguments: 403 . Xcoarse - vector on level one coarser than snes 404 405 Level: developer 406 407 .seealso: SNESFASSetRestriction(), SNESFASRestrict() 408 @*/ 409 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 410 { 411 PetscErrorCode ierr; 412 SNES_FAS *fas = (SNES_FAS*)snes->data; 413 414 PetscFunctionBegin; 415 if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 416 else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 417 else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 418 PetscFunctionReturn(0); 419 } 420 421 #undef __FUNCT__ 422 #define __FUNCT__ "SNESFASRestrict" 423 /*@ 424 SNESFASRestrict - restrict a Vec to the next coarser level 425 426 Collective 427 428 Input Arguments: 429 + fine - SNES from which to restrict 430 - Xfine - vector to restrict 431 432 Output Arguments: 433 . Xcoarse - result of restriction 434 435 Level: developer 436 437 .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 438 @*/ 439 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 440 { 441 PetscErrorCode ierr; 442 SNES_FAS *fas = (SNES_FAS*)fine->data; 443 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 446 PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 447 PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 448 if (fas->inject) { 449 ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 450 } else { 451 ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 452 ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 453 } 454 PetscFunctionReturn(0); 455 } 456 457 #undef __FUNCT__ 458 #define __FUNCT__ "FASCoarseCorrection" 459 /* 460 461 Performs the FAS coarse correction as: 462 463 fine problem: F(x) = 0 464 coarse problem: F^c(x) = b^c 465 466 b^c = F^c(I^c_fx^f - I^c_fF(x)) 467 468 */ 469 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 470 PetscErrorCode ierr; 471 Vec X_c, Xo_c, F_c, B_c; 472 SNESConvergedReason reason; 473 SNES next; 474 Mat restrct, interpolate; 475 PetscFunctionBegin; 476 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 477 if (next) { 478 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 479 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 480 481 X_c = next->vec_sol; 482 Xo_c = next->work[0]; 483 F_c = next->vec_func; 484 B_c = next->vec_rhs; 485 486 ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 487 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 488 489 /* restrict the defect */ 490 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 491 492 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 493 next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 494 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 495 496 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 497 498 /* set initial guess of the coarse problem to the projected fine solution */ 499 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 500 501 /* recurse to the next level */ 502 next->vec_rhs = B_c; 503 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 504 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 505 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 506 snes->reason = SNES_DIVERGED_INNER; 507 PetscFunctionReturn(0); 508 } 509 510 /* correct as x <- x + I(x^c - Rx)*/ 511 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 512 ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 513 } 514 PetscFunctionReturn(0); 515 } 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "FASCycle_Additive" 519 /* 520 521 The additive cycle looks like: 522 523 xhat = x 524 xhat = dS(x, b) 525 x = coarsecorrection(xhat, b_d) 526 x = x + nu*(xhat - x); 527 (optional) x = uS(x, b) 528 529 With the coarse RHS (defect correction) as below. 530 531 */ 532 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 533 Vec F, B, Xhat; 534 Vec X_c, Xo_c, F_c, B_c; 535 PetscErrorCode ierr; 536 SNESConvergedReason reason; 537 PetscReal xnorm, fnorm, ynorm; 538 PetscBool lssuccess; 539 SNES next; 540 Mat restrct, interpolate; 541 PetscFunctionBegin; 542 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 543 F = snes->vec_func; 544 B = snes->vec_rhs; 545 Xhat = snes->work[3]; 546 ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 547 /* recurse first */ 548 if (next) { 549 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 550 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 551 ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 552 553 X_c = next->vec_sol; 554 Xo_c = next->work[0]; 555 F_c = next->vec_func; 556 B_c = next->vec_rhs; 557 558 ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 559 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 560 561 /* restrict the defect */ 562 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 563 564 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 565 next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 566 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 567 568 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 569 570 /* set initial guess of the coarse problem to the projected fine solution */ 571 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 572 573 /* recurse */ 574 next->vec_rhs = B_c; 575 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 576 577 /* smooth on this level */ 578 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 579 580 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 581 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 582 snes->reason = SNES_DIVERGED_INNER; 583 PetscFunctionReturn(0); 584 } 585 586 /* correct as x <- x + I(x^c - Rx)*/ 587 ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 588 ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 589 590 /* additive correction of the coarse direction*/ 591 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 592 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 593 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 594 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 595 if (!lssuccess) { 596 if (++snes->numFailures >= snes->maxFailures) { 597 snes->reason = SNES_DIVERGED_LINE_SEARCH; 598 PetscFunctionReturn(0); 599 } 600 } 601 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 602 } else { 603 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 604 } 605 PetscFunctionReturn(0); 606 } 607 608 #undef __FUNCT__ 609 #define __FUNCT__ "FASCycle_Multiplicative" 610 /* 611 612 Defines the FAS cycle as: 613 614 fine problem: F(x) = 0 615 coarse problem: F^c(x) = b^c 616 617 b^c = F^c(I^c_fx^f - I^c_fF(x)) 618 619 correction: 620 621 x = x + I(x^c - Rx) 622 623 */ 624 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 625 626 PetscErrorCode ierr; 627 Vec F,B; 628 SNES_FAS *fas = (SNES_FAS *)snes->data; 629 630 PetscFunctionBegin; 631 F = snes->vec_func; 632 B = snes->vec_rhs; 633 /* pre-smooth -- just update using the pre-smoother */ 634 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 635 636 if (fas->level != 0) { 637 ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 638 ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 639 } 640 PetscFunctionReturn(0); 641 } 642 643 #undef __FUNCT__ 644 #define __FUNCT__ "SNESSolve_FAS" 645 646 PetscErrorCode SNESSolve_FAS(SNES snes) 647 { 648 PetscErrorCode ierr; 649 PetscInt i, maxits; 650 Vec X, F; 651 PetscReal fnorm; 652 SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 653 DM dm; 654 PetscBool isFine; 655 656 PetscFunctionBegin; 657 maxits = snes->max_its; /* maximum number of iterations */ 658 snes->reason = SNES_CONVERGED_ITERATING; 659 X = snes->vec_sol; 660 F = snes->vec_func; 661 662 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 663 /*norm setup */ 664 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 665 snes->iter = 0; 666 snes->norm = 0.; 667 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 668 if (isFine || fas->monitor) { 669 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 670 if (snes->domainerror) { 671 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 672 PetscFunctionReturn(0); 673 } 674 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 675 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 676 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 677 snes->norm = fnorm; 678 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 679 SNESLogConvHistory(snes,fnorm,0); 680 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 681 682 /* set parameter for default relative tolerance convergence test */ 683 snes->ttol = fnorm*snes->rtol; 684 /* test convergence */ 685 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 686 if (snes->reason) PetscFunctionReturn(0); 687 } 688 689 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 690 for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 691 DM dmcoarse; 692 ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 693 ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 694 dm = dmcoarse; 695 } 696 697 for (i = 0; i < maxits; i++) { 698 /* Call general purpose update function */ 699 700 if (snes->ops->update) { 701 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 702 } 703 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 704 ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 705 } else { 706 ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 707 } 708 709 /* check for FAS cycle divergence */ 710 if (snes->reason != SNES_CONVERGED_ITERATING) { 711 PetscFunctionReturn(0); 712 } 713 if (isFine || fas->monitor) { 714 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 715 } 716 /* Monitor convergence */ 717 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 718 snes->iter = i+1; 719 snes->norm = fnorm; 720 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 721 SNESLogConvHistory(snes,snes->norm,0); 722 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 723 /* Test for convergence */ 724 if (isFine) { 725 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 726 if (snes->reason) break; 727 } 728 } 729 if (i == maxits) { 730 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 731 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 732 } 733 PetscFunctionReturn(0); 734 } 735