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 36421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 37421d9b32SPeter Brune snes->data = (void*) fas; 38421d9b32SPeter Brune fas->level = 0; 39*293a7e31SPeter Brune fas->levels = 1; 40ee78dd50SPeter Brune fas->n_cycles = 1; 41ee78dd50SPeter Brune fas->max_up_it = 1; 42ee78dd50SPeter Brune fas->max_down_it = 1; 43ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 44ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 45421d9b32SPeter Brune fas->next = PETSC_NULL; 46421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 47421d9b32SPeter Brune fas->restrct = PETSC_NULL; 48421d9b32SPeter Brune 49421d9b32SPeter Brune PetscFunctionReturn(0); 50421d9b32SPeter Brune } 51421d9b32SPeter Brune EXTERN_C_END 52421d9b32SPeter Brune 53421d9b32SPeter Brune #undef __FUNCT__ 54421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 55421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 56421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 57421d9b32SPeter Brune PetscFunctionBegin; 58421d9b32SPeter Brune *levels = fas->level; 59421d9b32SPeter Brune PetscFunctionReturn(0); 60421d9b32SPeter Brune } 61421d9b32SPeter Brune 62421d9b32SPeter Brune #undef __FUNCT__ 63421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 64421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 65421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 66421d9b32SPeter Brune PetscInt levels = fas->level; 67421d9b32SPeter Brune PetscInt i; 68421d9b32SPeter Brune PetscFunctionBegin; 69421d9b32SPeter Brune *lsnes = snes; 70421d9b32SPeter Brune if (fas->level < level) { 71421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 72421d9b32SPeter Brune } 73421d9b32SPeter Brune if (level > levels - 1) { 74421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 75421d9b32SPeter Brune } 76421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 77421d9b32SPeter Brune *lsnes = fas->next; 78421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 79421d9b32SPeter Brune } 80421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 81421d9b32SPeter Brune PetscFunctionReturn(0); 82421d9b32SPeter Brune } 83421d9b32SPeter Brune 84421d9b32SPeter Brune #undef __FUNCT__ 85421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 86421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 87421d9b32SPeter Brune PetscErrorCode ierr; 88421d9b32SPeter Brune PetscInt i; 89421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 90421d9b32SPeter Brune MPI_Comm comm; 91421d9b32SPeter Brune PetscFunctionBegin; 92421d9b32SPeter Brune comm = PETSC_COMM_WORLD; 93421d9b32SPeter Brune /* user has changed the number of levels; reset */ 94421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 95421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 96421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 97421d9b32SPeter Brune /* setup the finest level */ 98421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 99421d9b32SPeter Brune if (comms) comm = comms[i]; 100421d9b32SPeter Brune fas->level = i; 101421d9b32SPeter Brune fas->levels = levels; 102421d9b32SPeter Brune if (i > 0) { 103421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 104*293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 105421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 106421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 107421d9b32SPeter Brune } 108421d9b32SPeter Brune } 109421d9b32SPeter Brune PetscFunctionReturn(0); 110421d9b32SPeter Brune } 111421d9b32SPeter Brune 112421d9b32SPeter Brune #undef __FUNCT__ 113421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 114421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 115421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 116d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 117421d9b32SPeter Brune 118421d9b32SPeter Brune PetscFunctionBegin; 119421d9b32SPeter Brune if (level > top_level) 120421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 121421d9b32SPeter Brune /* get to the correct level */ 122d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 123421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 124421d9b32SPeter Brune } 125421d9b32SPeter Brune if (fas->level != level) 126421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 127421d9b32SPeter Brune fas->interpolate = mat; 128421d9b32SPeter Brune PetscFunctionReturn(0); 129421d9b32SPeter Brune } 130421d9b32SPeter Brune 131421d9b32SPeter Brune #undef __FUNCT__ 132421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 133421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 134421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 135d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 136421d9b32SPeter Brune 137421d9b32SPeter Brune PetscFunctionBegin; 138421d9b32SPeter Brune if (level > top_level) 139421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 140421d9b32SPeter Brune /* get to the correct level */ 141d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 142421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 143421d9b32SPeter Brune } 144421d9b32SPeter Brune if (fas->level != level) 145421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 146421d9b32SPeter Brune fas->restrct = mat; 147421d9b32SPeter Brune PetscFunctionReturn(0); 148421d9b32SPeter Brune } 149421d9b32SPeter Brune 150421d9b32SPeter Brune #undef __FUNCT__ 151421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 152421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 153421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 154d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 155421d9b32SPeter Brune 156421d9b32SPeter Brune PetscFunctionBegin; 157421d9b32SPeter Brune if (level > top_level) 158421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 159421d9b32SPeter Brune /* get to the correct level */ 160d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 161421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 162421d9b32SPeter Brune } 163421d9b32SPeter Brune if (fas->level != level) 164421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 165421d9b32SPeter Brune fas->rscale = rscale; 166421d9b32SPeter Brune PetscFunctionReturn(0); 167421d9b32SPeter Brune } 168421d9b32SPeter Brune 169421d9b32SPeter Brune #undef __FUNCT__ 170421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 171421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 172421d9b32SPeter Brune { 17377df8cc4SPeter Brune PetscErrorCode ierr = 0; 174421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 175421d9b32SPeter Brune 176421d9b32SPeter Brune PetscFunctionBegin; 177421d9b32SPeter Brune /* destroy local data created in SNESSetup_FAS */ 178ee78dd50SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 179ee78dd50SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 180*293a7e31SPeter Brune if (fas->interpolate == fas->restrct) { 181*293a7e31SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 182*293a7e31SPeter Brune fas->restrct = PETSC_NULL; 183*293a7e31SPeter Brune } else { 184421d9b32SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 185421d9b32SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 186*293a7e31SPeter Brune } 187421d9b32SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 188421d9b32SPeter Brune 189*293a7e31SPeter Brune /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 190*293a7e31SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 191*293a7e31SPeter Brune 192fa9694d7SPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 193421d9b32SPeter Brune PetscFunctionReturn(0); 194421d9b32SPeter Brune } 195421d9b32SPeter Brune 196421d9b32SPeter Brune #undef __FUNCT__ 197421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 198421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 199421d9b32SPeter Brune { 200421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 201421d9b32SPeter Brune PetscErrorCode ierr; 202421d9b32SPeter Brune 203421d9b32SPeter Brune PetscFunctionBegin; 204421d9b32SPeter Brune /* recursively resets and then destroys */ 205421d9b32SPeter Brune ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 206421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 207421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 208421d9b32SPeter Brune PetscFunctionReturn(0); 209421d9b32SPeter Brune } 210421d9b32SPeter Brune 211421d9b32SPeter Brune #undef __FUNCT__ 212421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 213421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 214421d9b32SPeter Brune { 215421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 216421d9b32SPeter Brune PetscErrorCode ierr; 217421d9b32SPeter Brune 218421d9b32SPeter Brune PetscFunctionBegin; 219*293a7e31SPeter Brune ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 220421d9b32SPeter Brune /* gets the solver ready for solution */ 221421d9b32SPeter Brune if (snes->dm) { 222421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 223421d9b32SPeter Brune if (fas->next) { 224421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 225ee78dd50SPeter Brune if (!fas->next->dm) { 226ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 227ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 228ee78dd50SPeter Brune } 229421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 230ee78dd50SPeter Brune if (!fas->interpolate) { 231421d9b32SPeter Brune ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 232421d9b32SPeter Brune fas->restrct = fas->interpolate; 233421d9b32SPeter Brune } 234421d9b32SPeter Brune } 235ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 236ee78dd50SPeter Brune if (fas->upsmooth) SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr); 237ee78dd50SPeter Brune if (fas->downsmooth)SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr); 238421d9b32SPeter Brune } 239fa9694d7SPeter Brune if (fas->next) { 240fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 241fa9694d7SPeter Brune ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr); 242fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 243fa9694d7SPeter Brune } 244fa9694d7SPeter Brune /* got to set them all up at once */ 245421d9b32SPeter Brune PetscFunctionReturn(0); 246421d9b32SPeter Brune } 247421d9b32SPeter Brune 248421d9b32SPeter Brune #undef __FUNCT__ 249421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 250421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 251421d9b32SPeter Brune { 252ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 253ee78dd50SPeter Brune PetscInt levels = 1; 254ee78dd50SPeter Brune PetscBool flg, monflg; 255421d9b32SPeter Brune PetscErrorCode ierr; 256ee78dd50SPeter Brune const char * def_smooth = SNESNRICHARDSON; 257ee78dd50SPeter Brune char pre_type[256]; 258ee78dd50SPeter Brune char post_type[256]; 259ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 260421d9b32SPeter Brune 261421d9b32SPeter Brune PetscFunctionBegin; 262c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 263ee78dd50SPeter Brune 264ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 265ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 266ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 267ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 268ee78dd50SPeter Brune } 269ee78dd50SPeter Brune 270ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 271ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 272ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 273ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 274ee78dd50SPeter Brune if (flg) { 275ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 276ee78dd50SPeter Brune } else { 277ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 278ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 279ee78dd50SPeter Brune } 280ee78dd50SPeter Brune 281ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 282ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 283ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 284ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 285ee78dd50SPeter Brune 286ee78dd50SPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 287ee78dd50SPeter Brune 288ee78dd50SPeter Brune /* other options for the coarsest level */ 289ee78dd50SPeter Brune if (fas->level == 0) { 290ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 291ee78dd50SPeter Brune if (flg) { 292ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 293ee78dd50SPeter Brune } else { 294ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 295ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 296ee78dd50SPeter Brune } 297ee78dd50SPeter Brune } 298ee78dd50SPeter Brune 299421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 300ee78dd50SPeter Brune /* setup from the determined types if the smoothers don't exist */ 301ee78dd50SPeter Brune if (!fas->upsmooth) { 302ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 303*293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 304ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 305ee78dd50SPeter Brune } 306*293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 307*293a7e31SPeter Brune if (!fas->downsmooth && fas->level != 0) { 308ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 309*293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 310ee78dd50SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 311ee78dd50SPeter Brune } 312*293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 313ee78dd50SPeter Brune 314ee78dd50SPeter Brune if (monflg) { 315*293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 316*293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 317ee78dd50SPeter Brune } 318ee78dd50SPeter Brune 319ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 320ee78dd50SPeter Brune if (fas->next)ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr); 321421d9b32SPeter Brune PetscFunctionReturn(0); 322421d9b32SPeter Brune } 323421d9b32SPeter Brune 324421d9b32SPeter Brune #undef __FUNCT__ 325421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 326421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 327421d9b32SPeter Brune { 328421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 329421d9b32SPeter Brune PetscBool iascii; 330421d9b32SPeter Brune PetscErrorCode ierr; 331421d9b32SPeter Brune PetscInt levels = fas->levels; 332421d9b32SPeter Brune PetscInt i; 333421d9b32SPeter Brune 334421d9b32SPeter Brune PetscFunctionBegin; 335421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 336421d9b32SPeter Brune if (iascii) { 337421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 338421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 339421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 340421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 341ee78dd50SPeter Brune if (fas->upsmooth) { 342ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 343421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 344ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 345421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 346421d9b32SPeter Brune } else { 347ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 348421d9b32SPeter Brune } 349ee78dd50SPeter Brune if (fas->downsmooth) { 350ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 351421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 352ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 353421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 354421d9b32SPeter Brune } else { 355ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 356421d9b32SPeter Brune } 357421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 358421d9b32SPeter Brune } 359421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 360421d9b32SPeter Brune } else { 361421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 362421d9b32SPeter Brune } 363421d9b32SPeter Brune PetscFunctionReturn(0); 364421d9b32SPeter Brune } 365421d9b32SPeter Brune 366421d9b32SPeter Brune /* 367fa9694d7SPeter Brune 368fa9694d7SPeter Brune Defines the FAS cycle as: 369fa9694d7SPeter Brune 370fa9694d7SPeter Brune fine problem: F(x) = 0 371fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 372fa9694d7SPeter Brune 373fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 374fa9694d7SPeter Brune 375fa9694d7SPeter Brune correction: 376fa9694d7SPeter Brune 377fa9694d7SPeter Brune x = x + I(x^c - Rx) 378fa9694d7SPeter Brune 379421d9b32SPeter Brune */ 380fa9694d7SPeter Brune 381421d9b32SPeter Brune #undef __FUNCT__ 382421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 383c90fad12SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) { 384fa9694d7SPeter Brune 385fa9694d7SPeter Brune PetscErrorCode ierr; 386fa9694d7SPeter Brune Vec X_c, Xo_c, F_c, B_c; 387fa9694d7SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 388*293a7e31SPeter Brune PetscInt i; 389421d9b32SPeter Brune 390421d9b32SPeter Brune PetscFunctionBegin; 391fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 392ee78dd50SPeter Brune if (fas->upsmooth) { 393ee78dd50SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 394c90fad12SPeter Brune } else if (snes->pc) { 395c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 396fe6f9142SPeter Brune } 397fa9694d7SPeter Brune if (fas->next) { 398*293a7e31SPeter Brune for (i = 0; i < fas->n_cycles; i++) { 399ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 400c90fad12SPeter Brune X_c = fas->next->vec_sol; 401*293a7e31SPeter Brune Xo_c = fas->next->work[0]; 402c90fad12SPeter Brune F_c = fas->next->vec_func; 403*293a7e31SPeter Brune B_c = fas->next->work[1]; 404*293a7e31SPeter Brune /* inject the solution to coarse */ 405c90fad12SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 406fa9694d7SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 407c90fad12SPeter Brune if (B) { 408*293a7e31SPeter Brune ierr = VecAYPX(F, -1.0, B);CHKERRQ(ierr); /* F = B - F (defect) */ 409c90fad12SPeter Brune } else { 410*293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 411c90fad12SPeter Brune } 412*293a7e31SPeter Brune 413*293a7e31SPeter Brune /* restrict the defect */ 414*293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 415*293a7e31SPeter Brune ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr); 416*293a7e31SPeter Brune 417c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 418ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 419c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 420*293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 421c90fad12SPeter Brune 422ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 423ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 424c90fad12SPeter Brune 425c90fad12SPeter Brune /* recurse to the next level */ 426c90fad12SPeter Brune ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr); 427ee78dd50SPeter Brune 428fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 429fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 430ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 431*293a7e31SPeter Brune } 432fa9694d7SPeter Brune } 433ee78dd50SPeter Brune /* down-smooth -- just update using the down-smoother */ 434c90fad12SPeter Brune if (fas->level != 0) { 435ee78dd50SPeter Brune if (fas->downsmooth) { 436ee78dd50SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 437c90fad12SPeter Brune } else if (snes->pc) { 438c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 439c90fad12SPeter Brune } 440fe6f9142SPeter Brune } 441fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 442fa9694d7SPeter Brune 443fa9694d7SPeter Brune PetscFunctionReturn(0); 444421d9b32SPeter Brune } 445421d9b32SPeter Brune 446421d9b32SPeter Brune #undef __FUNCT__ 447c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private" 448*293a7e31SPeter Brune PetscErrorCode FASInitialGuess_Private(SNES snes, Vec B, Vec X) { 449c90fad12SPeter Brune 450*293a7e31SPeter Brune PetscErrorCode ierr; 451*293a7e31SPeter Brune Vec X_c, B_c; 452*293a7e31SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 453*293a7e31SPeter Brune 454*293a7e31SPeter Brune PetscFunctionBegin; 455*293a7e31SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 456*293a7e31SPeter Brune if (fas->level == 0) { 457*293a7e31SPeter Brune if (fas->upsmooth) { 458*293a7e31SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 459*293a7e31SPeter Brune } else if (snes->pc) { 460*293a7e31SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 461*293a7e31SPeter Brune } 462*293a7e31SPeter Brune } 463*293a7e31SPeter Brune if (fas->next) { 464*293a7e31SPeter Brune X_c = fas->next->vec_sol; 465*293a7e31SPeter Brune B_c = fas->next->work[0]; 466*293a7e31SPeter Brune /* inject the solution to coarse */ 467*293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr); 468*293a7e31SPeter Brune ierr = VecPointwiseMult(X_c, fas->rscale, X_c);CHKERRQ(ierr); 469*293a7e31SPeter Brune if (B) { 470*293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr); 471*293a7e31SPeter Brune ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr); 472*293a7e31SPeter Brune } else { 473*293a7e31SPeter Brune B_c = PETSC_NULL; 474*293a7e31SPeter Brune } 475*293a7e31SPeter Brune /* recurse to the next level */ 476*293a7e31SPeter Brune ierr = FASInitialGuess_Private(fas->next, B_c, X_c);CHKERRQ(ierr); 477*293a7e31SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, X);CHKERRQ(ierr); 478*293a7e31SPeter Brune } 479*293a7e31SPeter Brune /* down-smooth -- just update using the down-smoother */ 480*293a7e31SPeter Brune if (fas->level != 0) { 481*293a7e31SPeter Brune if (fas->downsmooth) { 482*293a7e31SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 483*293a7e31SPeter Brune } else if (snes->pc) { 484*293a7e31SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 485*293a7e31SPeter Brune } 486*293a7e31SPeter Brune } 487*293a7e31SPeter Brune PetscFunctionReturn(0); 488*293a7e31SPeter Brune } 489c90fad12SPeter Brune 490c90fad12SPeter Brune #undef __FUNCT__ 491421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 492421d9b32SPeter Brune 493421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 494421d9b32SPeter Brune { 495fa9694d7SPeter Brune PetscErrorCode ierr; 496fe6f9142SPeter Brune PetscInt i, maxits; 497fe6f9142SPeter Brune Vec X, F, B; 498fe6f9142SPeter Brune PetscReal fnorm; 499421d9b32SPeter Brune PetscFunctionBegin; 500fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 501fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 502fa9694d7SPeter Brune X = snes->vec_sol; 503fa9694d7SPeter Brune F = snes->vec_func; 504fe6f9142SPeter Brune B = snes->vec_rhs; 505fe6f9142SPeter Brune /* initial iteration */ 506*293a7e31SPeter Brune ierr = FASInitialGuess_Private(snes, B, X);CHKERRQ(ierr); 507*293a7e31SPeter Brune 508*293a7e31SPeter Brune /*norm setup */ 509fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 510fe6f9142SPeter Brune snes->iter = 0; 511fe6f9142SPeter Brune snes->norm = 0.; 512fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 513fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 514fe6f9142SPeter Brune if (snes->domainerror) { 515fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 516fe6f9142SPeter Brune PetscFunctionReturn(0); 517fe6f9142SPeter Brune } 518fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 519fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 520fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 521fe6f9142SPeter Brune snes->norm = fnorm; 522fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 523fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 524fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 525fe6f9142SPeter Brune 526fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 527fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 528fe6f9142SPeter Brune /* test convergence */ 529fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 530fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 531fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 532fe6f9142SPeter Brune /* Call general purpose update function */ 533fe6f9142SPeter Brune if (snes->ops->update) { 534fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 535fe6f9142SPeter Brune } 536c90fad12SPeter Brune ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr); 537c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 538c90fad12SPeter Brune /* Monitor convergence */ 539c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 540c90fad12SPeter Brune snes->iter = i+1; 541c90fad12SPeter Brune snes->norm = fnorm; 542c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 543c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 544c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 545c90fad12SPeter Brune /* Test for convergence */ 546c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 547c90fad12SPeter Brune if (snes->reason) break; 548fe6f9142SPeter Brune } 549fe6f9142SPeter Brune if (i == maxits) { 550fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 551fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 552fe6f9142SPeter Brune } 553421d9b32SPeter Brune PetscFunctionReturn(0); 554421d9b32SPeter Brune } 555