1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 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 fas->monitor = PETSC_NULL; 53 fas->usedmfornumberoflevels = PETSC_FALSE; 54 fas->useGS = PETSC_FALSE; 55 PetscFunctionReturn(0); 56 } 57 EXTERN_C_END 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "SNESFASGetLevels" 61 /*@ 62 SNESFASGetLevels - Gets the number of levels in a FAS. 63 64 Input Parameter: 65 . snes - the preconditioner context 66 67 Output parameter: 68 . levels - the number of levels 69 70 Level: advanced 71 72 .keywords: MG, get, levels, multigrid 73 74 .seealso: SNESFASSetLevels(), PCMGGetLevels() 75 @*/ 76 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 77 SNES_FAS * fas = (SNES_FAS *)snes->data; 78 PetscFunctionBegin; 79 *levels = fas->levels; 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "SNESFASSetCycles" 85 /*@ 86 SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 87 complicated cycling. 88 89 Logically Collective on SNES 90 91 Input Parameters: 92 + snes - the multigrid context 93 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 94 95 Options Database Key: 96 $ -snes_fas_cycles 1 or 2 97 98 Level: advanced 99 100 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 101 102 .seealso: SNESFASSetCyclesOnLevel() 103 @*/ 104 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 105 SNES_FAS * fas = (SNES_FAS *)snes->data; 106 PetscErrorCode ierr; 107 PetscFunctionBegin; 108 fas->n_cycles = cycles; 109 if (fas->next) { 110 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 111 } 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "SNESFASSetCyclesOnLevel" 117 /*@ 118 SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 119 120 Logically Collective on SNES 121 122 Input Parameters: 123 + snes - the multigrid context 124 . level - the level to set the number of cycles on 125 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 126 127 Level: advanced 128 129 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 130 131 .seealso: SNESFASSetCycles() 132 @*/ 133 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 134 SNES_FAS * fas = (SNES_FAS *)snes->data; 135 PetscInt top_level = fas->level,i; 136 137 PetscFunctionBegin; 138 if (level > top_level) 139 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 140 /* get to the correct level */ 141 for (i = fas->level; i > level; i--) { 142 fas = (SNES_FAS *)fas->next->data; 143 } 144 if (fas->level != level) 145 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 146 fas->n_cycles = cycles; 147 PetscFunctionReturn(0); 148 } 149 150 #undef __FUNCT__ 151 #define __FUNCT__ "SNESFASSetGS" 152 /*@C 153 SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 154 Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 155 and nonlinear preconditioners. 156 157 Logically Collective on SNES 158 159 Input Parameters: 160 + snes - the multigrid context 161 . gsfunc - the nonlinear smoother function 162 . ctx - the user context for the nonlinear smoother 163 - use_gs - whether to use the nonlinear smoother or not 164 165 Level: advanced 166 167 .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 168 169 .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 170 @*/ 171 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 172 PetscErrorCode ierr = 0; 173 SNES_FAS *fas = (SNES_FAS *)snes->data; 174 PetscFunctionBegin; 175 176 /* use or don't use it according to user wishes*/ 177 fas->useGS = use_gs; 178 if (gsfunc) { 179 ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 180 /* push the provided GS up the tree */ 181 if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 182 } else if (snes->ops->computegs) { 183 /* assume that the user has set the GS solver at this level */ 184 if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 185 } else if (use_gs) { 186 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 187 } 188 PetscFunctionReturn(0); 189 } 190 191 #undef __FUNCT__ 192 #define __FUNCT__ "SNESFASSetGSOnLevel" 193 /*@C 194 SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 195 196 Logically Collective on SNES 197 198 Input Parameters: 199 + snes - the multigrid context 200 . level - the level to set the nonlinear smoother on 201 . gsfunc - the nonlinear smoother function 202 . ctx - the user context for the nonlinear smoother 203 - use_gs - whether to use the nonlinear smoother or not 204 205 Level: advanced 206 207 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 208 209 .seealso: SNESSetGS(), SNESFASSetGS() 210 @*/ 211 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 212 SNES_FAS *fas = (SNES_FAS *)snes->data; 213 PetscErrorCode ierr; 214 PetscInt top_level = fas->level,i; 215 SNES cur_snes = snes; 216 PetscFunctionBegin; 217 if (level > top_level) 218 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 219 /* get to the correct level */ 220 for (i = fas->level; i > level; i--) { 221 fas = (SNES_FAS *)fas->next->data; 222 cur_snes = fas->next; 223 } 224 if (fas->level != level) 225 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 226 fas->useGS = use_gs; 227 if (gsfunc) { 228 ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 229 } 230 PetscFunctionReturn(0); 231 } 232 233 #undef __FUNCT__ 234 #define __FUNCT__ "SNESFASGetSNES" 235 /*@ 236 SNESFASGetSNES - Gets the SNES corresponding to a particular 237 level of the FAS hierarchy. 238 239 Input Parameters: 240 + snes - the multigrid context 241 level - the level to get 242 - lsnes - whether to use the nonlinear smoother or not 243 244 Level: advanced 245 246 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 247 248 .seealso: SNESFASSetLevels(), SNESFASGetLevels() 249 @*/ 250 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 251 SNES_FAS * fas = (SNES_FAS *)snes->data; 252 PetscInt levels = fas->level; 253 PetscInt i; 254 PetscFunctionBegin; 255 *lsnes = snes; 256 if (fas->level < level) { 257 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 258 } 259 if (level > levels - 1) { 260 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 261 } 262 for (i = fas->level; i > level; i--) { 263 *lsnes = fas->next; 264 fas = (SNES_FAS *)(*lsnes)->data; 265 } 266 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 267 PetscFunctionReturn(0); 268 } 269 270 #undef __FUNCT__ 271 #define __FUNCT__ "SNESFASSetLevels" 272 /*@C 273 SNESFASSetLevels - Sets the number of levels to use with FAS. 274 Must be called before any other FAS routine. 275 276 Input Parameters: 277 + snes - the snes context 278 . levels - the number of levels 279 - comms - optional communicators for each level; this is to allow solving the coarser 280 problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 281 Fortran. 282 283 Level: intermediate 284 285 Notes: 286 If the number of levels is one then the multigrid uses the -fas_levels prefix 287 for setting the level options rather than the -fas_coarse prefix. 288 289 .keywords: FAS, MG, set, levels, multigrid 290 291 .seealso: SNESFASGetLevels() 292 @*/ 293 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 294 PetscErrorCode ierr; 295 PetscInt i; 296 SNES_FAS * fas = (SNES_FAS *)snes->data; 297 MPI_Comm comm; 298 PetscFunctionBegin; 299 comm = ((PetscObject)snes)->comm; 300 if (levels == fas->levels) { 301 if (!comms) 302 PetscFunctionReturn(0); 303 } 304 /* user has changed the number of levels; reset */ 305 ierr = SNESReset(snes);CHKERRQ(ierr); 306 /* destroy any coarser levels if necessary */ 307 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 308 fas->next = PETSC_NULL; 309 /* setup the finest level */ 310 for (i = levels - 1; i >= 0; i--) { 311 if (comms) comm = comms[i]; 312 fas->level = i; 313 fas->levels = levels; 314 fas->next = PETSC_NULL; 315 if (i > 0) { 316 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 317 ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 318 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 319 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 320 fas = (SNES_FAS *)fas->next->data; 321 } 322 } 323 PetscFunctionReturn(0); 324 } 325 326 #undef __FUNCT__ 327 #define __FUNCT__ "SNESFASSetNumberSmoothUp" 328 /*@ 329 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 330 use on all levels. 331 332 Logically Collective on PC 333 334 Input Parameters: 335 + snes - the multigrid context 336 - n - the number of smoothing steps 337 338 Options Database Key: 339 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 340 341 Level: advanced 342 343 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 344 345 .seealso: SNESFASSetNumberSmoothDown() 346 @*/ 347 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 348 SNES_FAS * fas = (SNES_FAS *)snes->data; 349 PetscErrorCode ierr = 0; 350 PetscFunctionBegin; 351 fas->max_up_it = n; 352 if (fas->next) { 353 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 354 } 355 PetscFunctionReturn(0); 356 } 357 358 #undef __FUNCT__ 359 #define __FUNCT__ "SNESFASSetNumberSmoothDown" 360 /*@ 361 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 362 use on all levels. 363 364 Logically Collective on PC 365 366 Input Parameters: 367 + snes - the multigrid context 368 - n - the number of smoothing steps 369 370 Options Database Key: 371 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 372 373 Level: advanced 374 375 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 376 377 .seealso: SNESFASSetNumberSmoothUp() 378 @*/ 379 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 380 SNES_FAS * fas = (SNES_FAS *)snes->data; 381 PetscErrorCode ierr = 0; 382 PetscFunctionBegin; 383 fas->max_down_it = n; 384 if (fas->next) { 385 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 386 } 387 PetscFunctionReturn(0); 388 } 389 390 #undef __FUNCT__ 391 #define __FUNCT__ "SNESFASSetInterpolation" 392 /*@ 393 SNESFASSetInterpolation - Sets the function to be used to calculate the 394 interpolation from l-1 to the lth level 395 396 Input Parameters: 397 + snes - the multigrid context 398 . mat - the interpolation operator 399 - level - the level (0 is coarsest) to supply [do not supply 0] 400 401 Level: advanced 402 403 Notes: 404 Usually this is the same matrix used also to set the restriction 405 for the same level. 406 407 One can pass in the interpolation matrix or its transpose; PETSc figures 408 out from the matrix size which one it is. 409 410 .keywords: FAS, multigrid, set, interpolate, level 411 412 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 413 @*/ 414 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 415 SNES_FAS * fas = (SNES_FAS *)snes->data; 416 PetscInt top_level = fas->level,i; 417 418 PetscFunctionBegin; 419 if (level > top_level) 420 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 421 /* get to the correct level */ 422 for (i = fas->level; i > level; i--) { 423 fas = (SNES_FAS *)fas->next->data; 424 } 425 if (fas->level != level) 426 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 427 fas->interpolate = mat; 428 PetscFunctionReturn(0); 429 } 430 431 #undef __FUNCT__ 432 #define __FUNCT__ "SNESFASSetRestriction" 433 /*@ 434 SNESFASSetRestriction - Sets the function to be used to restrict the defect 435 from level l to l-1. 436 437 Input Parameters: 438 + snes - the multigrid context 439 . mat - the restriction matrix 440 - level - the level (0 is coarsest) to supply [Do not supply 0] 441 442 Level: advanced 443 444 Notes: 445 Usually this is the same matrix used also to set the interpolation 446 for the same level. 447 448 One can pass in the interpolation matrix or its transpose; PETSc figures 449 out from the matrix size which one it is. 450 451 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 452 is used. 453 454 .keywords: FAS, MG, set, multigrid, restriction, level 455 456 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 457 @*/ 458 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 459 SNES_FAS * fas = (SNES_FAS *)snes->data; 460 PetscInt top_level = fas->level,i; 461 462 PetscFunctionBegin; 463 if (level > top_level) 464 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 465 /* get to the correct level */ 466 for (i = fas->level; i > level; i--) { 467 fas = (SNES_FAS *)fas->next->data; 468 } 469 if (fas->level != level) 470 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 471 fas->restrct = mat; 472 PetscFunctionReturn(0); 473 } 474 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "SNESFASSetRScale" 478 /*@ 479 SNESFASSetRScale - Sets the scaling factor of the restriction 480 operator from level l to l-1. 481 482 Input Parameters: 483 + snes - the multigrid context 484 . rscale - the restriction scaling 485 - level - the level (0 is coarsest) to supply [Do not supply 0] 486 487 Level: advanced 488 489 Notes: 490 This is only used in the case that the injection is not set. 491 492 .keywords: FAS, MG, set, multigrid, restriction, level 493 494 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 495 @*/ 496 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 497 SNES_FAS * fas = (SNES_FAS *)snes->data; 498 PetscInt top_level = fas->level,i; 499 500 PetscFunctionBegin; 501 if (level > top_level) 502 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 503 /* get to the correct level */ 504 for (i = fas->level; i > level; i--) { 505 fas = (SNES_FAS *)fas->next->data; 506 } 507 if (fas->level != level) 508 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 509 fas->rscale = rscale; 510 PetscFunctionReturn(0); 511 } 512 513 514 #undef __FUNCT__ 515 #define __FUNCT__ "SNESFASSetInjection" 516 /*@ 517 SNESFASSetInjection - Sets the function to be used to inject the solution 518 from level l to l-1. 519 520 Input Parameters: 521 + snes - the multigrid context 522 . mat - the restriction matrix 523 - level - the level (0 is coarsest) to supply [Do not supply 0] 524 525 Level: advanced 526 527 Notes: 528 If you do not set this, the restriction and rscale is used to 529 project the solution instead. 530 531 .keywords: FAS, MG, set, multigrid, restriction, level 532 533 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 534 @*/ 535 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 536 SNES_FAS * fas = (SNES_FAS *)snes->data; 537 PetscInt top_level = fas->level,i; 538 539 PetscFunctionBegin; 540 if (level > top_level) 541 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 542 /* get to the correct level */ 543 for (i = fas->level; i > level; i--) { 544 fas = (SNES_FAS *)fas->next->data; 545 } 546 if (fas->level != level) 547 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 548 fas->inject = mat; 549 PetscFunctionReturn(0); 550 } 551 552 #undef __FUNCT__ 553 #define __FUNCT__ "SNESReset_FAS" 554 PetscErrorCode SNESReset_FAS(SNES snes) 555 { 556 PetscErrorCode ierr = 0; 557 SNES_FAS * fas = (SNES_FAS *)snes->data; 558 559 PetscFunctionBegin; 560 /* destroy local data created in SNESSetup_FAS */ 561 #if 0 562 /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 563 #endif 564 if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr); 565 #if 0 566 #endif 567 PetscFunctionReturn(0); 568 } 569 570 #undef __FUNCT__ 571 #define __FUNCT__ "SNESDestroy_FAS" 572 PetscErrorCode SNESDestroy_FAS(SNES snes) 573 { 574 SNES_FAS * fas = (SNES_FAS *)snes->data; 575 PetscErrorCode ierr; 576 577 PetscFunctionBegin; 578 /* recursively resets and then destroys */ 579 ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 580 if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 581 if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 582 if (fas->inject) { 583 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 584 } 585 if (fas->interpolate == fas->restrct) { 586 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 587 fas->restrct = PETSC_NULL; 588 } else { 589 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 590 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 591 } 592 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 593 if (snes->work) ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 594 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 595 ierr = PetscFree(fas);CHKERRQ(ierr); 596 597 PetscFunctionReturn(0); 598 } 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "SNESSetUp_FAS" 602 PetscErrorCode SNESSetUp_FAS(SNES snes) 603 { 604 SNES_FAS *fas = (SNES_FAS *) snes->data; 605 PetscErrorCode ierr; 606 VecScatter injscatter; 607 PetscInt dm_levels; 608 609 610 PetscFunctionBegin; 611 612 if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 613 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 614 dm_levels++; 615 if (dm_levels > fas->levels) { 616 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 617 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 618 } 619 } 620 621 if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */ 622 623 /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 624 if (fas->upsmooth) { 625 ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 626 if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 627 ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 628 } 629 if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 630 ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 631 } 632 if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 633 ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 634 } 635 } 636 if (fas->downsmooth) { 637 ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 638 if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 639 ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 640 } 641 if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 642 ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 643 } 644 if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 645 ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 646 } 647 } 648 649 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 650 if (fas->next) { 651 if (snes->ops->computefunction && !fas->next->ops->computefunction) { 652 ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 653 } 654 if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 655 ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 656 } 657 if (snes->ops->computegs && !fas->next->ops->computegs) { 658 ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 659 } 660 } 661 if (snes->dm) { 662 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 663 if (fas->next) { 664 /* for now -- assume the DM and the evaluation functions have been set externally */ 665 if (!fas->next->dm) { 666 ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 667 ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 668 } 669 /* set the interpolation and restriction from the DM */ 670 if (!fas->interpolate) { 671 ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 672 fas->restrct = fas->interpolate; 673 } 674 /* set the injection from the DM */ 675 if (!fas->inject) { 676 ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 677 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 678 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 679 } 680 } 681 /* set the DMs of the pre and post-smoothers here */ 682 if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 683 if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 684 } 685 686 if (fas->next) { 687 /* gotta set up the solution vector for this to work */ 688 if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 689 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 690 } 691 /* got to set them all up at once */ 692 PetscFunctionReturn(0); 693 } 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "SNESSetFromOptions_FAS" 697 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 698 { 699 SNES_FAS *fas = (SNES_FAS *) snes->data; 700 PetscInt levels = 1; 701 PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 702 PetscErrorCode ierr; 703 const char * def_smooth = SNESNRICHARDSON; 704 char pre_type[256]; 705 char post_type[256]; 706 char monfilename[PETSC_MAX_PATH_LEN]; 707 708 PetscFunctionBegin; 709 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 710 711 /* number of levels -- only process on the finest level */ 712 if (fas->levels - 1 == fas->level) { 713 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 714 if (!flg && snes->dm) { 715 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 716 levels++; 717 fas->usedmfornumberoflevels = PETSC_TRUE; 718 } 719 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 720 } 721 722 /* type of pre and/or post smoothers -- set both at once */ 723 ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 724 ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 725 ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 726 if (smoothflg) { 727 ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 728 } else { 729 ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 730 ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 731 } 732 733 /* options for the number of preconditioning cycles and cycle type */ 734 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 735 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 736 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 737 738 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 739 740 ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 741 742 /* other options for the coarsest level */ 743 if (fas->level == 0) { 744 ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 745 ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 746 } 747 748 ierr = PetscOptionsTail();CHKERRQ(ierr); 749 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 750 751 if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) { 752 const char *prefix; 753 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 754 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 755 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 756 if (fas->level || (fas->levels == 1)) { 757 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 758 } else { 759 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 760 } 761 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 762 ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 763 } 764 765 if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) { 766 const char *prefix; 767 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 768 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 769 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 770 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 771 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 772 ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 773 } 774 if (fas->upsmooth) { 775 ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 776 } 777 778 if (fas->downsmooth) { 779 ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 780 } 781 782 if (monflg) { 783 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 784 /* set the monitors for the upsmoother and downsmoother */ 785 if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 786 if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 787 } 788 789 /* recursive option setting for the smoothers */ 790 if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 791 PetscFunctionReturn(0); 792 } 793 794 #undef __FUNCT__ 795 #define __FUNCT__ "SNESView_FAS" 796 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 797 { 798 SNES_FAS *fas = (SNES_FAS *) snes->data; 799 PetscBool iascii; 800 PetscErrorCode ierr; 801 PetscInt levels = fas->levels; 802 PetscInt i; 803 804 PetscFunctionBegin; 805 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 806 if (iascii) { 807 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 808 ierr = PetscViewerASCIIPushTab(viewer); 809 for (i = levels - 1; i >= 0; i--) { 810 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 811 if (fas->upsmooth) { 812 ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 813 ierr = PetscViewerASCIIPushTab(viewer); 814 ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 815 ierr = PetscViewerASCIIPopTab(viewer); 816 } else { 817 ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 818 } 819 if (fas->downsmooth) { 820 ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 821 ierr = PetscViewerASCIIPushTab(viewer); 822 ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 823 ierr = PetscViewerASCIIPopTab(viewer); 824 } else { 825 ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 826 } 827 if (fas->useGS) { 828 ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 829 } 830 if (fas->next) fas = (SNES_FAS *)fas->next->data; 831 } 832 ierr = PetscViewerASCIIPopTab(viewer); 833 } else { 834 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 835 } 836 PetscFunctionReturn(0); 837 } 838 839 /* 840 841 Defines the FAS cycle as: 842 843 fine problem: F(x) = 0 844 coarse problem: F^c(x) = b^c 845 846 b^c = F^c(I^c_fx^f - I^c_fF(x)) 847 848 correction: 849 850 x = x + I(x^c - Rx) 851 852 */ 853 854 #undef __FUNCT__ 855 #define __FUNCT__ "FASCycle_Private" 856 PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 857 858 PetscErrorCode ierr; 859 Vec X_c, Xo_c, F_c, B_c,F,B; 860 SNES_FAS * fas = (SNES_FAS *)snes->data; 861 PetscInt i, k; 862 PetscReal fnorm; 863 864 PetscFunctionBegin; 865 F = snes->vec_func; 866 B = snes->vec_rhs; 867 /* pre-smooth -- just update using the pre-smoother */ 868 if (fas->downsmooth) { 869 ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 870 } else if (snes->ops->computegs) { 871 if (fas->monitor) { 872 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 873 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 874 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 875 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 876 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 877 } 878 for (k = 0; k < fas->max_up_it; k++) { 879 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 880 if (fas->monitor) { 881 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 882 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 883 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 884 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 885 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 886 } 887 } 888 } else if (snes->pc) { 889 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 890 } 891 if (fas->next) { 892 for (i = 0; i < fas->n_cycles; i++) { 893 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 894 895 X_c = fas->next->vec_sol; 896 Xo_c = fas->next->work[0]; 897 F_c = fas->next->vec_func; 898 B_c = fas->next->work[1]; 899 900 /* inject the solution */ 901 if (fas->inject) { 902 ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 903 } else { 904 ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 905 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 906 } 907 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 908 909 /* restrict the defect */ 910 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 911 912 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 913 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 914 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 915 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 916 917 /* set initial guess of the coarse problem to the projected fine solution */ 918 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 919 920 /* recurse to the next level */ 921 fas->next->vec_rhs = B_c; 922 ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); 923 fas->next->vec_rhs = PETSC_NULL; 924 925 /* correct as x <- x + I(x^c - Rx)*/ 926 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 927 ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 928 } 929 } 930 /* up-smooth -- just update using the down-smoother */ 931 if (fas->level != 0) { 932 if (fas->upsmooth) { 933 ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 934 } else if (snes->ops->computegs) { 935 if (fas->monitor) { 936 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 937 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 938 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 939 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 940 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 941 } 942 for (k = 0; k < fas->max_down_it; k++) { 943 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 944 if (fas->monitor) { 945 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 946 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 947 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 948 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 949 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 950 } 951 } 952 } else if (snes->pc) { 953 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 954 } 955 } 956 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 957 958 PetscFunctionReturn(0); 959 } 960 961 #undef __FUNCT__ 962 #define __FUNCT__ "SNESSolve_FAS" 963 964 PetscErrorCode SNESSolve_FAS(SNES snes) 965 { 966 PetscErrorCode ierr; 967 PetscInt i, maxits; 968 Vec X, B,F; 969 PetscReal fnorm; 970 PetscFunctionBegin; 971 maxits = snes->max_its; /* maximum number of iterations */ 972 snes->reason = SNES_CONVERGED_ITERATING; 973 X = snes->vec_sol; 974 B = snes->vec_rhs; 975 F = snes->vec_func; 976 977 /*norm setup */ 978 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 979 snes->iter = 0; 980 snes->norm = 0.; 981 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 982 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 983 if (snes->domainerror) { 984 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 985 PetscFunctionReturn(0); 986 } 987 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 988 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 989 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 990 snes->norm = fnorm; 991 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 992 SNESLogConvHistory(snes,fnorm,0); 993 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 994 995 /* set parameter for default relative tolerance convergence test */ 996 snes->ttol = fnorm*snes->rtol; 997 /* test convergence */ 998 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 999 if (snes->reason) PetscFunctionReturn(0); 1000 for (i = 0; i < maxits; i++) { 1001 /* Call general purpose update function */ 1002 1003 if (snes->ops->update) { 1004 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1005 } 1006 ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 1007 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1008 /* Monitor convergence */ 1009 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1010 snes->iter = i+1; 1011 snes->norm = fnorm; 1012 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1013 SNESLogConvHistory(snes,snes->norm,0); 1014 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1015 /* Test for convergence */ 1016 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1017 if (snes->reason) break; 1018 } 1019 if (i == maxits) { 1020 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1021 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1022 } 1023 PetscFunctionReturn(0); 1024 } 1025