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