1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 507144faaSPeter Brune 6421d9b32SPeter Brune /*MC 7ef536925SPeter Brune 8ef536925SPeter Brune SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 9421d9b32SPeter Brune 10421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 11421d9b32SPeter Brune 12421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 13421d9b32SPeter Brune M*/ 14421d9b32SPeter Brune 15421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 16421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 18421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 19421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 20421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 216273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 22421d9b32SPeter Brune 23421d9b32SPeter Brune EXTERN_C_BEGIN 24421d9b32SPeter Brune 25421d9b32SPeter Brune #undef __FUNCT__ 26421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 27421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 28421d9b32SPeter Brune { 29421d9b32SPeter Brune SNES_FAS * fas; 30421d9b32SPeter Brune PetscErrorCode ierr; 31421d9b32SPeter Brune 32421d9b32SPeter Brune PetscFunctionBegin; 33421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 34421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 35421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 36421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 37421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 38421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 39421d9b32SPeter Brune 40ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 41ed020824SBarry Smith snes->usespc = PETSC_FALSE; 42ed020824SBarry Smith 430e444f03SPeter Brune snes->max_funcs = 30000; 440e444f03SPeter Brune snes->max_its = 10000; 450e444f03SPeter Brune 46421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 47421d9b32SPeter Brune snes->data = (void*) fas; 48421d9b32SPeter Brune fas->level = 0; 49293a7e31SPeter Brune fas->levels = 1; 50ee78dd50SPeter Brune fas->n_cycles = 1; 51ee78dd50SPeter Brune fas->max_up_it = 1; 52ee78dd50SPeter Brune fas->max_down_it = 1; 53ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 54ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 55421d9b32SPeter Brune fas->next = PETSC_NULL; 566273346dSPeter Brune fas->previous = PETSC_NULL; 57421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 58421d9b32SPeter Brune fas->restrct = PETSC_NULL; 59efe1f98aSPeter Brune fas->inject = PETSC_NULL; 60cc05f883SPeter Brune fas->monitor = PETSC_NULL; 61cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 62ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 6322c1e704SPeter Brune fas->linesearch_smooth = PETSC_NULL; 64ddebd997SPeter Brune 65421d9b32SPeter Brune PetscFunctionReturn(0); 66421d9b32SPeter Brune } 67421d9b32SPeter Brune EXTERN_C_END 68421d9b32SPeter Brune 69421d9b32SPeter Brune #undef __FUNCT__ 70421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 71c79ef259SPeter Brune /*@ 722e8ce248SJed Brown SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 73c79ef259SPeter Brune 74c79ef259SPeter Brune Input Parameter: 752e8ce248SJed Brown . snes - the nonlinear solver context 76c79ef259SPeter Brune 77c79ef259SPeter Brune Output parameter: 78c79ef259SPeter Brune . levels - the number of levels 79c79ef259SPeter Brune 80c79ef259SPeter Brune Level: advanced 81c79ef259SPeter Brune 82c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 83c79ef259SPeter Brune 84c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 85c79ef259SPeter Brune @*/ 86421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 87421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 88421d9b32SPeter Brune PetscFunctionBegin; 89ee1fd11aSPeter Brune *levels = fas->levels; 90421d9b32SPeter Brune PetscFunctionReturn(0); 91421d9b32SPeter Brune } 92421d9b32SPeter Brune 93421d9b32SPeter Brune #undef __FUNCT__ 94646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 95c79ef259SPeter Brune /*@ 962e8ce248SJed Brown SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 97c79ef259SPeter Brune complicated cycling. 98c79ef259SPeter Brune 99c79ef259SPeter Brune Logically Collective on SNES 100c79ef259SPeter Brune 101c79ef259SPeter Brune Input Parameters: 102c79ef259SPeter Brune + snes - the multigrid context 103c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 104c79ef259SPeter Brune 105c79ef259SPeter Brune Options Database Key: 106c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 107c79ef259SPeter Brune 108c79ef259SPeter Brune Level: advanced 109c79ef259SPeter Brune 110c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 111c79ef259SPeter Brune 112c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 113c79ef259SPeter Brune @*/ 114646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 115646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 116646217ecSPeter Brune PetscErrorCode ierr; 117646217ecSPeter Brune PetscFunctionBegin; 118646217ecSPeter Brune fas->n_cycles = cycles; 119646217ecSPeter Brune if (fas->next) { 120646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 121646217ecSPeter Brune } 122646217ecSPeter Brune PetscFunctionReturn(0); 123646217ecSPeter Brune } 124646217ecSPeter Brune 125eff52c0eSPeter Brune #undef __FUNCT__ 126c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 127c79ef259SPeter Brune /*@ 128722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 129c79ef259SPeter Brune 130c79ef259SPeter Brune Logically Collective on SNES 131c79ef259SPeter Brune 132c79ef259SPeter Brune Input Parameters: 133c79ef259SPeter Brune + snes - the multigrid context 134c79ef259SPeter Brune . level - the level to set the number of cycles on 135c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 136c79ef259SPeter Brune 137c79ef259SPeter Brune Level: advanced 138c79ef259SPeter Brune 139c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 140c79ef259SPeter Brune 141c79ef259SPeter Brune .seealso: SNESFASSetCycles() 142c79ef259SPeter Brune @*/ 143c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 144c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 145c79ef259SPeter Brune PetscInt top_level = fas->level,i; 146c79ef259SPeter Brune 147c79ef259SPeter Brune PetscFunctionBegin; 1482e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D cycle count from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 149c79ef259SPeter Brune /* get to the correct level */ 150c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 151c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 152c79ef259SPeter Brune } 1532e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 154c79ef259SPeter Brune fas->n_cycles = cycles; 155c79ef259SPeter Brune PetscFunctionReturn(0); 156c79ef259SPeter Brune } 157c79ef259SPeter Brune 158c79ef259SPeter Brune #undef __FUNCT__ 159eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 160aeed3662SMatthew G Knepley /*@C 161c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 162c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 163c79ef259SPeter Brune and nonlinear preconditioners. 164c79ef259SPeter Brune 165c79ef259SPeter Brune Logically Collective on SNES 166c79ef259SPeter Brune 167c79ef259SPeter Brune Input Parameters: 168c79ef259SPeter Brune + snes - the multigrid context 169c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 170c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 171c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 172c79ef259SPeter Brune 173c79ef259SPeter Brune Level: advanced 174c79ef259SPeter Brune 175c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 176c79ef259SPeter Brune 177c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 178c79ef259SPeter Brune @*/ 179eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 180eff52c0eSPeter Brune PetscErrorCode ierr = 0; 181eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 182eff52c0eSPeter Brune PetscFunctionBegin; 183eff52c0eSPeter Brune 184eff52c0eSPeter Brune if (gsfunc) { 185eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 186eff52c0eSPeter Brune /* push the provided GS up the tree */ 187eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 1886cab3a1bSJed Brown } else { 1896cab3a1bSJed Brown ierr = SNESGetGS(snes,&gsfunc,&ctx);CHKERRQ(ierr); 1906cab3a1bSJed Brown if (gsfunc) { 191eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 192eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 193eff52c0eSPeter Brune } else if (use_gs) { 1942e8ce248SJed Brown SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level); 195eff52c0eSPeter Brune } 1966cab3a1bSJed Brown } 197eff52c0eSPeter Brune PetscFunctionReturn(0); 198eff52c0eSPeter Brune } 199eff52c0eSPeter Brune 200eff52c0eSPeter Brune #undef __FUNCT__ 201eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 202aeed3662SMatthew G Knepley /*@C 203c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 204c79ef259SPeter Brune 205c79ef259SPeter Brune Logically Collective on SNES 206c79ef259SPeter Brune 207c79ef259SPeter Brune Input Parameters: 208c79ef259SPeter Brune + snes - the multigrid context 209c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 210c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 211c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 212c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 213c79ef259SPeter Brune 214c79ef259SPeter Brune Level: advanced 215c79ef259SPeter Brune 216c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 217c79ef259SPeter Brune 218c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 219c79ef259SPeter Brune @*/ 220eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 221eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 222eff52c0eSPeter Brune PetscErrorCode ierr; 223eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 224eff52c0eSPeter Brune SNES cur_snes = snes; 225eff52c0eSPeter Brune PetscFunctionBegin; 2262e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D GS from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 227eff52c0eSPeter Brune /* get to the correct level */ 228eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 229eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 230eff52c0eSPeter Brune cur_snes = fas->next; 231eff52c0eSPeter Brune } 2322e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 233eff52c0eSPeter Brune if (gsfunc) { 2346273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 235eff52c0eSPeter Brune } 236eff52c0eSPeter Brune PetscFunctionReturn(0); 237eff52c0eSPeter Brune } 238eff52c0eSPeter Brune 239646217ecSPeter Brune #undef __FUNCT__ 240421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 241c79ef259SPeter Brune /*@ 242c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 243c79ef259SPeter Brune level of the FAS hierarchy. 244c79ef259SPeter Brune 245c79ef259SPeter Brune Input Parameters: 246c79ef259SPeter Brune + snes - the multigrid context 247c79ef259SPeter Brune level - the level to get 248c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 249c79ef259SPeter Brune 250c79ef259SPeter Brune Level: advanced 251c79ef259SPeter Brune 252c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 253c79ef259SPeter Brune 254c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 255c79ef259SPeter Brune @*/ 256421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) { 257421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 258421d9b32SPeter Brune PetscInt i; 2592e8ce248SJed Brown 260421d9b32SPeter Brune PetscFunctionBegin; 2612e8ce248SJed Brown if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 2622e8ce248SJed Brown if (fas->level < level) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS on level %D, must call from SNES at least as fine as level",level,fas->level); 263421d9b32SPeter Brune *lsnes = snes; 264421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 265421d9b32SPeter Brune *lsnes = fas->next; 266421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 267421d9b32SPeter Brune } 2682e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 269421d9b32SPeter Brune PetscFunctionReturn(0); 270421d9b32SPeter Brune } 271421d9b32SPeter Brune 272421d9b32SPeter Brune #undef __FUNCT__ 27307144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 27407144faaSPeter Brune /*@ 27507144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 27607144faaSPeter Brune e 27707144faaSPeter Brune 27807144faaSPeter Brune 27907144faaSPeter Brune @*/ 28007144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 28107144faaSPeter Brune { 28207144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 28307144faaSPeter Brune 28407144faaSPeter Brune PetscFunctionBegin; 28507144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 28607144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 28707144faaSPeter Brune fas->fastype = fastype; 28807144faaSPeter Brune PetscFunctionReturn(0); 28907144faaSPeter Brune } 29007144faaSPeter Brune 29107144faaSPeter Brune #undef __FUNCT__ 292421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 293aeed3662SMatthew G Knepley /*@C 294c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 295c79ef259SPeter Brune Must be called before any other FAS routine. 296c79ef259SPeter Brune 297c79ef259SPeter Brune Input Parameters: 298c79ef259SPeter Brune + snes - the snes context 299c79ef259SPeter Brune . levels - the number of levels 300c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 301c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 302c79ef259SPeter Brune Fortran. 303c79ef259SPeter Brune 304c79ef259SPeter Brune Level: intermediate 305c79ef259SPeter Brune 306c79ef259SPeter Brune Notes: 307c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 308c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 309c79ef259SPeter Brune 310c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 311c79ef259SPeter Brune 312c79ef259SPeter Brune .seealso: SNESFASGetLevels() 313c79ef259SPeter Brune @*/ 314421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 315421d9b32SPeter Brune PetscErrorCode ierr; 316421d9b32SPeter Brune PetscInt i; 317421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3186273346dSPeter Brune SNES prevsnes; 319421d9b32SPeter Brune MPI_Comm comm; 320421d9b32SPeter Brune PetscFunctionBegin; 321ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 322ee1fd11aSPeter Brune if (levels == fas->levels) { 323ee1fd11aSPeter Brune if (!comms) 324ee1fd11aSPeter Brune PetscFunctionReturn(0); 325ee1fd11aSPeter Brune } 326421d9b32SPeter Brune /* user has changed the number of levels; reset */ 327421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 328421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 329421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 330ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3316273346dSPeter Brune fas->previous = PETSC_NULL; 3326273346dSPeter Brune prevsnes = snes; 333421d9b32SPeter Brune /* setup the finest level */ 334421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 335421d9b32SPeter Brune if (comms) comm = comms[i]; 336421d9b32SPeter Brune fas->level = i; 337421d9b32SPeter Brune fas->levels = levels; 338ee1fd11aSPeter Brune fas->next = PETSC_NULL; 339421d9b32SPeter Brune if (i > 0) { 340421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3414a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 342421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 343794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3446273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3456273346dSPeter Brune prevsnes = fas->next; 3466273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 347421d9b32SPeter Brune } 348421d9b32SPeter Brune } 349421d9b32SPeter Brune PetscFunctionReturn(0); 350421d9b32SPeter Brune } 351421d9b32SPeter Brune 352421d9b32SPeter Brune #undef __FUNCT__ 353c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 354c79ef259SPeter Brune /*@ 355c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 356c79ef259SPeter Brune use on all levels. 357c79ef259SPeter Brune 358fde0ff24SPeter Brune Logically Collective on SNES 359c79ef259SPeter Brune 360c79ef259SPeter Brune Input Parameters: 361c79ef259SPeter Brune + snes - the multigrid context 362c79ef259SPeter Brune - n - the number of smoothing steps 363c79ef259SPeter Brune 364c79ef259SPeter Brune Options Database Key: 365d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 366c79ef259SPeter Brune 367c79ef259SPeter Brune Level: advanced 368c79ef259SPeter Brune 369c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 370c79ef259SPeter Brune 371c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 372c79ef259SPeter Brune @*/ 373c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 374c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 375c79ef259SPeter Brune PetscErrorCode ierr = 0; 376c79ef259SPeter Brune PetscFunctionBegin; 377c79ef259SPeter Brune fas->max_up_it = n; 378c79ef259SPeter Brune if (fas->next) { 379c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 380c79ef259SPeter Brune } 381c79ef259SPeter Brune PetscFunctionReturn(0); 382c79ef259SPeter Brune } 383c79ef259SPeter Brune 384c79ef259SPeter Brune #undef __FUNCT__ 385c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 386c79ef259SPeter Brune /*@ 387c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 388c79ef259SPeter Brune use on all levels. 389c79ef259SPeter Brune 390fde0ff24SPeter Brune Logically Collective on SNES 391c79ef259SPeter Brune 392c79ef259SPeter Brune Input Parameters: 393c79ef259SPeter Brune + snes - the multigrid context 394c79ef259SPeter Brune - n - the number of smoothing steps 395c79ef259SPeter Brune 396c79ef259SPeter Brune Options Database Key: 397d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 398c79ef259SPeter Brune 399c79ef259SPeter Brune Level: advanced 400c79ef259SPeter Brune 401c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 402c79ef259SPeter Brune 403c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 404c79ef259SPeter Brune @*/ 405c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 406c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 407c79ef259SPeter Brune PetscErrorCode ierr = 0; 408c79ef259SPeter Brune PetscFunctionBegin; 409c79ef259SPeter Brune fas->max_down_it = n; 410c79ef259SPeter Brune if (fas->next) { 411c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 412c79ef259SPeter Brune } 413c79ef259SPeter Brune PetscFunctionReturn(0); 414c79ef259SPeter Brune } 415c79ef259SPeter Brune 416c79ef259SPeter Brune #undef __FUNCT__ 417421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 418c79ef259SPeter Brune /*@ 419c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 420c79ef259SPeter Brune interpolation from l-1 to the lth level 421c79ef259SPeter Brune 422c79ef259SPeter Brune Input Parameters: 423c79ef259SPeter Brune + snes - the multigrid context 424c79ef259SPeter Brune . mat - the interpolation operator 425c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 426c79ef259SPeter Brune 427c79ef259SPeter Brune Level: advanced 428c79ef259SPeter Brune 429c79ef259SPeter Brune Notes: 430c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 431c79ef259SPeter Brune for the same level. 432c79ef259SPeter Brune 433c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 434c79ef259SPeter Brune out from the matrix size which one it is. 435c79ef259SPeter Brune 436c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 437c79ef259SPeter Brune 438bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 439c79ef259SPeter Brune @*/ 440421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 441421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 442d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 443bccf9bb3SJed Brown PetscErrorCode ierr; 444421d9b32SPeter Brune 445421d9b32SPeter Brune PetscFunctionBegin; 4462e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D interpolation from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 447421d9b32SPeter Brune /* get to the correct level */ 448d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 449421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 450421d9b32SPeter Brune } 4512e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 452bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 453bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 454421d9b32SPeter Brune fas->interpolate = mat; 455421d9b32SPeter Brune PetscFunctionReturn(0); 456421d9b32SPeter Brune } 457421d9b32SPeter Brune 458421d9b32SPeter Brune #undef __FUNCT__ 459421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 460c79ef259SPeter Brune /*@ 461c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 462c79ef259SPeter Brune from level l to l-1. 463c79ef259SPeter Brune 464c79ef259SPeter Brune Input Parameters: 465c79ef259SPeter Brune + snes - the multigrid context 466c79ef259SPeter Brune . mat - the restriction matrix 467c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 468c79ef259SPeter Brune 469c79ef259SPeter Brune Level: advanced 470c79ef259SPeter Brune 471c79ef259SPeter Brune Notes: 472c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 473c79ef259SPeter Brune for the same level. 474c79ef259SPeter Brune 475c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 476c79ef259SPeter Brune out from the matrix size which one it is. 477c79ef259SPeter Brune 478fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 479c79ef259SPeter Brune is used. 480c79ef259SPeter Brune 481c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 482c79ef259SPeter Brune 483c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 484c79ef259SPeter Brune @*/ 485421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 486421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 487d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 488bccf9bb3SJed Brown PetscErrorCode ierr; 489421d9b32SPeter Brune 490421d9b32SPeter Brune PetscFunctionBegin; 4912e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D restriction from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 492421d9b32SPeter Brune /* get to the correct level */ 493d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 494421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 495421d9b32SPeter Brune } 4962e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 497bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 498bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 499421d9b32SPeter Brune fas->restrct = mat; 500421d9b32SPeter Brune PetscFunctionReturn(0); 501421d9b32SPeter Brune } 502421d9b32SPeter Brune 503421d9b32SPeter Brune #undef __FUNCT__ 504421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 505c79ef259SPeter Brune /*@ 506c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 507c79ef259SPeter Brune operator from level l to l-1. 508c79ef259SPeter Brune 509c79ef259SPeter Brune Input Parameters: 510c79ef259SPeter Brune + snes - the multigrid context 511c79ef259SPeter Brune . rscale - the restriction scaling 512c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 513c79ef259SPeter Brune 514c79ef259SPeter Brune Level: advanced 515c79ef259SPeter Brune 516c79ef259SPeter Brune Notes: 517c79ef259SPeter Brune This is only used in the case that the injection is not set. 518c79ef259SPeter Brune 519c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 520c79ef259SPeter Brune 521c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 522c79ef259SPeter Brune @*/ 523421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 524421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 525d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 526bccf9bb3SJed Brown PetscErrorCode ierr; 527421d9b32SPeter Brune 528421d9b32SPeter Brune PetscFunctionBegin; 5292e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D rscale from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 530421d9b32SPeter Brune /* get to the correct level */ 531d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 532421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 533421d9b32SPeter Brune } 5342e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 535bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 536bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 537421d9b32SPeter Brune fas->rscale = rscale; 538421d9b32SPeter Brune PetscFunctionReturn(0); 539421d9b32SPeter Brune } 540421d9b32SPeter Brune 541efe1f98aSPeter Brune 542efe1f98aSPeter Brune #undef __FUNCT__ 543efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 544c79ef259SPeter Brune /*@ 545c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 546c79ef259SPeter Brune from level l to l-1. 547c79ef259SPeter Brune 548c79ef259SPeter Brune Input Parameters: 549c79ef259SPeter Brune + snes - the multigrid context 550c79ef259SPeter Brune . mat - the restriction matrix 551c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 552c79ef259SPeter Brune 553c79ef259SPeter Brune Level: advanced 554c79ef259SPeter Brune 555c79ef259SPeter Brune Notes: 556c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 557c79ef259SPeter Brune project the solution instead. 558c79ef259SPeter Brune 559c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 560c79ef259SPeter Brune 561c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 562c79ef259SPeter Brune @*/ 563efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 564efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 565efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 566bccf9bb3SJed Brown PetscErrorCode ierr; 567efe1f98aSPeter Brune 568efe1f98aSPeter Brune PetscFunctionBegin; 5692e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D injection from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 570efe1f98aSPeter Brune /* get to the correct level */ 571efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 572efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 573efe1f98aSPeter Brune } 5742e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 575bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 576bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 577efe1f98aSPeter Brune fas->inject = mat; 578efe1f98aSPeter Brune PetscFunctionReturn(0); 579efe1f98aSPeter Brune } 580efe1f98aSPeter Brune 581421d9b32SPeter Brune #undef __FUNCT__ 582421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 583421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 584421d9b32SPeter Brune { 58577df8cc4SPeter Brune PetscErrorCode ierr = 0; 586421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 587421d9b32SPeter Brune 588421d9b32SPeter Brune PetscFunctionBegin; 589bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 590bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 5913dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 592bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 593bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 594bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 5953dccd265SPeter Brune 5963dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 5973dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 598742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 59922c1e704SPeter Brune 6006188f407SPeter Brune ierr = PetscLineSearchDestroy(&fas->linesearch_smooth);CHKERRQ(ierr); 60122c1e704SPeter Brune 602421d9b32SPeter Brune PetscFunctionReturn(0); 603421d9b32SPeter Brune } 604421d9b32SPeter Brune 605421d9b32SPeter Brune #undef __FUNCT__ 606421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 607421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 608421d9b32SPeter Brune { 609421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 610742fe5e2SPeter Brune PetscErrorCode ierr = 0; 611421d9b32SPeter Brune 612421d9b32SPeter Brune PetscFunctionBegin; 613421d9b32SPeter Brune /* recursively resets and then destroys */ 61479d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 615421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 616421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 617ee1fd11aSPeter Brune 618421d9b32SPeter Brune PetscFunctionReturn(0); 619421d9b32SPeter Brune } 620421d9b32SPeter Brune 621421d9b32SPeter Brune #undef __FUNCT__ 622421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 623421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 624421d9b32SPeter Brune { 62548bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 626421d9b32SPeter Brune PetscErrorCode ierr; 62722c1e704SPeter Brune const char *optionsprefix; 628efe1f98aSPeter Brune VecScatter injscatter; 629d1adcc6fSPeter Brune PetscInt dm_levels; 6303dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 631d1adcc6fSPeter Brune 632421d9b32SPeter Brune PetscFunctionBegin; 633eff52c0eSPeter Brune 634cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 635d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 636d1adcc6fSPeter Brune dm_levels++; 637cc05f883SPeter Brune if (dm_levels > fas->levels) { 6383dccd265SPeter Brune 6392e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6403dccd265SPeter Brune vec_sol = snes->vec_sol; 6413dccd265SPeter Brune vec_func = snes->vec_func; 6423dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6433dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6443dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6453dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6463dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6473dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6483dccd265SPeter Brune 6493dccd265SPeter Brune /* reset the number of levels */ 650d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 651cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6523dccd265SPeter Brune 6533dccd265SPeter Brune snes->vec_sol = vec_sol; 6543dccd265SPeter Brune snes->vec_func = vec_func; 6553dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6563dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 657d1adcc6fSPeter Brune } 658d1adcc6fSPeter Brune } 659d1adcc6fSPeter Brune 6603dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6613dccd265SPeter Brune 66207144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 663fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 66407144faaSPeter Brune } else { 665ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 66607144faaSPeter Brune } 667cc05f883SPeter Brune 66879d9a41aSPeter Brune if (snes->dm) { 66979d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 67079d9a41aSPeter Brune if (fas->next) { 67179d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 67279d9a41aSPeter Brune if (!fas->next->dm) { 67379d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 67479d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 67579d9a41aSPeter Brune } 67679d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 67779d9a41aSPeter Brune if (!fas->interpolate) { 67879d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 679bccf9bb3SJed Brown if (!fas->restrct) { 680bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 68179d9a41aSPeter Brune fas->restrct = fas->interpolate; 68279d9a41aSPeter Brune } 683bccf9bb3SJed Brown } 68479d9a41aSPeter Brune /* set the injection from the DM */ 68579d9a41aSPeter Brune if (!fas->inject) { 68679d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 68779d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 68879d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 68979d9a41aSPeter Brune } 69079d9a41aSPeter Brune } 69179d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 69279d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 69379d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 69479d9a41aSPeter Brune } 69579d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 69679d9a41aSPeter Brune 69779d9a41aSPeter Brune if (fas->next) { 69879d9a41aSPeter Brune if (fas->galerkin) { 69979d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 70079d9a41aSPeter Brune } 70179d9a41aSPeter Brune } 70279d9a41aSPeter Brune 70379d9a41aSPeter Brune if (fas->next) { 70479d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 705938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 706938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 70779d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 70879d9a41aSPeter Brune } 70979d9a41aSPeter Brune 7106cab3a1bSJed Brown /* setup the pre and post smoothers */ 7116cab3a1bSJed Brown if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);} 7126cab3a1bSJed Brown if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);} 713cc05f883SPeter Brune 71422c1e704SPeter Brune /* if the pre and post smoothers don't exist, set up line searches in their place */ 71522c1e704SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 71622c1e704SPeter Brune if (!fas->upsmooth || !fas->downsmooth) { 7176188f407SPeter Brune ierr = PetscLineSearchCreate(((PetscObject)snes)->comm, &fas->linesearch_smooth);CHKERRQ(ierr); 7186188f407SPeter Brune ierr = PetscLineSearchSetSNES(fas->linesearch_smooth, snes);CHKERRQ(ierr); 7196188f407SPeter Brune ierr = PetscLineSearchSetType(fas->linesearch_smooth, PETSCLINESEARCHL2);CHKERRQ(ierr); 7206188f407SPeter Brune ierr = PetscLineSearchAppendOptionsPrefix(fas->linesearch_smooth, "fas_");CHKERRQ(ierr); 7216188f407SPeter Brune ierr = PetscLineSearchAppendOptionsPrefix(fas->linesearch_smooth, optionsprefix);CHKERRQ(ierr); 7226188f407SPeter Brune ierr = PetscLineSearchSetFromOptions(fas->linesearch_smooth);CHKERRQ(ierr); 72322c1e704SPeter Brune } 72422c1e704SPeter Brune 7256273346dSPeter Brune /* setup FAS work vectors */ 7266273346dSPeter Brune if (fas->galerkin) { 7276273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7286273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7296273346dSPeter Brune } 730421d9b32SPeter Brune PetscFunctionReturn(0); 731421d9b32SPeter Brune } 732421d9b32SPeter Brune 733421d9b32SPeter Brune #undef __FUNCT__ 734421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 735421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 736421d9b32SPeter Brune { 737ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 738ee78dd50SPeter Brune PetscInt levels = 1; 739fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 740421d9b32SPeter Brune PetscErrorCode ierr; 741ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 74207144faaSPeter Brune SNESFASType fastype; 743fde0ff24SPeter Brune const char *optionsprefix; 744fde0ff24SPeter Brune const char *prefix; 745*9e764e56SPeter Brune PetscLineSearch linesearch; 746421d9b32SPeter Brune 747421d9b32SPeter Brune PetscFunctionBegin; 748c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 749ee78dd50SPeter Brune 750ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 751ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 752ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 753c732cbdbSBarry Smith if (!flg && snes->dm) { 754c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 755c732cbdbSBarry Smith levels++; 756d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 757c732cbdbSBarry Smith } 758ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 759ee78dd50SPeter Brune } 76007144faaSPeter Brune fastype = fas->fastype; 76107144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 76207144faaSPeter Brune if (flg) { 76307144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 76407144faaSPeter Brune } 765ee78dd50SPeter Brune 766fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 767fde0ff24SPeter Brune 768fde0ff24SPeter Brune /* smoother setup options */ 769fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 770fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 771fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 772fde0ff24SPeter Brune if (fas->level == 0) { 773fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 774fde0ff24SPeter Brune } 775ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 776fde0ff24SPeter Brune 777d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 778ee78dd50SPeter Brune 7796273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 7806273346dSPeter Brune 781646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 782ee78dd50SPeter Brune 783ee78dd50SPeter Brune 784421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7858cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 786eff52c0eSPeter Brune 787fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 788fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 789fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 790fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 791fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 792fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 793fde0ff24SPeter Brune } 794fde0ff24SPeter Brune 795fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 7968cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 7978cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 7988cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 799fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8008cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8018cc86e31SPeter Brune } 8028cc86e31SPeter Brune 803fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 80467339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 805ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 80667339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 807fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 808293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 809ee78dd50SPeter Brune } 810fde0ff24SPeter Brune 811fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 812fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 813fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 814fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 815fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 816fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 817fde0ff24SPeter Brune } 818fde0ff24SPeter Brune 819fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 820fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 821fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 822fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 823fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 824fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 825fde0ff24SPeter Brune } 826fde0ff24SPeter Brune 8271a266240SBarry Smith if (fas->upsmooth) { 828fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8291a266240SBarry Smith } 8301a266240SBarry Smith 8311a266240SBarry Smith if (fas->downsmooth) { 832fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8331a266240SBarry Smith } 834ee78dd50SPeter Brune 835742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 836fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 837742fe5e2SPeter Brune } 838742fe5e2SPeter Brune 839ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 840ce11c761SPeter Brune if (!fas->downsmooth) { 84193dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 8425d115551SPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 843ce11c761SPeter Brune if (fas->level == 0) { 84479d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 845ce11c761SPeter Brune } 846ce11c761SPeter Brune } 847ce11c761SPeter Brune 848ce11c761SPeter Brune if (!fas->upsmooth) { 84993dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 8505d115551SPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 851ce11c761SPeter Brune } 852ce11c761SPeter Brune 853ee78dd50SPeter Brune if (monflg) { 854646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 855794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8562f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 857742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 858293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 859293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 860d28543b3SPeter Brune } else { 861d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 862d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 863d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 864d28543b3SPeter Brune } 865ee78dd50SPeter Brune } 866ee78dd50SPeter Brune 867*9e764e56SPeter Brune 868*9e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 869*9e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 870*9e764e56SPeter Brune if (!snes->linesearch) { 871*9e764e56SPeter Brune ierr = SNESGetPetscLineSearch(snes, &linesearch);CHKERRQ(ierr); 872*9e764e56SPeter Brune ierr = PetscLineSearchSetType(linesearch, PETSCLINESEARCHL2);CHKERRQ(ierr); 873*9e764e56SPeter Brune } 874*9e764e56SPeter Brune } 875*9e764e56SPeter Brune 876ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 877d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 878421d9b32SPeter Brune PetscFunctionReturn(0); 879421d9b32SPeter Brune } 880421d9b32SPeter Brune 881421d9b32SPeter Brune #undef __FUNCT__ 882421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 883421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 884421d9b32SPeter Brune { 885421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 886421d9b32SPeter Brune PetscBool iascii; 887421d9b32SPeter Brune PetscErrorCode ierr; 888421d9b32SPeter Brune 889421d9b32SPeter Brune PetscFunctionBegin; 890421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 891421d9b32SPeter Brune if (iascii) { 8922e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 893421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 894bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 895ee78dd50SPeter Brune if (fas->upsmooth) { 89639a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 897421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 898ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 899421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 900421d9b32SPeter Brune } else { 901ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 902421d9b32SPeter Brune } 903ee78dd50SPeter Brune if (fas->downsmooth) { 90439a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 905421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 906ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 907421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 908421d9b32SPeter Brune } else { 909ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 910421d9b32SPeter Brune } 911421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 912421d9b32SPeter Brune } else { 913421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 914421d9b32SPeter Brune } 915421d9b32SPeter Brune PetscFunctionReturn(0); 916421d9b32SPeter Brune } 917421d9b32SPeter Brune 918421d9b32SPeter Brune #undef __FUNCT__ 91939bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92039bd7f45SPeter Brune /* 92139bd7f45SPeter Brune Defines the action of the downsmoother 92239bd7f45SPeter Brune */ 92339bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 92439bd7f45SPeter Brune PetscErrorCode ierr = 0; 92522c1e704SPeter Brune PetscReal fnorm; 926742fe5e2SPeter Brune SNESConvergedReason reason; 92739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 92822c1e704SPeter Brune Vec Y, FPC; 929fde0ff24SPeter Brune PetscBool lssuccess; 930fde0ff24SPeter Brune PetscInt k; 931421d9b32SPeter Brune PetscFunctionBegin; 932fde0ff24SPeter Brune Y = snes->work[3]; 933d1adcc6fSPeter Brune if (fas->downsmooth) { 934d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 935742fe5e2SPeter Brune /* check convergence reason for the smoother */ 936742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 937742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 938742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 939742fe5e2SPeter Brune PetscFunctionReturn(0); 940742fe5e2SPeter Brune } 9414b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9424b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 943fde0ff24SPeter Brune } else { 944fde0ff24SPeter Brune /* basic richardson smoother */ 945fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 946794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 947794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 948fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 9496188f407SPeter Brune ierr = PetscLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 9506188f407SPeter Brune ierr = PetscLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 951fde0ff24SPeter Brune if (!lssuccess) { 952fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 953fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9542f7ea302SPeter Brune PetscFunctionReturn(0); 9552f7ea302SPeter Brune } 956fe6f9142SPeter Brune } 957fde0ff24SPeter Brune } 958fde0ff24SPeter Brune } 95939bd7f45SPeter Brune PetscFunctionReturn(0); 96039bd7f45SPeter Brune } 96139bd7f45SPeter Brune 96239bd7f45SPeter Brune 96339bd7f45SPeter Brune #undef __FUNCT__ 96439bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 96539bd7f45SPeter Brune /* 96607144faaSPeter Brune Defines the action of the upsmoother 96739bd7f45SPeter Brune */ 96839bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 96939bd7f45SPeter Brune PetscErrorCode ierr = 0; 97022c1e704SPeter Brune PetscReal fnorm; 97139bd7f45SPeter Brune SNESConvergedReason reason; 97239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 97322c1e704SPeter Brune Vec Y, FPC; 974fde0ff24SPeter Brune PetscBool lssuccess; 975fde0ff24SPeter Brune PetscInt k; 97639bd7f45SPeter Brune PetscFunctionBegin; 977fde0ff24SPeter Brune Y = snes->work[3]; 97839bd7f45SPeter Brune if (fas->upsmooth) { 979fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 98039bd7f45SPeter Brune /* check convergence reason for the smoother */ 981fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 98239bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 98339bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 98439bd7f45SPeter Brune PetscFunctionReturn(0); 98539bd7f45SPeter Brune } 9864b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9874b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 988fde0ff24SPeter Brune } else { 989fde0ff24SPeter Brune /* basic richardson smoother */ 990fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 99139bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 99239bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 993fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 9946188f407SPeter Brune ierr = PetscLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 9956188f407SPeter Brune ierr = PetscLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 996fde0ff24SPeter Brune if (!lssuccess) { 997fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 998fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 99939bd7f45SPeter Brune PetscFunctionReturn(0); 100039bd7f45SPeter Brune } 100139bd7f45SPeter Brune } 1002fde0ff24SPeter Brune } 1003fde0ff24SPeter Brune } 100439bd7f45SPeter Brune PetscFunctionReturn(0); 100539bd7f45SPeter Brune } 100639bd7f45SPeter Brune 100739bd7f45SPeter Brune #undef __FUNCT__ 1008938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1009938e4a01SJed Brown /*@ 1010938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1011938e4a01SJed Brown 1012938e4a01SJed Brown Collective 1013938e4a01SJed Brown 1014938e4a01SJed Brown Input Arguments: 1015938e4a01SJed Brown . snes - SNESFAS 1016938e4a01SJed Brown 1017938e4a01SJed Brown Output Arguments: 1018938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1019938e4a01SJed Brown 1020938e4a01SJed Brown Level: developer 1021938e4a01SJed Brown 1022938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1023938e4a01SJed Brown @*/ 1024938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1025938e4a01SJed Brown { 1026938e4a01SJed Brown PetscErrorCode ierr; 1027938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1028938e4a01SJed Brown 1029938e4a01SJed Brown PetscFunctionBegin; 1030938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1031938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1032938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1033938e4a01SJed Brown PetscFunctionReturn(0); 1034938e4a01SJed Brown } 1035938e4a01SJed Brown 1036e9923e8dSJed Brown #undef __FUNCT__ 1037e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1038e9923e8dSJed Brown /*@ 1039e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1040e9923e8dSJed Brown 1041e9923e8dSJed Brown Collective 1042e9923e8dSJed Brown 1043e9923e8dSJed Brown Input Arguments: 1044e9923e8dSJed Brown + fine - SNES from which to restrict 1045e9923e8dSJed Brown - Xfine - vector to restrict 1046e9923e8dSJed Brown 1047e9923e8dSJed Brown Output Arguments: 1048e9923e8dSJed Brown . Xcoarse - result of restriction 1049e9923e8dSJed Brown 1050e9923e8dSJed Brown Level: developer 1051e9923e8dSJed Brown 1052e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1053e9923e8dSJed Brown @*/ 1054e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1055e9923e8dSJed Brown { 1056e9923e8dSJed Brown PetscErrorCode ierr; 1057e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1058e9923e8dSJed Brown 1059e9923e8dSJed Brown PetscFunctionBegin; 1060e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1061e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1062e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1063e9923e8dSJed Brown if (fas->inject) { 1064e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1065e9923e8dSJed Brown } else { 1066e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1067e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1068e9923e8dSJed Brown } 1069e9923e8dSJed Brown PetscFunctionReturn(0); 1070e9923e8dSJed Brown } 1071e9923e8dSJed Brown 1072e9923e8dSJed Brown #undef __FUNCT__ 107339bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 107439bd7f45SPeter Brune /* 107539bd7f45SPeter Brune 107639bd7f45SPeter Brune Performs the FAS coarse correction as: 107739bd7f45SPeter Brune 107839bd7f45SPeter Brune fine problem: F(x) = 0 107939bd7f45SPeter Brune coarse problem: F^c(x) = b^c 108039bd7f45SPeter Brune 108139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 108239bd7f45SPeter Brune 108339bd7f45SPeter Brune with correction: 108439bd7f45SPeter Brune 108539bd7f45SPeter Brune 108639bd7f45SPeter Brune 108739bd7f45SPeter Brune */ 108839a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 108939bd7f45SPeter Brune PetscErrorCode ierr; 109039bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 109139bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 109239bd7f45SPeter Brune SNESConvergedReason reason; 109339bd7f45SPeter Brune PetscFunctionBegin; 1094fa9694d7SPeter Brune if (fas->next) { 1095c90fad12SPeter Brune X_c = fas->next->vec_sol; 1096293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1097c90fad12SPeter Brune F_c = fas->next->vec_func; 1098742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1099efe1f98aSPeter Brune 1100938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1101293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1102293a7e31SPeter Brune 1103293a7e31SPeter Brune /* restrict the defect */ 1104293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1105293a7e31SPeter Brune 1106c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1107ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1108c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1109742fe5e2SPeter Brune 1110293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1111c90fad12SPeter Brune 1112ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1113ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1114c90fad12SPeter Brune 1115c90fad12SPeter Brune /* recurse to the next level */ 1116f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1117742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1118742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1119742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1120742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1121742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1122742fe5e2SPeter Brune PetscFunctionReturn(0); 1123742fe5e2SPeter Brune } 1124742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1125ee78dd50SPeter Brune 1126fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1127fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 112839bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1129293a7e31SPeter Brune } 113039bd7f45SPeter Brune PetscFunctionReturn(0); 113139bd7f45SPeter Brune } 113239bd7f45SPeter Brune 113339bd7f45SPeter Brune #undef __FUNCT__ 113439bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 113539bd7f45SPeter Brune /* 113639bd7f45SPeter Brune 113739bd7f45SPeter Brune The additive cycle looks like: 113839bd7f45SPeter Brune 113907144faaSPeter Brune xhat = x 114007144faaSPeter Brune xhat = dS(x, b) 114107144faaSPeter Brune x = coarsecorrection(xhat, b_d) 114207144faaSPeter Brune x = x + nu*(xhat - x); 114339bd7f45SPeter Brune (optional) x = uS(x, b) 114439bd7f45SPeter Brune 114539bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 114639bd7f45SPeter Brune 114739bd7f45SPeter Brune */ 114839bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 114907144faaSPeter Brune Vec F, B, Xhat; 115022c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 115139bd7f45SPeter Brune PetscErrorCode ierr; 115207144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 115307144faaSPeter Brune SNESConvergedReason reason; 115422c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 115522c1e704SPeter Brune PetscBool lssuccess; 115639bd7f45SPeter Brune PetscFunctionBegin; 115739bd7f45SPeter Brune 115839bd7f45SPeter Brune F = snes->vec_func; 115939bd7f45SPeter Brune B = snes->vec_rhs; 1160fde0ff24SPeter Brune Xhat = snes->work[3]; 116107144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 116207144faaSPeter Brune /* recurse first */ 116307144faaSPeter Brune if (fas->next) { 116407144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 116539bd7f45SPeter Brune 116607144faaSPeter Brune X_c = fas->next->vec_sol; 116707144faaSPeter Brune Xo_c = fas->next->work[0]; 116807144faaSPeter Brune F_c = fas->next->vec_func; 116907144faaSPeter Brune B_c = fas->next->vec_rhs; 117039bd7f45SPeter Brune 1171938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 117207144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 117307144faaSPeter Brune 117407144faaSPeter Brune /* restrict the defect */ 117507144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 117607144faaSPeter Brune 117707144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 117807144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 117907144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 118007144faaSPeter Brune 118107144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 118207144faaSPeter Brune 118307144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 118407144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 118507144faaSPeter Brune 118607144faaSPeter Brune /* recurse */ 118707144faaSPeter Brune fas->next->vec_rhs = B_c; 118807144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 118907144faaSPeter Brune 119007144faaSPeter Brune /* smooth on this level */ 119107144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 119207144faaSPeter Brune 119307144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 119407144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 119507144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 119607144faaSPeter Brune PetscFunctionReturn(0); 119707144faaSPeter Brune } 119807144faaSPeter Brune 119907144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1200c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1201ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 120207144faaSPeter Brune 1203ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1204ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1205ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1206*9e764e56SPeter Brune ierr = PetscLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 1207*9e764e56SPeter Brune ierr = PetscLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 1208*9e764e56SPeter Brune if (!lssuccess) { 1209*9e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1210*9e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 1211*9e764e56SPeter Brune PetscFunctionReturn(0); 1212*9e764e56SPeter Brune } 1213*9e764e56SPeter Brune } 1214*9e764e56SPeter Brune ierr = PetscLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 121507144faaSPeter Brune } else { 121607144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 121707144faaSPeter Brune } 121839bd7f45SPeter Brune PetscFunctionReturn(0); 121939bd7f45SPeter Brune } 122039bd7f45SPeter Brune 122139bd7f45SPeter Brune #undef __FUNCT__ 122239bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 122339bd7f45SPeter Brune /* 122439bd7f45SPeter Brune 122539bd7f45SPeter Brune Defines the FAS cycle as: 122639bd7f45SPeter Brune 122739bd7f45SPeter Brune fine problem: F(x) = 0 122839bd7f45SPeter Brune coarse problem: F^c(x) = b^c 122939bd7f45SPeter Brune 123039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 123139bd7f45SPeter Brune 123239bd7f45SPeter Brune correction: 123339bd7f45SPeter Brune 123439bd7f45SPeter Brune x = x + I(x^c - Rx) 123539bd7f45SPeter Brune 123639bd7f45SPeter Brune */ 123739bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 123839bd7f45SPeter Brune 123939bd7f45SPeter Brune PetscErrorCode ierr; 124039bd7f45SPeter Brune Vec F,B; 124139bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 124239bd7f45SPeter Brune 124339bd7f45SPeter Brune PetscFunctionBegin; 124439bd7f45SPeter Brune F = snes->vec_func; 124539bd7f45SPeter Brune B = snes->vec_rhs; 124639bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 124739bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 124839bd7f45SPeter Brune 124939a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 125039bd7f45SPeter Brune 1251c90fad12SPeter Brune if (fas->level != 0) { 125239bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1253fe6f9142SPeter Brune } 1254fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1255fa9694d7SPeter Brune 1256fa9694d7SPeter Brune PetscFunctionReturn(0); 1257421d9b32SPeter Brune } 1258421d9b32SPeter Brune 1259421d9b32SPeter Brune #undef __FUNCT__ 1260421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1261421d9b32SPeter Brune 1262421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1263421d9b32SPeter Brune { 1264fa9694d7SPeter Brune PetscErrorCode ierr; 1265fe6f9142SPeter Brune PetscInt i, maxits; 1266ddb5aff1SPeter Brune Vec X, F; 1267fe6f9142SPeter Brune PetscReal fnorm; 1268b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1269b17ce1afSJed Brown DM dm; 1270b17ce1afSJed Brown 1271421d9b32SPeter Brune PetscFunctionBegin; 1272fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1273fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1274fa9694d7SPeter Brune X = snes->vec_sol; 1275f5a6d4f9SBarry Smith F = snes->vec_func; 1276293a7e31SPeter Brune 1277293a7e31SPeter Brune /*norm setup */ 1278fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1279fe6f9142SPeter Brune snes->iter = 0; 1280fe6f9142SPeter Brune snes->norm = 0.; 1281fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1282fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1283fe6f9142SPeter Brune if (snes->domainerror) { 1284fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1285fe6f9142SPeter Brune PetscFunctionReturn(0); 1286fe6f9142SPeter Brune } 1287fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1288fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1289fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1290fe6f9142SPeter Brune snes->norm = fnorm; 1291fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1292fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1293fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1294fe6f9142SPeter Brune 1295fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1296fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1297fe6f9142SPeter Brune /* test convergence */ 1298fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1299fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1300b17ce1afSJed Brown 1301b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1302b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1303b17ce1afSJed Brown DM dmcoarse; 1304b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1305b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1306b17ce1afSJed Brown dm = dmcoarse; 1307b17ce1afSJed Brown } 1308b17ce1afSJed Brown 1309fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1310fe6f9142SPeter Brune /* Call general purpose update function */ 1311646217ecSPeter Brune 1312fe6f9142SPeter Brune if (snes->ops->update) { 1313fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1314fe6f9142SPeter Brune } 131507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 131639bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 131707144faaSPeter Brune } else { 131807144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 131907144faaSPeter Brune } 1320742fe5e2SPeter Brune 1321742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1322742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1323742fe5e2SPeter Brune PetscFunctionReturn(0); 1324742fe5e2SPeter Brune } 1325c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1326c90fad12SPeter Brune /* Monitor convergence */ 1327c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1328c90fad12SPeter Brune snes->iter = i+1; 1329c90fad12SPeter Brune snes->norm = fnorm; 1330c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1331c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1332c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1333c90fad12SPeter Brune /* Test for convergence */ 1334c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1335c90fad12SPeter Brune if (snes->reason) break; 1336fe6f9142SPeter Brune } 1337fe6f9142SPeter Brune if (i == maxits) { 1338fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1339fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1340fe6f9142SPeter Brune } 1341421d9b32SPeter Brune PetscFunctionReturn(0); 1342421d9b32SPeter Brune } 1343