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*ee78dd50SPeter Brune fas->n_cycles = 1; 40*ee78dd50SPeter Brune fas->max_up_it = 1; 41*ee78dd50SPeter Brune fas->max_down_it = 1; 42*ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 43*ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 44421d9b32SPeter Brune fas->next = PETSC_NULL; 45421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 46421d9b32SPeter Brune fas->restrct = PETSC_NULL; 47421d9b32SPeter Brune 48421d9b32SPeter Brune PetscFunctionReturn(0); 49421d9b32SPeter Brune } 50421d9b32SPeter Brune EXTERN_C_END 51421d9b32SPeter Brune 52421d9b32SPeter Brune #undef __FUNCT__ 53421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 54421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 55421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 56421d9b32SPeter Brune PetscFunctionBegin; 57421d9b32SPeter Brune *levels = fas->level; 58421d9b32SPeter Brune PetscFunctionReturn(0); 59421d9b32SPeter Brune } 60421d9b32SPeter Brune 61421d9b32SPeter Brune #undef __FUNCT__ 62421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 63421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 64421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 65421d9b32SPeter Brune PetscInt levels = fas->level; 66421d9b32SPeter Brune PetscInt i; 67421d9b32SPeter Brune PetscFunctionBegin; 68421d9b32SPeter Brune *lsnes = snes; 69421d9b32SPeter Brune if (fas->level < level) { 70421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 71421d9b32SPeter Brune } 72421d9b32SPeter Brune if (level > levels - 1) { 73421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 74421d9b32SPeter Brune } 75421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 76421d9b32SPeter Brune *lsnes = fas->next; 77421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 78421d9b32SPeter Brune } 79421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 80421d9b32SPeter Brune PetscFunctionReturn(0); 81421d9b32SPeter Brune } 82421d9b32SPeter Brune 83421d9b32SPeter Brune #undef __FUNCT__ 84421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 85421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 86421d9b32SPeter Brune PetscErrorCode ierr; 87421d9b32SPeter Brune PetscInt i; 88421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 89421d9b32SPeter Brune MPI_Comm comm; 90*ee78dd50SPeter Brune SNES so; 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 */ 98*ee78dd50SPeter Brune so = snes; 99421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 100421d9b32SPeter Brune if (comms) comm = comms[i]; 101421d9b32SPeter Brune fas->level = i; 102421d9b32SPeter Brune fas->levels = levels; 103421d9b32SPeter Brune if (i > 0) { 104421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 105*ee78dd50SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)so, 1);CHKERRQ(ierr); 106421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 107421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 108*ee78dd50SPeter Brune so = fas->next; 109421d9b32SPeter Brune } 110421d9b32SPeter Brune } 111421d9b32SPeter Brune PetscFunctionReturn(0); 112421d9b32SPeter Brune } 113421d9b32SPeter Brune 114421d9b32SPeter Brune #undef __FUNCT__ 115421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 116421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 117421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 118d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 119421d9b32SPeter Brune 120421d9b32SPeter Brune PetscFunctionBegin; 121421d9b32SPeter Brune if (level > top_level) 122421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 123421d9b32SPeter Brune /* get to the correct level */ 124d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 125421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 126421d9b32SPeter Brune } 127421d9b32SPeter Brune if (fas->level != level) 128421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 129421d9b32SPeter Brune fas->interpolate = mat; 130421d9b32SPeter Brune PetscFunctionReturn(0); 131421d9b32SPeter Brune } 132421d9b32SPeter Brune 133421d9b32SPeter Brune #undef __FUNCT__ 134421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 135421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 136421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 137d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 138421d9b32SPeter Brune 139421d9b32SPeter Brune PetscFunctionBegin; 140421d9b32SPeter Brune if (level > top_level) 141421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 142421d9b32SPeter Brune /* get to the correct level */ 143d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 144421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 145421d9b32SPeter Brune } 146421d9b32SPeter Brune if (fas->level != level) 147421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 148421d9b32SPeter Brune fas->restrct = mat; 149421d9b32SPeter Brune PetscFunctionReturn(0); 150421d9b32SPeter Brune } 151421d9b32SPeter Brune 152421d9b32SPeter Brune #undef __FUNCT__ 153421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 154421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 155421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 156d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 157421d9b32SPeter Brune 158421d9b32SPeter Brune PetscFunctionBegin; 159421d9b32SPeter Brune if (level > top_level) 160421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 161421d9b32SPeter Brune /* get to the correct level */ 162d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 163421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 164421d9b32SPeter Brune } 165421d9b32SPeter Brune if (fas->level != level) 166421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 167421d9b32SPeter Brune fas->rscale = rscale; 168421d9b32SPeter Brune PetscFunctionReturn(0); 169421d9b32SPeter Brune } 170421d9b32SPeter Brune 171421d9b32SPeter Brune #undef __FUNCT__ 172421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 173421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 174421d9b32SPeter Brune { 17577df8cc4SPeter Brune PetscErrorCode ierr = 0; 176421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 177421d9b32SPeter Brune 178421d9b32SPeter Brune PetscFunctionBegin; 179421d9b32SPeter Brune /* destroy local data created in SNESSetup_FAS */ 180*ee78dd50SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 181*ee78dd50SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 182421d9b32SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 183421d9b32SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 184421d9b32SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 185421d9b32SPeter Brune 186421d9b32SPeter Brune /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */ 187421d9b32SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 188fa9694d7SPeter Brune if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 189421d9b32SPeter Brune PetscFunctionReturn(0); 190421d9b32SPeter Brune } 191421d9b32SPeter Brune 192421d9b32SPeter Brune #undef __FUNCT__ 193421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 194421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 195421d9b32SPeter Brune { 196421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 197421d9b32SPeter Brune PetscErrorCode ierr; 198421d9b32SPeter Brune 199421d9b32SPeter Brune PetscFunctionBegin; 200421d9b32SPeter Brune /* recursively resets and then destroys */ 201421d9b32SPeter Brune ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 202421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 203421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 204421d9b32SPeter Brune PetscFunctionReturn(0); 205421d9b32SPeter Brune } 206421d9b32SPeter Brune 207421d9b32SPeter Brune #undef __FUNCT__ 208421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 209421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 210421d9b32SPeter Brune { 211421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 212421d9b32SPeter Brune PetscErrorCode ierr; 213421d9b32SPeter Brune 214421d9b32SPeter Brune PetscFunctionBegin; 215fa9694d7SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* the four work vectors are used to transfer stuff BACK */ 216421d9b32SPeter Brune /* gets the solver ready for solution */ 217421d9b32SPeter Brune if (snes->dm) { 218421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 219421d9b32SPeter Brune if (fas->next) { 220421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 221*ee78dd50SPeter Brune if (!fas->next->dm) { 222*ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 223*ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 224*ee78dd50SPeter Brune } 225421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 226*ee78dd50SPeter Brune if (!fas->interpolate) { 227421d9b32SPeter Brune ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 228421d9b32SPeter Brune fas->restrct = fas->interpolate; 229421d9b32SPeter Brune } 230421d9b32SPeter Brune } 231*ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 232*ee78dd50SPeter Brune if (fas->upsmooth) SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr); 233*ee78dd50SPeter Brune if (fas->downsmooth)SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr); 234421d9b32SPeter Brune } 235fa9694d7SPeter Brune if (fas->next) { 236fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 237fa9694d7SPeter Brune ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr); 238fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 239fa9694d7SPeter Brune } 240fa9694d7SPeter Brune /* got to set them all up at once */ 241421d9b32SPeter Brune PetscFunctionReturn(0); 242421d9b32SPeter Brune } 243421d9b32SPeter Brune 244421d9b32SPeter Brune #undef __FUNCT__ 245421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 246421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 247421d9b32SPeter Brune { 248*ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 249*ee78dd50SPeter Brune PetscInt levels = 1; 250*ee78dd50SPeter Brune PetscBool flg, monflg; 251421d9b32SPeter Brune PetscErrorCode ierr; 252*ee78dd50SPeter Brune const char * def_smooth = SNESNRICHARDSON; 253*ee78dd50SPeter Brune char pre_type[256]; 254*ee78dd50SPeter Brune char post_type[256]; 255*ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 256*ee78dd50SPeter Brune PetscViewer monviewer; 257421d9b32SPeter Brune 258421d9b32SPeter Brune PetscFunctionBegin; 259c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 260*ee78dd50SPeter Brune 261*ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 262*ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 263*ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 264*ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 265*ee78dd50SPeter Brune } 266*ee78dd50SPeter Brune 267*ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 268*ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 269*ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 270*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 271*ee78dd50SPeter Brune if (flg) { 272*ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 273*ee78dd50SPeter Brune } else { 274*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 275*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 276*ee78dd50SPeter Brune } 277*ee78dd50SPeter Brune 278*ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 279*ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 280*ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 281*ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 282*ee78dd50SPeter Brune 283*ee78dd50SPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 284*ee78dd50SPeter Brune 285*ee78dd50SPeter Brune /* other options for the coarsest level */ 286*ee78dd50SPeter Brune if (fas->level == 0) { 287*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr); 288*ee78dd50SPeter Brune if (flg) { 289*ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 290*ee78dd50SPeter Brune } else { 291*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr); 292*ee78dd50SPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr); 293*ee78dd50SPeter Brune } 294*ee78dd50SPeter Brune } 295*ee78dd50SPeter Brune 296421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 297*ee78dd50SPeter Brune /* setup from the determined types if the smoothers don't exist */ 298*ee78dd50SPeter Brune if (!fas->upsmooth) { 299*ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 300*ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 301*ee78dd50SPeter Brune } 302*ee78dd50SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 303*ee78dd50SPeter Brune if (!fas->downsmooth) { 304*ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 305*ee78dd50SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 306*ee78dd50SPeter Brune } 307*ee78dd50SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 308*ee78dd50SPeter Brune 309*ee78dd50SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 310*ee78dd50SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 311*ee78dd50SPeter Brune if (monflg) { 312*ee78dd50SPeter Brune ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr); 313*ee78dd50SPeter Brune ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 314*ee78dd50SPeter Brune ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 315*ee78dd50SPeter Brune } 316*ee78dd50SPeter Brune 317*ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 318*ee78dd50SPeter Brune if (fas->next)ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr); 319421d9b32SPeter Brune PetscFunctionReturn(0); 320421d9b32SPeter Brune } 321421d9b32SPeter Brune 322421d9b32SPeter Brune #undef __FUNCT__ 323421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 324421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 325421d9b32SPeter Brune { 326421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 327421d9b32SPeter Brune PetscBool iascii; 328421d9b32SPeter Brune PetscErrorCode ierr; 329421d9b32SPeter Brune PetscInt levels = fas->levels; 330421d9b32SPeter Brune PetscInt i; 331421d9b32SPeter Brune 332421d9b32SPeter Brune PetscFunctionBegin; 333421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 334421d9b32SPeter Brune if (iascii) { 335421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 336421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 337421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 338421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 339*ee78dd50SPeter Brune if (fas->upsmooth) { 340*ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 341421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 342*ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 343421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 344421d9b32SPeter Brune } else { 345*ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 346421d9b32SPeter Brune } 347*ee78dd50SPeter Brune if (fas->downsmooth) { 348*ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 349421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 350*ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 351421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 352421d9b32SPeter Brune } else { 353*ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 354421d9b32SPeter Brune } 355421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 356421d9b32SPeter Brune } 357421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 358421d9b32SPeter Brune } else { 359421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 360421d9b32SPeter Brune } 361421d9b32SPeter Brune PetscFunctionReturn(0); 362421d9b32SPeter Brune } 363421d9b32SPeter Brune 364421d9b32SPeter Brune /* 365fa9694d7SPeter Brune 366fa9694d7SPeter Brune Defines the FAS cycle as: 367fa9694d7SPeter Brune 368fa9694d7SPeter Brune fine problem: F(x) = 0 369fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 370fa9694d7SPeter Brune 371fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 372fa9694d7SPeter Brune 373fa9694d7SPeter Brune correction: 374fa9694d7SPeter Brune 375fa9694d7SPeter Brune x = x + I(x^c - Rx) 376fa9694d7SPeter Brune 377421d9b32SPeter Brune */ 378fa9694d7SPeter Brune 379421d9b32SPeter Brune #undef __FUNCT__ 380421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 381c90fad12SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) { 382fa9694d7SPeter Brune 383fa9694d7SPeter Brune PetscErrorCode ierr; 384fa9694d7SPeter Brune Vec X_c, Xo_c, F_c, B_c; 385fa9694d7SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 386421d9b32SPeter Brune 387421d9b32SPeter Brune PetscFunctionBegin; 388fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 389*ee78dd50SPeter Brune if (fas->upsmooth) { 390*ee78dd50SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 391c90fad12SPeter Brune } else if (snes->pc) { 392c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 393fe6f9142SPeter Brune } 394fa9694d7SPeter Brune if (fas->next) { 395*ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 396c90fad12SPeter Brune X_c = fas->next->vec_sol; 397fa9694d7SPeter Brune Xo_c = fas->next->work[1]; 398c90fad12SPeter Brune F_c = fas->next->vec_func; 399c90fad12SPeter Brune B_c = fas->next->work[2]; 400fa9694d7SPeter Brune /* project to coarse */ 401c90fad12SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 402fa9694d7SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 403c90fad12SPeter Brune ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr); 404fa9694d7SPeter Brune ierr = VecPointwiseMult(F_c, fas->rscale, F_c);CHKERRQ(ierr); 405c90fad12SPeter Brune if (B) { 406c90fad12SPeter Brune ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr); 407c90fad12SPeter Brune ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr); 408c90fad12SPeter Brune } else { 409c90fad12SPeter Brune ierr = VecSet(B_c, 0.0);CHKERRQ(ierr); 410c90fad12SPeter Brune } 411c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 412*ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 413*ee78dd50SPeter Brune ierr = VecAXPY(B_c, -1.0, F_c);CHKERRQ(ierr); /* B_c = RB + F(X_c) - RF */ 414c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 415c90fad12SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); 416c90fad12SPeter Brune 417*ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 418*ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 419c90fad12SPeter Brune 420c90fad12SPeter Brune /* recurse to the next level */ 421c90fad12SPeter Brune ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr); 422*ee78dd50SPeter Brune 423fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 424fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 425*ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 426*ee78dd50SPeter Brune 427*ee78dd50SPeter Brune #if 0 428*ee78dd50SPeter Brune /* CHECK CHECK INTERPOLATION CHECK */ 429*ee78dd50SPeter Brune PetscScalar qnorm; 430*ee78dd50SPeter Brune ierr = VecNorm(Xo_c, &qnorm);CHKERRQ(ierr); 431*ee78dd50SPeter Brune ierr = MatInterpolate(fas->interpolate, Xo_c, F);CHKERRQ(ierr); 432*ee78dd50SPeter Brune ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr); 433*ee78dd50SPeter Brune ierr = VecPointwiseMult(F_c, fas->rscale, F_c);CHKERRQ(ierr); 434*ee78dd50SPeter Brune ierr = VecAXPY(F_c, -1.0, Xo_c);CHKERRQ(ierr); 435*ee78dd50SPeter Brune PetscScalar tnorm; 436*ee78dd50SPeter Brune ierr = VecNorm(F_c, &tnorm);CHKERRQ(ierr); 437*ee78dd50SPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD, "NORM OF PROJECTION ERROR %g, (%g)\n", tnorm, qnorm);CHKERRQ(ierr); 438*ee78dd50SPeter Brune #endif 439fa9694d7SPeter Brune } 440*ee78dd50SPeter Brune /* down-smooth -- just update using the down-smoother */ 441c90fad12SPeter Brune if (fas->level != 0) { 442*ee78dd50SPeter Brune if (fas->downsmooth) { 443*ee78dd50SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 444c90fad12SPeter Brune } else if (snes->pc) { 445c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 446c90fad12SPeter Brune } 447fe6f9142SPeter Brune } 448fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 449fa9694d7SPeter Brune 450fa9694d7SPeter Brune PetscFunctionReturn(0); 451421d9b32SPeter Brune } 452421d9b32SPeter Brune 453421d9b32SPeter Brune #undef __FUNCT__ 454c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private" 455c90fad12SPeter Brune 456c90fad12SPeter Brune 457c90fad12SPeter Brune #undef __FUNCT__ 458421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 459421d9b32SPeter Brune 460421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 461421d9b32SPeter Brune { 462fa9694d7SPeter Brune PetscErrorCode ierr; 463fe6f9142SPeter Brune PetscInt i, maxits; 464fe6f9142SPeter Brune Vec X, F, B; 465fe6f9142SPeter Brune PetscReal fnorm; 466421d9b32SPeter Brune PetscFunctionBegin; 467fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 468fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 469fa9694d7SPeter Brune X = snes->vec_sol; 470fa9694d7SPeter Brune F = snes->vec_func; 471fe6f9142SPeter Brune B = snes->vec_rhs; 472fe6f9142SPeter Brune /* initial iteration */ 473fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 474fe6f9142SPeter Brune snes->iter = 0; 475fe6f9142SPeter Brune snes->norm = 0.; 476fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 477fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 478fe6f9142SPeter Brune if (snes->domainerror) { 479fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 480fe6f9142SPeter Brune PetscFunctionReturn(0); 481fe6f9142SPeter Brune } 482fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 483fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 484fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 485fe6f9142SPeter Brune snes->norm = fnorm; 486fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 487fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 488fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 489fe6f9142SPeter Brune 490fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 491fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 492fe6f9142SPeter Brune /* test convergence */ 493fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 494fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 495fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 496fe6f9142SPeter Brune /* Call general purpose update function */ 497fe6f9142SPeter Brune if (snes->ops->update) { 498fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 499fe6f9142SPeter Brune } 500c90fad12SPeter Brune ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr); 501c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 502c90fad12SPeter Brune /* Monitor convergence */ 503c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 504c90fad12SPeter Brune snes->iter = i+1; 505c90fad12SPeter Brune snes->norm = fnorm; 506c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 507c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 508c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 509c90fad12SPeter Brune /* Test for convergence */ 510c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 511c90fad12SPeter Brune if (snes->reason) break; 512fe6f9142SPeter Brune } 513fe6f9142SPeter Brune if (i == maxits) { 514fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 515fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 516fe6f9142SPeter Brune } 517421d9b32SPeter Brune PetscFunctionReturn(0); 518421d9b32SPeter Brune } 519