1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> 3 4 /*MC 5 Full Approximation Scheme nonlinear multigrid solver. 6 7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8 9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10 M*/ 11 12 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16 extern PetscErrorCode SNESSolve_FAS(SNES snes); 17 extern PetscErrorCode SNESReset_FAS(SNES snes); 18 19 EXTERN_C_BEGIN 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "SNESCreate_FAS" 23 PetscErrorCode SNESCreate_FAS(SNES snes) 24 { 25 SNES_FAS * fas; 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 snes->ops->destroy = SNESDestroy_FAS; 30 snes->ops->setup = SNESSetUp_FAS; 31 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32 snes->ops->view = SNESView_FAS; 33 snes->ops->solve = SNESSolve_FAS; 34 snes->ops->reset = SNESReset_FAS; 35 36 snes->usesksp = PETSC_FALSE; 37 snes->usespc = PETSC_FALSE; 38 39 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 40 snes->data = (void*) fas; 41 fas->level = 0; 42 fas->levels = 1; 43 fas->n_cycles = 1; 44 fas->max_up_it = 1; 45 fas->max_down_it = 1; 46 fas->upsmooth = PETSC_NULL; 47 fas->downsmooth = PETSC_NULL; 48 fas->next = PETSC_NULL; 49 fas->interpolate = PETSC_NULL; 50 fas->restrct = PETSC_NULL; 51 fas->inject = PETSC_NULL; 52 53 PetscFunctionReturn(0); 54 } 55 EXTERN_C_END 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "SNESFASGetLevels" 59 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 60 SNES_FAS * fas = (SNES_FAS *)snes->data; 61 PetscFunctionBegin; 62 *levels = fas->levels; 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "SNESFASSetCycles" 68 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 69 SNES_FAS * fas = (SNES_FAS *)snes->data; 70 PetscErrorCode ierr; 71 PetscFunctionBegin; 72 fas->n_cycles = cycles; 73 if (fas->next) { 74 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 75 } 76 PetscFunctionReturn(0); 77 } 78 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "SNESFASSetCyclesOnLevel" 82 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, PetscInt cycles) { 83 SNES_FAS * fas = (SNES_FAS *)snes->data; 84 PetscInt top_level = fas->level,i; 85 86 PetscFunctionBegin; 87 if (level > top_level) 88 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 89 /* get to the correct level */ 90 for (i = fas->level; i > level; i--) { 91 fas = (SNES_FAS *)fas->next->data; 92 } 93 if (fas->level != level) 94 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 95 fas->n_cycles = cycles; 96 PetscFunctionReturn(0); 97 } 98 99 100 #undef __FUNCT__ 101 #define __FUNCT__ "SNESFASGetSNES" 102 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 103 SNES_FAS * fas = (SNES_FAS *)snes->data; 104 PetscInt levels = fas->level; 105 PetscInt i; 106 PetscFunctionBegin; 107 *lsnes = snes; 108 if (fas->level < level) { 109 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 110 } 111 if (level > levels - 1) { 112 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 113 } 114 for (i = fas->level; i > level; i--) { 115 *lsnes = fas->next; 116 fas = (SNES_FAS *)(*lsnes)->data; 117 } 118 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 119 PetscFunctionReturn(0); 120 } 121 122 #undef __FUNCT__ 123 #define __FUNCT__ "SNESFASSetLevels" 124 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 125 PetscErrorCode ierr; 126 PetscInt i; 127 SNES_FAS * fas = (SNES_FAS *)snes->data; 128 MPI_Comm comm; 129 PetscFunctionBegin; 130 comm = ((PetscObject)snes)->comm; 131 if (levels == fas->levels) { 132 if (!comms) 133 PetscFunctionReturn(0); 134 } 135 /* user has changed the number of levels; reset */ 136 ierr = SNESReset(snes);CHKERRQ(ierr); 137 /* destroy any coarser levels if necessary */ 138 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 139 fas->next = PETSC_NULL; 140 /* setup the finest level */ 141 for (i = levels - 1; i >= 0; i--) { 142 if (comms) comm = comms[i]; 143 fas->level = i; 144 fas->levels = levels; 145 fas->next = PETSC_NULL; 146 if (i > 0) { 147 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 148 ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 149 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 150 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 151 fas = (SNES_FAS *)fas->next->data; 152 } 153 } 154 PetscFunctionReturn(0); 155 } 156 157 #undef __FUNCT__ 158 #define __FUNCT__ "SNESFASSetInterpolation" 159 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 160 SNES_FAS * fas = (SNES_FAS *)snes->data; 161 PetscInt top_level = fas->level,i; 162 163 PetscFunctionBegin; 164 if (level > top_level) 165 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 166 /* get to the correct level */ 167 for (i = fas->level; i > level; i--) { 168 fas = (SNES_FAS *)fas->next->data; 169 } 170 if (fas->level != level) 171 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 172 fas->interpolate = mat; 173 PetscFunctionReturn(0); 174 } 175 176 #undef __FUNCT__ 177 #define __FUNCT__ "SNESFASSetRestriction" 178 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 179 SNES_FAS * fas = (SNES_FAS *)snes->data; 180 PetscInt top_level = fas->level,i; 181 182 PetscFunctionBegin; 183 if (level > top_level) 184 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 185 /* get to the correct level */ 186 for (i = fas->level; i > level; i--) { 187 fas = (SNES_FAS *)fas->next->data; 188 } 189 if (fas->level != level) 190 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 191 fas->restrct = mat; 192 PetscFunctionReturn(0); 193 } 194 195 196 #undef __FUNCT__ 197 #define __FUNCT__ "SNESFASSetRScale" 198 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 199 SNES_FAS * fas = (SNES_FAS *)snes->data; 200 PetscInt top_level = fas->level,i; 201 202 PetscFunctionBegin; 203 if (level > top_level) 204 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 205 /* get to the correct level */ 206 for (i = fas->level; i > level; i--) { 207 fas = (SNES_FAS *)fas->next->data; 208 } 209 if (fas->level != level) 210 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 211 fas->rscale = rscale; 212 PetscFunctionReturn(0); 213 } 214 215 216 #undef __FUNCT__ 217 #define __FUNCT__ "SNESFASSetInjection" 218 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 219 SNES_FAS * fas = (SNES_FAS *)snes->data; 220 PetscInt top_level = fas->level,i; 221 222 PetscFunctionBegin; 223 if (level > top_level) 224 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 225 /* get to the correct level */ 226 for (i = fas->level; i > level; i--) { 227 fas = (SNES_FAS *)fas->next->data; 228 } 229 if (fas->level != level) 230 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 231 fas->inject = mat; 232 PetscFunctionReturn(0); 233 } 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "SNESReset_FAS" 237 PetscErrorCode SNESReset_FAS(SNES snes) 238 { 239 PetscErrorCode ierr = 0; 240 SNES_FAS * fas = (SNES_FAS *)snes->data; 241 242 PetscFunctionBegin; 243 /* destroy local data created in SNESSetup_FAS */ 244 #if 0 245 /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 246 #endif 247 if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr); 248 #if 0 249 #endif 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "SNESDestroy_FAS" 255 PetscErrorCode SNESDestroy_FAS(SNES snes) 256 { 257 SNES_FAS * fas = (SNES_FAS *)snes->data; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 /* recursively resets and then destroys */ 262 ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 263 if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 264 if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 265 if (fas->inject) { 266 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 267 } 268 if (fas->interpolate == fas->restrct) { 269 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 270 fas->restrct = PETSC_NULL; 271 } else { 272 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 273 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 274 } 275 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 276 if (snes->work) ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 277 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 278 ierr = PetscFree(fas);CHKERRQ(ierr); 279 280 PetscFunctionReturn(0); 281 } 282 283 #undef __FUNCT__ 284 #define __FUNCT__ "SNESSetUp_FAS" 285 PetscErrorCode SNESSetUp_FAS(SNES snes) 286 { 287 SNES_FAS *fas = (SNES_FAS *) snes->data,*tmp; 288 PetscErrorCode ierr; 289 VecScatter injscatter; 290 PetscInt dm_levels; 291 292 293 PetscFunctionBegin; 294 /* reset to use the DM if appropriate */ 295 296 if (fas->usedmfornumberoflevels && !fas->level) { 297 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 298 dm_levels++; 299 if (dm_levels > fas->levels) { /* the problem is now being solved on a finer grid */ 300 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 301 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 302 } 303 } 304 305 /* should call the SNESSetFromOptions() only when appropriate */ 306 tmp = fas; 307 while (tmp) { 308 if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);} 309 if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);} 310 tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0; 311 } 312 313 if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */ 314 /* gets the solver ready for solution */ 315 if (fas->next) { 316 if (snes->ops->computegs) { 317 ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 318 } 319 } 320 if (snes->dm) { 321 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 322 if (fas->next) { 323 /* for now -- assume the DM and the evaluation functions have been set externally */ 324 if (!fas->next->dm) { 325 ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 326 ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 327 } 328 /* set the interpolation and restriction from the DM */ 329 if (!fas->interpolate) { 330 ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 331 fas->restrct = fas->interpolate; 332 } 333 /* set the injection from the DM */ 334 if (!fas->inject) { 335 ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 336 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 337 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 338 } 339 } 340 /* set the DMs of the pre and post-smoothers here */ 341 if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 342 if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 343 } 344 if (fas->next) { 345 /* gotta set up the solution vector for this to work */ 346 if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 347 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 348 } 349 /* got to set them all up at once */ 350 PetscFunctionReturn(0); 351 } 352 353 #undef __FUNCT__ 354 #define __FUNCT__ "SNESSetFromOptions_FAS" 355 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 356 { 357 SNES_FAS *fas = (SNES_FAS *) snes->data; 358 PetscInt levels = 1; 359 PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 360 PetscErrorCode ierr; 361 const char * def_smooth = SNESNRICHARDSON; 362 char pre_type[256]; 363 char post_type[256]; 364 char monfilename[PETSC_MAX_PATH_LEN]; 365 366 PetscFunctionBegin; 367 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 368 369 /* number of levels -- only process on the finest level */ 370 if (fas->levels - 1 == fas->level) { 371 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 372 if (!flg && snes->dm) { 373 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 374 levels++; 375 fas->usedmfornumberoflevels = PETSC_TRUE; 376 } 377 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 378 } 379 380 /* type of pre and/or post smoothers -- set both at once */ 381 ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 382 ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 383 ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 384 if (smoothflg) { 385 ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 386 } else { 387 ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 388 ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 389 } 390 391 /* options for the number of preconditioning cycles and cycle type */ 392 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 393 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 394 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 395 396 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 397 398 /* other options for the coarsest level */ 399 if (fas->level == 0) { 400 ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 401 } 402 403 ierr = PetscOptionsTail();CHKERRQ(ierr); 404 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 405 if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->ops->computegs)) { 406 const char *prefix; 407 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 408 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 409 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 410 if (fas->level || (fas->levels == 1)) { 411 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr); 412 } else { 413 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 414 } 415 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 416 ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 417 if (snes->ops->computefunction) { 418 ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 419 } 420 } 421 422 if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->ops->computegs)) { 423 const char *prefix; 424 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 425 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 426 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 427 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr); 428 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 429 ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 430 if (snes->ops->computefunction) { 431 ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 432 } 433 } 434 if (fas->upsmooth) { 435 ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 436 } 437 438 if (fas->downsmooth) { 439 ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 440 } 441 442 if (monflg) { 443 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 444 /* set the monitors for the upsmoother and downsmoother */ 445 if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 446 if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 447 } 448 449 /* recursive option setting for the smoothers */ 450 if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 451 PetscFunctionReturn(0); 452 } 453 454 #undef __FUNCT__ 455 #define __FUNCT__ "SNESView_FAS" 456 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 457 { 458 SNES_FAS *fas = (SNES_FAS *) snes->data; 459 PetscBool iascii; 460 PetscErrorCode ierr; 461 PetscInt levels = fas->levels; 462 PetscInt i; 463 464 PetscFunctionBegin; 465 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 466 if (iascii) { 467 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 468 ierr = PetscViewerASCIIPushTab(viewer); 469 for (i = levels - 1; i >= 0; i--) { 470 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 471 if (fas->upsmooth) { 472 ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIPushTab(viewer); 474 ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 475 ierr = PetscViewerASCIIPopTab(viewer); 476 } else { 477 ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 478 } 479 if (fas->downsmooth) { 480 ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 481 ierr = PetscViewerASCIIPushTab(viewer); 482 ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 483 ierr = PetscViewerASCIIPopTab(viewer); 484 } else { 485 ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 486 } 487 if (fas->next) fas = (SNES_FAS *)fas->next->data; 488 } 489 ierr = PetscViewerASCIIPopTab(viewer); 490 } else { 491 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 492 } 493 PetscFunctionReturn(0); 494 } 495 496 /* 497 498 Defines the FAS cycle as: 499 500 fine problem: F(x) = 0 501 coarse problem: F^c(x) = b^c 502 503 b^c = F^c(I^c_fx^f - I^c_fF(x)) 504 505 correction: 506 507 x = x + I(x^c - Rx) 508 509 */ 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "FASCycle_Private" 513 PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 514 515 PetscErrorCode ierr; 516 Vec X_c, Xo_c, F_c, B_c,F,B; 517 SNES_FAS * fas = (SNES_FAS *)snes->data; 518 PetscInt i, k; 519 PetscScalar fnorm; 520 521 PetscFunctionBegin; 522 F = snes->vec_func; 523 B = snes->vec_rhs; 524 /* pre-smooth -- just update using the pre-smoother */ 525 if (fas->downsmooth) { 526 ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 527 } else if (snes->ops->computegs) { 528 if (fas->monitor) { 529 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 530 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 531 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 532 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", 0, fnorm);CHKERRQ(ierr); 533 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 534 } 535 for (k = 0; k < fas->max_up_it; k++) { 536 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 537 if (snes->monitor) { 538 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 539 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 540 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 541 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", k+1, fnorm);CHKERRQ(ierr); 542 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 543 } 544 } 545 } else if (snes->pc) { 546 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 547 } 548 if (fas->next) { 549 for (i = 0; i < fas->n_cycles; i++) { 550 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 551 552 X_c = fas->next->vec_sol; 553 Xo_c = fas->next->work[0]; 554 F_c = fas->next->vec_func; 555 B_c = fas->next->work[1]; 556 557 /* inject the solution */ 558 if (fas->inject) { 559 ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 560 } else { 561 ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 562 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 563 } 564 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 565 566 /* restrict the defect */ 567 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 568 569 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 570 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 571 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 572 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 573 574 /* set initial guess of the coarse problem to the projected fine solution */ 575 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 576 577 /* recurse to the next level */ 578 fas->next->vec_rhs = B_c; 579 ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); 580 fas->next->vec_rhs = PETSC_NULL; 581 582 /* correct as x <- x + I(x^c - Rx)*/ 583 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 584 ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 585 } 586 } 587 /* up-smooth -- just update using the down-smoother */ 588 if (fas->level != 0) { 589 if (fas->upsmooth) { 590 ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 591 } else if (snes->ops->computegs) { 592 if (fas->monitor) { 593 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 594 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 595 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 596 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", 0, fnorm);CHKERRQ(ierr); 597 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 598 } 599 for (k = 0; k < fas->max_down_it; k++) { 600 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 601 if (fas->monitor) { 602 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 603 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 604 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 605 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", k+1, fnorm);CHKERRQ(ierr); 606 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 607 } 608 } 609 } else if (snes->pc) { 610 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 611 } 612 } 613 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 614 615 PetscFunctionReturn(0); 616 } 617 618 #undef __FUNCT__ 619 #define __FUNCT__ "SNESSolve_FAS" 620 621 PetscErrorCode SNESSolve_FAS(SNES snes) 622 { 623 PetscErrorCode ierr; 624 PetscInt i, maxits; 625 Vec X, B,F; 626 PetscReal fnorm; 627 PetscFunctionBegin; 628 maxits = snes->max_its; /* maximum number of iterations */ 629 snes->reason = SNES_CONVERGED_ITERATING; 630 X = snes->vec_sol; 631 B = snes->vec_rhs; 632 F = snes->vec_func; 633 634 /*norm setup */ 635 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 636 snes->iter = 0; 637 snes->norm = 0.; 638 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 639 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 640 if (snes->domainerror) { 641 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 642 PetscFunctionReturn(0); 643 } 644 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 645 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 646 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 647 snes->norm = fnorm; 648 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 649 SNESLogConvHistory(snes,fnorm,0); 650 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 651 652 /* set parameter for default relative tolerance convergence test */ 653 snes->ttol = fnorm*snes->rtol; 654 /* test convergence */ 655 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 656 if (snes->reason) PetscFunctionReturn(0); 657 for (i = 0; i < maxits; i++) { 658 /* Call general purpose update function */ 659 660 if (snes->ops->update) { 661 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 662 } 663 ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 664 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 665 /* Monitor convergence */ 666 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 667 snes->iter = i+1; 668 snes->norm = fnorm; 669 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 670 SNESLogConvHistory(snes,snes->norm,0); 671 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 672 /* Test for convergence */ 673 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 674 if (snes->reason) break; 675 } 676 if (i == maxits) { 677 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 678 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 679 } 680 PetscFunctionReturn(0); 681 } 682