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