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