1421d9b32SPeter Brune /* Defines the basic SNES object */ 2421d9b32SPeter Brune #include <../src/snes/impls/fas/fasimpls.h> 3421d9b32SPeter Brune 4421d9b32SPeter Brune /*MC 5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver. 6421d9b32SPeter Brune 7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8421d9b32SPeter Brune 9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10421d9b32SPeter Brune M*/ 11421d9b32SPeter Brune 12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 18421d9b32SPeter Brune 19421d9b32SPeter Brune EXTERN_C_BEGIN 20421d9b32SPeter Brune 21421d9b32SPeter Brune #undef __FUNCT__ 22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 24421d9b32SPeter Brune { 25421d9b32SPeter Brune SNES_FAS * fas; 26421d9b32SPeter Brune PetscErrorCode ierr; 27421d9b32SPeter Brune 28421d9b32SPeter Brune PetscFunctionBegin; 29421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 30421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 31421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 33421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 34421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 35421d9b32SPeter Brune 36ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 37ed020824SBarry Smith snes->usespc = PETSC_FALSE; 38ed020824SBarry Smith 39421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 40421d9b32SPeter Brune snes->data = (void*) fas; 41421d9b32SPeter Brune fas->level = 0; 42293a7e31SPeter Brune fas->levels = 1; 43ee78dd50SPeter Brune fas->n_cycles = 1; 44ee78dd50SPeter Brune fas->max_up_it = 1; 45ee78dd50SPeter Brune fas->max_down_it = 1; 46ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 47ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 48421d9b32SPeter Brune fas->next = PETSC_NULL; 49421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 50421d9b32SPeter Brune fas->restrct = PETSC_NULL; 51efe1f98aSPeter Brune fas->inject = PETSC_NULL; 52421d9b32SPeter Brune 53421d9b32SPeter Brune PetscFunctionReturn(0); 54421d9b32SPeter Brune } 55421d9b32SPeter Brune EXTERN_C_END 56421d9b32SPeter Brune 57421d9b32SPeter Brune #undef __FUNCT__ 58421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 59421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 60421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 61421d9b32SPeter Brune PetscFunctionBegin; 62ee1fd11aSPeter Brune *levels = fas->levels; 63421d9b32SPeter Brune PetscFunctionReturn(0); 64421d9b32SPeter Brune } 65421d9b32SPeter Brune 66421d9b32SPeter Brune #undef __FUNCT__ 67421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 68421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 69421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 70421d9b32SPeter Brune PetscInt levels = fas->level; 71421d9b32SPeter Brune PetscInt i; 72421d9b32SPeter Brune PetscFunctionBegin; 73421d9b32SPeter Brune *lsnes = snes; 74421d9b32SPeter Brune if (fas->level < level) { 75421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 76421d9b32SPeter Brune } 77421d9b32SPeter Brune if (level > levels - 1) { 78421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 79421d9b32SPeter Brune } 80421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 81421d9b32SPeter Brune *lsnes = fas->next; 82421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 83421d9b32SPeter Brune } 84421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 85421d9b32SPeter Brune PetscFunctionReturn(0); 86421d9b32SPeter Brune } 87421d9b32SPeter Brune 88421d9b32SPeter Brune #undef __FUNCT__ 89421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 90421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 91421d9b32SPeter Brune PetscErrorCode ierr; 92421d9b32SPeter Brune PetscInt i; 93421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 94421d9b32SPeter Brune MPI_Comm comm; 95421d9b32SPeter Brune PetscFunctionBegin; 96ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 97ee1fd11aSPeter Brune if (levels == fas->levels) { 98ee1fd11aSPeter Brune if (!comms) 99ee1fd11aSPeter Brune PetscFunctionReturn(0); 100ee1fd11aSPeter Brune } 101421d9b32SPeter Brune /* user has changed the number of levels; reset */ 102421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 103421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 104421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 105ee1fd11aSPeter Brune fas->next = PETSC_NULL; 106421d9b32SPeter Brune /* setup the finest level */ 107421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 108421d9b32SPeter Brune if (comms) comm = comms[i]; 109421d9b32SPeter Brune fas->level = i; 110421d9b32SPeter Brune fas->levels = levels; 111ee1fd11aSPeter Brune fas->next = PETSC_NULL; 112421d9b32SPeter Brune if (i > 0) { 113421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 1144a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 115293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 116421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 117421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 118421d9b32SPeter Brune } 119421d9b32SPeter Brune } 120421d9b32SPeter Brune PetscFunctionReturn(0); 121421d9b32SPeter Brune } 122421d9b32SPeter Brune 123421d9b32SPeter Brune #undef __FUNCT__ 124421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 125421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 126421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 127d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 128421d9b32SPeter Brune 129421d9b32SPeter Brune PetscFunctionBegin; 130421d9b32SPeter Brune if (level > top_level) 131421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 132421d9b32SPeter Brune /* get to the correct level */ 133d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 134421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 135421d9b32SPeter Brune } 136421d9b32SPeter Brune if (fas->level != level) 137421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 138421d9b32SPeter Brune fas->interpolate = mat; 139421d9b32SPeter Brune PetscFunctionReturn(0); 140421d9b32SPeter Brune } 141421d9b32SPeter Brune 142421d9b32SPeter Brune #undef __FUNCT__ 143421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 144421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 145421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 146d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 147421d9b32SPeter Brune 148421d9b32SPeter Brune PetscFunctionBegin; 149421d9b32SPeter Brune if (level > top_level) 150421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 151421d9b32SPeter Brune /* get to the correct level */ 152d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 153421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 154421d9b32SPeter Brune } 155421d9b32SPeter Brune if (fas->level != level) 156421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 157421d9b32SPeter Brune fas->restrct = mat; 158421d9b32SPeter Brune PetscFunctionReturn(0); 159421d9b32SPeter Brune } 160421d9b32SPeter Brune 161efe1f98aSPeter Brune 162421d9b32SPeter Brune #undef __FUNCT__ 163421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 164421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 165421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 166d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 167421d9b32SPeter Brune 168421d9b32SPeter Brune PetscFunctionBegin; 169421d9b32SPeter Brune if (level > top_level) 170421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 171421d9b32SPeter Brune /* get to the correct level */ 172d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 173421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 174421d9b32SPeter Brune } 175421d9b32SPeter Brune if (fas->level != level) 176421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 177421d9b32SPeter Brune fas->rscale = rscale; 178421d9b32SPeter Brune PetscFunctionReturn(0); 179421d9b32SPeter Brune } 180421d9b32SPeter Brune 181efe1f98aSPeter Brune 182efe1f98aSPeter Brune #undef __FUNCT__ 183efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 184efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 185efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 186efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 187efe1f98aSPeter Brune 188efe1f98aSPeter Brune PetscFunctionBegin; 189efe1f98aSPeter Brune if (level > top_level) 190efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 191efe1f98aSPeter Brune /* get to the correct level */ 192efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 193efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 194efe1f98aSPeter Brune } 195efe1f98aSPeter Brune if (fas->level != level) 196efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 197efe1f98aSPeter Brune fas->inject = mat; 198efe1f98aSPeter Brune PetscFunctionReturn(0); 199efe1f98aSPeter Brune } 200efe1f98aSPeter Brune 201421d9b32SPeter Brune #undef __FUNCT__ 202421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 203421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 204421d9b32SPeter Brune { 20577df8cc4SPeter Brune PetscErrorCode ierr = 0; 206421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 207421d9b32SPeter Brune 208421d9b32SPeter Brune PetscFunctionBegin; 209421d9b32SPeter Brune /* destroy local data created in SNESSetup_FAS */ 21051e86f29SPeter Brune #if 0 211293a7e31SPeter Brune /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 21251e86f29SPeter Brune #endif 213ee1fd11aSPeter Brune if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr); 21451e86f29SPeter Brune #if 0 21551e86f29SPeter Brune #endif 216421d9b32SPeter Brune PetscFunctionReturn(0); 217421d9b32SPeter Brune } 218421d9b32SPeter Brune 219421d9b32SPeter Brune #undef __FUNCT__ 220421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 221421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 222421d9b32SPeter Brune { 223421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 224421d9b32SPeter Brune PetscErrorCode ierr; 225421d9b32SPeter Brune 226421d9b32SPeter Brune PetscFunctionBegin; 227421d9b32SPeter Brune /* recursively resets and then destroys */ 228421d9b32SPeter Brune ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 22951e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 23051e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 231efe1f98aSPeter Brune if (fas->inject) { 232efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 233efe1f98aSPeter Brune } 23451e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 23551e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 23651e86f29SPeter Brune fas->restrct = PETSC_NULL; 23751e86f29SPeter Brune } else { 23851e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 23951e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 24051e86f29SPeter Brune } 24151e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 24251e86f29SPeter Brune if (snes->work) ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 243421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 244421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 245ee1fd11aSPeter Brune 246421d9b32SPeter Brune PetscFunctionReturn(0); 247421d9b32SPeter Brune } 248421d9b32SPeter Brune 249421d9b32SPeter Brune #undef __FUNCT__ 250421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 251421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 252421d9b32SPeter Brune { 2531a266240SBarry Smith SNES_FAS *fas = (SNES_FAS *) snes->data,*tmp; 254421d9b32SPeter Brune PetscErrorCode ierr; 255efe1f98aSPeter Brune VecScatter injscatter; 256421d9b32SPeter Brune 257421d9b32SPeter Brune PetscFunctionBegin; 2581a266240SBarry Smith /* should call the SNESSetFromOptions() only when approriate */ 2591a266240SBarry Smith tmp = fas; 2601a266240SBarry Smith while (tmp) { 261*6bed07a3SBarry Smith if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);} 262*6bed07a3SBarry Smith if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);} 2631a266240SBarry Smith tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0; 2641a266240SBarry Smith } 2651a266240SBarry Smith 2661a266240SBarry Smith if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */ 267421d9b32SPeter Brune /* gets the solver ready for solution */ 268421d9b32SPeter Brune if (snes->dm) { 269421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 270421d9b32SPeter Brune if (fas->next) { 271421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 272ee78dd50SPeter Brune if (!fas->next->dm) { 273ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 274ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 275ee78dd50SPeter Brune } 276421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 277ee78dd50SPeter Brune if (!fas->interpolate) { 278421d9b32SPeter Brune ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 279421d9b32SPeter Brune fas->restrct = fas->interpolate; 280421d9b32SPeter Brune } 281efe1f98aSPeter Brune /* set the injection from the DM */ 282efe1f98aSPeter Brune if (!fas->inject) { 283efe1f98aSPeter Brune ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 284efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 285efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 286efe1f98aSPeter Brune } 287421d9b32SPeter Brune } 288ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 289*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 290*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 291421d9b32SPeter Brune } 292fa9694d7SPeter Brune if (fas->next) { 293fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 294*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 295fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 296fa9694d7SPeter Brune } 297fa9694d7SPeter Brune /* got to set them all up at once */ 298421d9b32SPeter Brune PetscFunctionReturn(0); 299421d9b32SPeter Brune } 300421d9b32SPeter Brune 301421d9b32SPeter Brune #undef __FUNCT__ 302421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 303421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 304421d9b32SPeter Brune { 305ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 306ee78dd50SPeter Brune PetscInt levels = 1; 307ee78dd50SPeter Brune PetscBool flg, monflg; 308421d9b32SPeter Brune PetscErrorCode ierr; 309ee78dd50SPeter Brune const char * def_smooth = SNESNRICHARDSON; 310ee78dd50SPeter Brune char pre_type[256]; 311ee78dd50SPeter Brune char post_type[256]; 312ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 313421d9b32SPeter Brune 314421d9b32SPeter Brune PetscFunctionBegin; 315c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 316ee78dd50SPeter Brune 317ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 318ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 319ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 320c732cbdbSBarry Smith if (!flg && snes->dm) { 321c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 322c732cbdbSBarry Smith levels++; 323c732cbdbSBarry Smith } 324ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 325ee78dd50SPeter Brune } 326ee78dd50SPeter Brune 327ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 328ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 329ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 330ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 331ee78dd50SPeter Brune if (flg) { 332ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 333ee78dd50SPeter Brune } else { 334ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 335ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 336ee78dd50SPeter Brune } 337ee78dd50SPeter Brune 338ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 339ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 340ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 341ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 342ee78dd50SPeter Brune 343ee78dd50SPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 344ee78dd50SPeter Brune 345ee78dd50SPeter Brune /* other options for the coarsest level */ 346ee78dd50SPeter Brune if (fas->level == 0) { 347ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 348ee78dd50SPeter Brune if (flg) { 349ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 350ee78dd50SPeter Brune } else { 351ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 352ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 353ee78dd50SPeter Brune } 354ee78dd50SPeter Brune } 355ee78dd50SPeter Brune 356421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 357ee78dd50SPeter Brune /* setup from the determined types if the smoothers don't exist */ 358ee78dd50SPeter Brune if (!fas->upsmooth) { 35967339d5cSBarry Smith const char *prefix; 36067339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 361ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 36267339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 3631a266240SBarry Smith if (fas->level || (fas->levels == 1)) { 36467339d5cSBarry Smith ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr); 36567339d5cSBarry Smith } else { 36667339d5cSBarry Smith ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_coarse_");CHKERRQ(ierr); 36767339d5cSBarry Smith } 368293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 369ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 37069bed9afSBarry Smith if (snes->ops->computefunction) { 37169bed9afSBarry Smith ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 37269bed9afSBarry Smith } 373ee78dd50SPeter Brune } 3741a266240SBarry Smith if (fas->upsmooth) { 3751a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 3761a266240SBarry Smith } 3771a266240SBarry Smith 378293a7e31SPeter Brune if (!fas->downsmooth && fas->level != 0) { 37967339d5cSBarry Smith const char *prefix; 38067339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 381ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 38267339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 38367339d5cSBarry Smith ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr); 384293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 385ee78dd50SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 38669bed9afSBarry Smith if (snes->ops->computefunction) { 38769bed9afSBarry Smith ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 38869bed9afSBarry Smith } 389ee78dd50SPeter Brune } 3901a266240SBarry Smith if (fas->downsmooth) { 3911a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 3921a266240SBarry Smith } 393ee78dd50SPeter Brune 394ee78dd50SPeter Brune if (monflg) { 395293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 396293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 397ee78dd50SPeter Brune } 398ee78dd50SPeter Brune 399ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 4001a266240SBarry Smith if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 401421d9b32SPeter Brune PetscFunctionReturn(0); 402421d9b32SPeter Brune } 403421d9b32SPeter Brune 404421d9b32SPeter Brune #undef __FUNCT__ 405421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 406421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 407421d9b32SPeter Brune { 408421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 409421d9b32SPeter Brune PetscBool iascii; 410421d9b32SPeter Brune PetscErrorCode ierr; 411421d9b32SPeter Brune PetscInt levels = fas->levels; 412421d9b32SPeter Brune PetscInt i; 413421d9b32SPeter Brune 414421d9b32SPeter Brune PetscFunctionBegin; 415421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 416421d9b32SPeter Brune if (iascii) { 417421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 418421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 419421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 420421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 421ee78dd50SPeter Brune if (fas->upsmooth) { 422ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 423421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 424ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 425421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 426421d9b32SPeter Brune } else { 427ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 428421d9b32SPeter Brune } 429ee78dd50SPeter Brune if (fas->downsmooth) { 430ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 431421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 432ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 433421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 434421d9b32SPeter Brune } else { 435ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 436421d9b32SPeter Brune } 437421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 438421d9b32SPeter Brune } 439421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 440421d9b32SPeter Brune } else { 441421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 442421d9b32SPeter Brune } 443421d9b32SPeter Brune PetscFunctionReturn(0); 444421d9b32SPeter Brune } 445421d9b32SPeter Brune 446421d9b32SPeter Brune /* 447fa9694d7SPeter Brune 448fa9694d7SPeter Brune Defines the FAS cycle as: 449fa9694d7SPeter Brune 450fa9694d7SPeter Brune fine problem: F(x) = 0 451fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 452fa9694d7SPeter Brune 453fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 454fa9694d7SPeter Brune 455fa9694d7SPeter Brune correction: 456fa9694d7SPeter Brune 457fa9694d7SPeter Brune x = x + I(x^c - Rx) 458fa9694d7SPeter Brune 459421d9b32SPeter Brune */ 460fa9694d7SPeter Brune 461421d9b32SPeter Brune #undef __FUNCT__ 462421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 463f5a6d4f9SBarry Smith PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X) { 464fa9694d7SPeter Brune 465fa9694d7SPeter Brune PetscErrorCode ierr; 466f5a6d4f9SBarry Smith Vec X_c, Xo_c, F_c, B_c,F; 467fa9694d7SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 468293a7e31SPeter Brune PetscInt i; 469421d9b32SPeter Brune 470421d9b32SPeter Brune PetscFunctionBegin; 471f5a6d4f9SBarry Smith F = snes->vec_func; 472fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 473ee78dd50SPeter Brune if (fas->upsmooth) { 474ee78dd50SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 475c90fad12SPeter Brune } else if (snes->pc) { 476c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 477fe6f9142SPeter Brune } 478fa9694d7SPeter Brune if (fas->next) { 479293a7e31SPeter Brune for (i = 0; i < fas->n_cycles; i++) { 480ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 481c90fad12SPeter Brune X_c = fas->next->vec_sol; 482293a7e31SPeter Brune Xo_c = fas->next->work[0]; 483c90fad12SPeter Brune F_c = fas->next->vec_func; 484293a7e31SPeter Brune B_c = fas->next->work[1]; 485efe1f98aSPeter Brune 486efe1f98aSPeter Brune /* inject the solution */ 487efe1f98aSPeter Brune if (fas->inject) { 488a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 489efe1f98aSPeter Brune } else { 490a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 491a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 492efe1f98aSPeter Brune } 493293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 494293a7e31SPeter Brune 495293a7e31SPeter Brune /* restrict the defect */ 496293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 497293a7e31SPeter Brune 498c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 499ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 500c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 501293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 502c90fad12SPeter Brune 503ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 504ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 505c90fad12SPeter Brune 506c90fad12SPeter Brune /* recurse to the next level */ 507f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 508f5a6d4f9SBarry Smith ierr = FASCycle_Private(fas->next, B_c, X_c);CHKERRQ(ierr); 509f5a6d4f9SBarry Smith fas->next->vec_rhs = PETSC_NULL; 510ee78dd50SPeter Brune 511fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 512fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 513ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 514293a7e31SPeter Brune } 515fa9694d7SPeter Brune } 516ee78dd50SPeter Brune /* down-smooth -- just update using the down-smoother */ 517c90fad12SPeter Brune if (fas->level != 0) { 518ee78dd50SPeter Brune if (fas->downsmooth) { 519ee78dd50SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 520c90fad12SPeter Brune } else if (snes->pc) { 521c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 522c90fad12SPeter Brune } 523fe6f9142SPeter Brune } 524fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 525fa9694d7SPeter Brune 526fa9694d7SPeter Brune PetscFunctionReturn(0); 527421d9b32SPeter Brune } 528421d9b32SPeter Brune 529421d9b32SPeter Brune #undef __FUNCT__ 530c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private" 531293a7e31SPeter Brune PetscErrorCode FASInitialGuess_Private(SNES snes, Vec B, Vec X) { 532c90fad12SPeter Brune 533293a7e31SPeter Brune PetscErrorCode ierr; 534293a7e31SPeter Brune Vec X_c, B_c; 535ee1fd11aSPeter Brune SNES_FAS * fas; 536293a7e31SPeter Brune 537293a7e31SPeter Brune PetscFunctionBegin; 538ee1fd11aSPeter Brune fas = (SNES_FAS *)snes->data; 539293a7e31SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 540293a7e31SPeter Brune if (fas->level == 0) { 541293a7e31SPeter Brune if (fas->upsmooth) { 542293a7e31SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 543293a7e31SPeter Brune } else if (snes->pc) { 544293a7e31SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 545293a7e31SPeter Brune } 546293a7e31SPeter Brune } 547293a7e31SPeter Brune if (fas->next) { 548293a7e31SPeter Brune X_c = fas->next->vec_sol; 549293a7e31SPeter Brune B_c = fas->next->work[0]; 550293a7e31SPeter Brune /* inject the solution to coarse */ 551efe1f98aSPeter Brune if (fas->inject) { 552efe1f98aSPeter Brune ierr = MatRestrict(fas->inject, X, X_c);CHKERRQ(ierr); 553efe1f98aSPeter Brune } else { 554293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr); 555293a7e31SPeter Brune ierr = VecPointwiseMult(X_c, fas->rscale, X_c);CHKERRQ(ierr); 556efe1f98aSPeter Brune } 557293a7e31SPeter Brune if (B) { 558293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr); 559293a7e31SPeter Brune ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr); 560293a7e31SPeter Brune } else { 561293a7e31SPeter Brune B_c = PETSC_NULL; 562293a7e31SPeter Brune } 563293a7e31SPeter Brune /* recurse to the next level */ 564293a7e31SPeter Brune ierr = FASInitialGuess_Private(fas->next, B_c, X_c);CHKERRQ(ierr); 565293a7e31SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, X);CHKERRQ(ierr); 566293a7e31SPeter Brune } 567293a7e31SPeter Brune /* down-smooth -- just update using the down-smoother */ 568293a7e31SPeter Brune if (fas->level != 0) { 569293a7e31SPeter Brune if (fas->downsmooth) { 570293a7e31SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 571293a7e31SPeter Brune } else if (snes->pc) { 572293a7e31SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 573293a7e31SPeter Brune } 574293a7e31SPeter Brune } 575293a7e31SPeter Brune PetscFunctionReturn(0); 576293a7e31SPeter Brune } 577c90fad12SPeter Brune 578c90fad12SPeter Brune #undef __FUNCT__ 579421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 580421d9b32SPeter Brune 581421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 582421d9b32SPeter Brune { 583fa9694d7SPeter Brune PetscErrorCode ierr; 584fe6f9142SPeter Brune PetscInt i, maxits; 585f5a6d4f9SBarry Smith Vec X, B,F; 586fe6f9142SPeter Brune PetscReal fnorm; 587421d9b32SPeter Brune PetscFunctionBegin; 588fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 589fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 590fa9694d7SPeter Brune X = snes->vec_sol; 591fe6f9142SPeter Brune B = snes->vec_rhs; 592f5a6d4f9SBarry Smith F = snes->vec_func; 593293a7e31SPeter Brune 594293a7e31SPeter Brune /*norm setup */ 595fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 596fe6f9142SPeter Brune snes->iter = 0; 597fe6f9142SPeter Brune snes->norm = 0.; 598fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 599fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 600fe6f9142SPeter Brune if (snes->domainerror) { 601fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 602fe6f9142SPeter Brune PetscFunctionReturn(0); 603fe6f9142SPeter Brune } 604fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 605fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 606fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 607fe6f9142SPeter Brune snes->norm = fnorm; 608fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 609fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 610fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 611fe6f9142SPeter Brune 612fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 613fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 614fe6f9142SPeter Brune /* test convergence */ 615fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 616fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 617fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 618fe6f9142SPeter Brune /* Call general purpose update function */ 619fe6f9142SPeter Brune if (snes->ops->update) { 620fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 621fe6f9142SPeter Brune } 622f5a6d4f9SBarry Smith ierr = FASCycle_Private(snes, B, X);CHKERRQ(ierr); 623c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 624c90fad12SPeter Brune /* Monitor convergence */ 625c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 626c90fad12SPeter Brune snes->iter = i+1; 627c90fad12SPeter Brune snes->norm = fnorm; 628c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 629c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 630c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 631c90fad12SPeter Brune /* Test for convergence */ 632c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 633c90fad12SPeter Brune if (snes->reason) break; 634fe6f9142SPeter Brune } 635fe6f9142SPeter Brune if (i == maxits) { 636fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 637fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 638fe6f9142SPeter Brune } 639421d9b32SPeter Brune PetscFunctionReturn(0); 640421d9b32SPeter Brune } 641