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