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 if (next) { 199 /* gotta set up the solution vector for this to work */ 200 if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 201 if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 202 ierr = SNESSetUp(next);CHKERRQ(ierr); 203 } 204 /* setup FAS work vectors */ 205 if (fas->galerkin) { 206 ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 207 ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 208 } 209 PetscFunctionReturn(0); 210 } 211 212 #undef __FUNCT__ 213 #define __FUNCT__ "SNESSetFromOptions_FAS" 214 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 215 { 216 SNES_FAS *fas = (SNES_FAS *) snes->data; 217 PetscInt levels = 1; 218 PetscBool flg, upflg, downflg, monflg, galerkinflg; 219 PetscErrorCode ierr; 220 char monfilename[PETSC_MAX_PATH_LEN]; 221 SNESFASType fastype; 222 const char *optionsprefix; 223 SNESLineSearch linesearch; 224 PetscInt m, n_up, n_down; 225 SNES next; 226 PetscBool isFine; 227 228 PetscFunctionBegin; 229 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 230 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 231 232 /* number of levels -- only process most options on the finest level */ 233 if (isFine) { 234 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 235 if (!flg && snes->dm) { 236 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 237 levels++; 238 fas->usedmfornumberoflevels = PETSC_TRUE; 239 } 240 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 241 fastype = fas->fastype; 242 ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 243 if (flg) { 244 ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 245 } 246 247 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 248 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 249 if (flg) { 250 ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 251 } 252 253 ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 254 if (flg) { 255 ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 256 } 257 258 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 259 260 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 261 262 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 263 if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 264 } 265 266 ierr = PetscOptionsTail();CHKERRQ(ierr); 267 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 268 if (upflg) { 269 ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 270 } 271 if (downflg) { 272 ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 273 } 274 275 /* set up the default line search for coarse grid corrections */ 276 if (fas->fastype == SNES_FAS_ADDITIVE) { 277 if (!snes->linesearch) { 278 ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 279 ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr); 280 } 281 } 282 283 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 284 /* recursive option setting for the smoothers */ 285 if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 286 PetscFunctionReturn(0); 287 } 288 289 #undef __FUNCT__ 290 #define __FUNCT__ "SNESView_FAS" 291 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 292 { 293 SNES_FAS *fas = (SNES_FAS *) snes->data; 294 PetscBool isFine, iascii; 295 PetscInt i; 296 PetscErrorCode ierr; 297 SNES smoothu, smoothd, levelsnes; 298 299 PetscFunctionBegin; 300 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 301 if (isFine) { 302 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 303 if (iascii) { 304 ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 305 if (fas->galerkin) { 306 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 307 } else { 308 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 309 } 310 for (i=0; i<fas->levels; i++) { 311 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 312 ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 313 ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 314 if (!i) { 315 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 316 } else { 317 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 318 } 319 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 320 ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 321 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 322 if (i && (smoothd == smoothu)) { 323 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 324 } else if (i) { 325 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 326 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 327 ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 328 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 329 } 330 } 331 } else { 332 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 333 } 334 } 335 PetscFunctionReturn(0); 336 } 337 338 #undef __FUNCT__ 339 #define __FUNCT__ "FASDownSmooth" 340 /* 341 Defines the action of the downsmoother 342 */ 343 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 344 PetscErrorCode ierr = 0; 345 SNESConvergedReason reason; 346 Vec FPC; 347 SNES smoothd; 348 PetscFunctionBegin; 349 ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 350 ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 351 /* check convergence reason for the smoother */ 352 ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 353 if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 354 snes->reason = SNES_DIVERGED_INNER; 355 PetscFunctionReturn(0); 356 } 357 ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 358 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 359 PetscFunctionReturn(0); 360 } 361 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "FASUpSmooth" 365 /* 366 Defines the action of the upsmoother 367 */ 368 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 369 PetscErrorCode ierr = 0; 370 SNESConvergedReason reason; 371 Vec FPC; 372 SNES smoothu; 373 PetscFunctionBegin; 374 375 ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 376 ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 377 /* check convergence reason for the smoother */ 378 ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 379 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 380 snes->reason = SNES_DIVERGED_INNER; 381 PetscFunctionReturn(0); 382 } 383 ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 384 ierr = VecCopy(FPC, F);CHKERRQ(ierr); 385 PetscFunctionReturn(0); 386 } 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "SNESFASCreateCoarseVec" 390 /*@ 391 SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 392 393 Collective 394 395 Input Arguments: 396 . snes - SNESFAS 397 398 Output Arguments: 399 . Xcoarse - vector on level one coarser than snes 400 401 Level: developer 402 403 .seealso: SNESFASSetRestriction(), SNESFASRestrict() 404 @*/ 405 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 406 { 407 PetscErrorCode ierr; 408 SNES_FAS *fas = (SNES_FAS*)snes->data; 409 410 PetscFunctionBegin; 411 if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 412 else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 413 else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "SNESFASRestrict" 419 /*@ 420 SNESFASRestrict - restrict a Vec to the next coarser level 421 422 Collective 423 424 Input Arguments: 425 + fine - SNES from which to restrict 426 - Xfine - vector to restrict 427 428 Output Arguments: 429 . Xcoarse - result of restriction 430 431 Level: developer 432 433 .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 434 @*/ 435 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 436 { 437 PetscErrorCode ierr; 438 SNES_FAS *fas = (SNES_FAS*)fine->data; 439 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 442 PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 443 PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 444 if (fas->inject) { 445 ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 446 } else { 447 ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 448 ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 449 } 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "FASCoarseCorrection" 455 /* 456 457 Performs the FAS coarse correction as: 458 459 fine problem: F(x) = 0 460 coarse problem: F^c(x) = b^c 461 462 b^c = F^c(I^c_fx^f - I^c_fF(x)) 463 464 */ 465 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 466 PetscErrorCode ierr; 467 Vec X_c, Xo_c, F_c, B_c; 468 SNESConvergedReason reason; 469 SNES next; 470 Mat restrct, interpolate; 471 PetscFunctionBegin; 472 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 473 if (next) { 474 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 475 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 476 477 X_c = next->vec_sol; 478 Xo_c = next->work[0]; 479 F_c = next->vec_func; 480 B_c = next->vec_rhs; 481 482 ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 483 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 484 485 /* restrict the defect */ 486 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 487 488 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 489 next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 490 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 491 492 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 493 494 /* set initial guess of the coarse problem to the projected fine solution */ 495 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 496 497 /* recurse to the next level */ 498 next->vec_rhs = B_c; 499 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 500 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 501 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 502 snes->reason = SNES_DIVERGED_INNER; 503 PetscFunctionReturn(0); 504 } 505 506 /* correct as x <- x + I(x^c - Rx)*/ 507 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 508 ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 509 } 510 PetscFunctionReturn(0); 511 } 512 513 #undef __FUNCT__ 514 #define __FUNCT__ "FASCycle_Additive" 515 /* 516 517 The additive cycle looks like: 518 519 xhat = x 520 xhat = dS(x, b) 521 x = coarsecorrection(xhat, b_d) 522 x = x + nu*(xhat - x); 523 (optional) x = uS(x, b) 524 525 With the coarse RHS (defect correction) as below. 526 527 */ 528 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 529 Vec F, B, Xhat; 530 Vec X_c, Xo_c, F_c, B_c; 531 PetscErrorCode ierr; 532 SNESConvergedReason reason; 533 PetscReal xnorm, fnorm, ynorm; 534 PetscBool lssuccess; 535 SNES next; 536 Mat restrct, interpolate; 537 PetscFunctionBegin; 538 ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 539 F = snes->vec_func; 540 B = snes->vec_rhs; 541 Xhat = snes->work[3]; 542 ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 543 /* recurse first */ 544 if (next) { 545 ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 546 ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 547 ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 548 549 X_c = next->vec_sol; 550 Xo_c = next->work[0]; 551 F_c = next->vec_func; 552 B_c = next->vec_rhs; 553 554 ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 555 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 556 557 /* restrict the defect */ 558 ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 559 560 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 561 next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 562 ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 563 564 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 565 566 /* set initial guess of the coarse problem to the projected fine solution */ 567 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 568 569 /* recurse */ 570 next->vec_rhs = B_c; 571 ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 572 573 /* smooth on this level */ 574 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 575 576 ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 577 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 578 snes->reason = SNES_DIVERGED_INNER; 579 PetscFunctionReturn(0); 580 } 581 582 /* correct as x <- x + I(x^c - Rx)*/ 583 ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 584 ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 585 586 /* additive correction of the coarse direction*/ 587 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 588 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 589 ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 590 ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 591 if (!lssuccess) { 592 if (++snes->numFailures >= snes->maxFailures) { 593 snes->reason = SNES_DIVERGED_LINE_SEARCH; 594 PetscFunctionReturn(0); 595 } 596 } 597 ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 598 } else { 599 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 600 } 601 PetscFunctionReturn(0); 602 } 603 604 #undef __FUNCT__ 605 #define __FUNCT__ "FASCycle_Multiplicative" 606 /* 607 608 Defines the FAS cycle as: 609 610 fine problem: F(x) = 0 611 coarse problem: F^c(x) = b^c 612 613 b^c = F^c(I^c_fx^f - I^c_fF(x)) 614 615 correction: 616 617 x = x + I(x^c - Rx) 618 619 */ 620 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 621 622 PetscErrorCode ierr; 623 Vec F,B; 624 SNES_FAS *fas = (SNES_FAS *)snes->data; 625 626 PetscFunctionBegin; 627 F = snes->vec_func; 628 B = snes->vec_rhs; 629 /* pre-smooth -- just update using the pre-smoother */ 630 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 631 632 if (fas->level != 0) { 633 ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 634 ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 635 } 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "SNESSolve_FAS" 641 642 PetscErrorCode SNESSolve_FAS(SNES snes) 643 { 644 PetscErrorCode ierr; 645 PetscInt i, maxits; 646 Vec X, F; 647 PetscReal fnorm; 648 SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 649 DM dm; 650 PetscBool isFine; 651 652 PetscFunctionBegin; 653 maxits = snes->max_its; /* maximum number of iterations */ 654 snes->reason = SNES_CONVERGED_ITERATING; 655 X = snes->vec_sol; 656 F = snes->vec_func; 657 658 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 659 /*norm setup */ 660 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 661 snes->iter = 0; 662 snes->norm = 0.; 663 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 664 if (isFine || fas->monitor) { 665 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 666 if (snes->domainerror) { 667 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 668 PetscFunctionReturn(0); 669 } 670 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 671 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 672 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 673 snes->norm = fnorm; 674 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 675 SNESLogConvHistory(snes,fnorm,0); 676 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 677 678 /* set parameter for default relative tolerance convergence test */ 679 snes->ttol = fnorm*snes->rtol; 680 /* test convergence */ 681 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 682 if (snes->reason) PetscFunctionReturn(0); 683 } 684 685 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 686 for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 687 DM dmcoarse; 688 ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 689 ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 690 dm = dmcoarse; 691 } 692 693 for (i = 0; i < maxits; i++) { 694 /* Call general purpose update function */ 695 696 if (snes->ops->update) { 697 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 698 } 699 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 700 ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 701 } else { 702 ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 703 } 704 705 /* check for FAS cycle divergence */ 706 if (snes->reason != SNES_CONVERGED_ITERATING) { 707 PetscFunctionReturn(0); 708 } 709 if (isFine || fas->monitor) { 710 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 711 } 712 /* Monitor convergence */ 713 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 714 snes->iter = i+1; 715 snes->norm = fnorm; 716 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 717 SNESLogConvHistory(snes,snes->norm,0); 718 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 719 /* Test for convergence */ 720 if (isFine) { 721 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 722 if (snes->reason) break; 723 } 724 } 725 if (i == maxits) { 726 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 727 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 728 } 729 PetscFunctionReturn(0); 730 } 731