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; 6422c1e704SPeter Brune fas->linesearch = PETSC_NULL; 65ddebd997SPeter Brune 66421d9b32SPeter Brune PetscFunctionReturn(0); 67421d9b32SPeter Brune } 68421d9b32SPeter Brune EXTERN_C_END 69421d9b32SPeter Brune 70421d9b32SPeter Brune #undef __FUNCT__ 71421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 72c79ef259SPeter Brune /*@ 732e8ce248SJed Brown SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 74c79ef259SPeter Brune 75c79ef259SPeter Brune Input Parameter: 762e8ce248SJed Brown . snes - the nonlinear solver context 77c79ef259SPeter Brune 78c79ef259SPeter Brune Output parameter: 79c79ef259SPeter Brune . levels - the number of levels 80c79ef259SPeter Brune 81c79ef259SPeter Brune Level: advanced 82c79ef259SPeter Brune 83c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 84c79ef259SPeter Brune 85c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 86c79ef259SPeter Brune @*/ 87421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 88421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 89421d9b32SPeter Brune PetscFunctionBegin; 90ee1fd11aSPeter Brune *levels = fas->levels; 91421d9b32SPeter Brune PetscFunctionReturn(0); 92421d9b32SPeter Brune } 93421d9b32SPeter Brune 94421d9b32SPeter Brune #undef __FUNCT__ 95646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 96c79ef259SPeter Brune /*@ 972e8ce248SJed Brown SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 98c79ef259SPeter Brune complicated cycling. 99c79ef259SPeter Brune 100c79ef259SPeter Brune Logically Collective on SNES 101c79ef259SPeter Brune 102c79ef259SPeter Brune Input Parameters: 103c79ef259SPeter Brune + snes - the multigrid context 104c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 105c79ef259SPeter Brune 106c79ef259SPeter Brune Options Database Key: 107c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 108c79ef259SPeter Brune 109c79ef259SPeter Brune Level: advanced 110c79ef259SPeter Brune 111c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 112c79ef259SPeter Brune 113c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 114c79ef259SPeter Brune @*/ 115646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 116646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 117646217ecSPeter Brune PetscErrorCode ierr; 118646217ecSPeter Brune PetscFunctionBegin; 119646217ecSPeter Brune fas->n_cycles = cycles; 120646217ecSPeter Brune if (fas->next) { 121646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 122646217ecSPeter Brune } 123646217ecSPeter Brune PetscFunctionReturn(0); 124646217ecSPeter Brune } 125646217ecSPeter Brune 126eff52c0eSPeter Brune #undef __FUNCT__ 127c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 128c79ef259SPeter Brune /*@ 129722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 130c79ef259SPeter Brune 131c79ef259SPeter Brune Logically Collective on SNES 132c79ef259SPeter Brune 133c79ef259SPeter Brune Input Parameters: 134c79ef259SPeter Brune + snes - the multigrid context 135c79ef259SPeter Brune . level - the level to set the number of cycles on 136c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 137c79ef259SPeter Brune 138c79ef259SPeter Brune Level: advanced 139c79ef259SPeter Brune 140c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 141c79ef259SPeter Brune 142c79ef259SPeter Brune .seealso: SNESFASSetCycles() 143c79ef259SPeter Brune @*/ 144c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 145c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 146c79ef259SPeter Brune PetscInt top_level = fas->level,i; 147c79ef259SPeter Brune 148c79ef259SPeter Brune PetscFunctionBegin; 1492e8ce248SJed 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); 150c79ef259SPeter Brune /* get to the correct level */ 151c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 152c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 153c79ef259SPeter Brune } 1542e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 155c79ef259SPeter Brune fas->n_cycles = cycles; 156c79ef259SPeter Brune PetscFunctionReturn(0); 157c79ef259SPeter Brune } 158c79ef259SPeter Brune 159c79ef259SPeter Brune #undef __FUNCT__ 160eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 161aeed3662SMatthew G Knepley /*@C 162c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 163c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 164c79ef259SPeter Brune and nonlinear preconditioners. 165c79ef259SPeter Brune 166c79ef259SPeter Brune Logically Collective on SNES 167c79ef259SPeter Brune 168c79ef259SPeter Brune Input Parameters: 169c79ef259SPeter Brune + snes - the multigrid context 170c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 171c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 172c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 173c79ef259SPeter Brune 174c79ef259SPeter Brune Level: advanced 175c79ef259SPeter Brune 176c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 177c79ef259SPeter Brune 178c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 179c79ef259SPeter Brune @*/ 180eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 181eff52c0eSPeter Brune PetscErrorCode ierr = 0; 182eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 183eff52c0eSPeter Brune PetscFunctionBegin; 184eff52c0eSPeter Brune 185eff52c0eSPeter Brune if (gsfunc) { 186eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 187eff52c0eSPeter Brune /* push the provided GS up the tree */ 188eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 1896cab3a1bSJed Brown } else { 1906cab3a1bSJed Brown ierr = SNESGetGS(snes,&gsfunc,&ctx);CHKERRQ(ierr); 1916cab3a1bSJed Brown if (gsfunc) { 192eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 193eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 194eff52c0eSPeter Brune } else if (use_gs) { 1952e8ce248SJed Brown SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level); 196eff52c0eSPeter Brune } 1976cab3a1bSJed Brown } 198eff52c0eSPeter Brune PetscFunctionReturn(0); 199eff52c0eSPeter Brune } 200eff52c0eSPeter Brune 201eff52c0eSPeter Brune #undef __FUNCT__ 202eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 203aeed3662SMatthew G Knepley /*@C 204c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 205c79ef259SPeter Brune 206c79ef259SPeter Brune Logically Collective on SNES 207c79ef259SPeter Brune 208c79ef259SPeter Brune Input Parameters: 209c79ef259SPeter Brune + snes - the multigrid context 210c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 211c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 212c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 213c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 214c79ef259SPeter Brune 215c79ef259SPeter Brune Level: advanced 216c79ef259SPeter Brune 217c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 218c79ef259SPeter Brune 219c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 220c79ef259SPeter Brune @*/ 221eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 222eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 223eff52c0eSPeter Brune PetscErrorCode ierr; 224eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 225eff52c0eSPeter Brune SNES cur_snes = snes; 226eff52c0eSPeter Brune PetscFunctionBegin; 2272e8ce248SJed 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); 228eff52c0eSPeter Brune /* get to the correct level */ 229eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 230eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 231eff52c0eSPeter Brune cur_snes = fas->next; 232eff52c0eSPeter Brune } 2332e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 234eff52c0eSPeter Brune if (gsfunc) { 2356273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 236eff52c0eSPeter Brune } 237eff52c0eSPeter Brune PetscFunctionReturn(0); 238eff52c0eSPeter Brune } 239eff52c0eSPeter Brune 240646217ecSPeter Brune #undef __FUNCT__ 241421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 242c79ef259SPeter Brune /*@ 243c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 244c79ef259SPeter Brune level of the FAS hierarchy. 245c79ef259SPeter Brune 246c79ef259SPeter Brune Input Parameters: 247c79ef259SPeter Brune + snes - the multigrid context 248c79ef259SPeter Brune level - the level to get 249c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 250c79ef259SPeter Brune 251c79ef259SPeter Brune Level: advanced 252c79ef259SPeter Brune 253c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 254c79ef259SPeter Brune 255c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 256c79ef259SPeter Brune @*/ 257421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) { 258421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 259421d9b32SPeter Brune PetscInt i; 2602e8ce248SJed Brown 261421d9b32SPeter Brune PetscFunctionBegin; 2622e8ce248SJed 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); 2632e8ce248SJed 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); 264421d9b32SPeter Brune *lsnes = snes; 265421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 266421d9b32SPeter Brune *lsnes = fas->next; 267421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 268421d9b32SPeter Brune } 2692e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 270421d9b32SPeter Brune PetscFunctionReturn(0); 271421d9b32SPeter Brune } 272421d9b32SPeter Brune 273421d9b32SPeter Brune #undef __FUNCT__ 27407144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 27507144faaSPeter Brune /*@ 27607144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 27707144faaSPeter Brune e 27807144faaSPeter Brune 27907144faaSPeter Brune 28007144faaSPeter Brune @*/ 28107144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 28207144faaSPeter Brune { 28307144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 28407144faaSPeter Brune 28507144faaSPeter Brune PetscFunctionBegin; 28607144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 28707144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 28807144faaSPeter Brune fas->fastype = fastype; 28907144faaSPeter Brune PetscFunctionReturn(0); 29007144faaSPeter Brune } 29107144faaSPeter Brune 29207144faaSPeter Brune #undef __FUNCT__ 293421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 294aeed3662SMatthew G Knepley /*@C 295c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 296c79ef259SPeter Brune Must be called before any other FAS routine. 297c79ef259SPeter Brune 298c79ef259SPeter Brune Input Parameters: 299c79ef259SPeter Brune + snes - the snes context 300c79ef259SPeter Brune . levels - the number of levels 301c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 302c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 303c79ef259SPeter Brune Fortran. 304c79ef259SPeter Brune 305c79ef259SPeter Brune Level: intermediate 306c79ef259SPeter Brune 307c79ef259SPeter Brune Notes: 308c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 309c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 310c79ef259SPeter Brune 311c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 312c79ef259SPeter Brune 313c79ef259SPeter Brune .seealso: SNESFASGetLevels() 314c79ef259SPeter Brune @*/ 315421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 316421d9b32SPeter Brune PetscErrorCode ierr; 317421d9b32SPeter Brune PetscInt i; 318421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3196273346dSPeter Brune SNES prevsnes; 320421d9b32SPeter Brune MPI_Comm comm; 321421d9b32SPeter Brune PetscFunctionBegin; 322ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 323ee1fd11aSPeter Brune if (levels == fas->levels) { 324ee1fd11aSPeter Brune if (!comms) 325ee1fd11aSPeter Brune PetscFunctionReturn(0); 326ee1fd11aSPeter Brune } 327421d9b32SPeter Brune /* user has changed the number of levels; reset */ 328421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 329421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 330421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 331ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3326273346dSPeter Brune fas->previous = PETSC_NULL; 3336273346dSPeter Brune prevsnes = snes; 334421d9b32SPeter Brune /* setup the finest level */ 335421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 336421d9b32SPeter Brune if (comms) comm = comms[i]; 337421d9b32SPeter Brune fas->level = i; 338421d9b32SPeter Brune fas->levels = levels; 339ee1fd11aSPeter Brune fas->next = PETSC_NULL; 340421d9b32SPeter Brune if (i > 0) { 341421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3424a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 343421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 344794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3456273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3466273346dSPeter Brune prevsnes = fas->next; 3476273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 348421d9b32SPeter Brune } 349421d9b32SPeter Brune } 350421d9b32SPeter Brune PetscFunctionReturn(0); 351421d9b32SPeter Brune } 352421d9b32SPeter Brune 353421d9b32SPeter Brune #undef __FUNCT__ 354c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 355c79ef259SPeter Brune /*@ 356c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 357c79ef259SPeter Brune use on all levels. 358c79ef259SPeter Brune 359fde0ff24SPeter Brune Logically Collective on SNES 360c79ef259SPeter Brune 361c79ef259SPeter Brune Input Parameters: 362c79ef259SPeter Brune + snes - the multigrid context 363c79ef259SPeter Brune - n - the number of smoothing steps 364c79ef259SPeter Brune 365c79ef259SPeter Brune Options Database Key: 366d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 367c79ef259SPeter Brune 368c79ef259SPeter Brune Level: advanced 369c79ef259SPeter Brune 370c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 371c79ef259SPeter Brune 372c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 373c79ef259SPeter Brune @*/ 374c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 375c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 376c79ef259SPeter Brune PetscErrorCode ierr = 0; 377c79ef259SPeter Brune PetscFunctionBegin; 378c79ef259SPeter Brune fas->max_up_it = n; 379c79ef259SPeter Brune if (fas->next) { 380c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 381c79ef259SPeter Brune } 382c79ef259SPeter Brune PetscFunctionReturn(0); 383c79ef259SPeter Brune } 384c79ef259SPeter Brune 385c79ef259SPeter Brune #undef __FUNCT__ 386c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 387c79ef259SPeter Brune /*@ 388c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 389c79ef259SPeter Brune use on all levels. 390c79ef259SPeter Brune 391fde0ff24SPeter Brune Logically Collective on SNES 392c79ef259SPeter Brune 393c79ef259SPeter Brune Input Parameters: 394c79ef259SPeter Brune + snes - the multigrid context 395c79ef259SPeter Brune - n - the number of smoothing steps 396c79ef259SPeter Brune 397c79ef259SPeter Brune Options Database Key: 398d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 399c79ef259SPeter Brune 400c79ef259SPeter Brune Level: advanced 401c79ef259SPeter Brune 402c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 403c79ef259SPeter Brune 404c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 405c79ef259SPeter Brune @*/ 406c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 407c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 408c79ef259SPeter Brune PetscErrorCode ierr = 0; 409c79ef259SPeter Brune PetscFunctionBegin; 410c79ef259SPeter Brune fas->max_down_it = n; 411c79ef259SPeter Brune if (fas->next) { 412c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 413c79ef259SPeter Brune } 414c79ef259SPeter Brune PetscFunctionReturn(0); 415c79ef259SPeter Brune } 416c79ef259SPeter Brune 417c79ef259SPeter Brune #undef __FUNCT__ 418421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 419c79ef259SPeter Brune /*@ 420c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 421c79ef259SPeter Brune interpolation from l-1 to the lth level 422c79ef259SPeter Brune 423c79ef259SPeter Brune Input Parameters: 424c79ef259SPeter Brune + snes - the multigrid context 425c79ef259SPeter Brune . mat - the interpolation operator 426c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 427c79ef259SPeter Brune 428c79ef259SPeter Brune Level: advanced 429c79ef259SPeter Brune 430c79ef259SPeter Brune Notes: 431c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 432c79ef259SPeter Brune for the same level. 433c79ef259SPeter Brune 434c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 435c79ef259SPeter Brune out from the matrix size which one it is. 436c79ef259SPeter Brune 437c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 438c79ef259SPeter Brune 439bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 440c79ef259SPeter Brune @*/ 441421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 442421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 443d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 444bccf9bb3SJed Brown PetscErrorCode ierr; 445421d9b32SPeter Brune 446421d9b32SPeter Brune PetscFunctionBegin; 4472e8ce248SJed 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); 448421d9b32SPeter Brune /* get to the correct level */ 449d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 450421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 451421d9b32SPeter Brune } 4522e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 453bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 454bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 455421d9b32SPeter Brune fas->interpolate = mat; 456421d9b32SPeter Brune PetscFunctionReturn(0); 457421d9b32SPeter Brune } 458421d9b32SPeter Brune 459421d9b32SPeter Brune #undef __FUNCT__ 460421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 461c79ef259SPeter Brune /*@ 462c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 463c79ef259SPeter Brune from level l to l-1. 464c79ef259SPeter Brune 465c79ef259SPeter Brune Input Parameters: 466c79ef259SPeter Brune + snes - the multigrid context 467c79ef259SPeter Brune . mat - the restriction matrix 468c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 469c79ef259SPeter Brune 470c79ef259SPeter Brune Level: advanced 471c79ef259SPeter Brune 472c79ef259SPeter Brune Notes: 473c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 474c79ef259SPeter Brune for the same level. 475c79ef259SPeter Brune 476c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 477c79ef259SPeter Brune out from the matrix size which one it is. 478c79ef259SPeter Brune 479fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 480c79ef259SPeter Brune is used. 481c79ef259SPeter Brune 482c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 483c79ef259SPeter Brune 484c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 485c79ef259SPeter Brune @*/ 486421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 487421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 488d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 489bccf9bb3SJed Brown PetscErrorCode ierr; 490421d9b32SPeter Brune 491421d9b32SPeter Brune PetscFunctionBegin; 4922e8ce248SJed 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); 493421d9b32SPeter Brune /* get to the correct level */ 494d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 495421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 496421d9b32SPeter Brune } 4972e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 498bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 499bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 500421d9b32SPeter Brune fas->restrct = mat; 501421d9b32SPeter Brune PetscFunctionReturn(0); 502421d9b32SPeter Brune } 503421d9b32SPeter Brune 504421d9b32SPeter Brune #undef __FUNCT__ 505421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 506c79ef259SPeter Brune /*@ 507c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 508c79ef259SPeter Brune operator from level l to l-1. 509c79ef259SPeter Brune 510c79ef259SPeter Brune Input Parameters: 511c79ef259SPeter Brune + snes - the multigrid context 512c79ef259SPeter Brune . rscale - the restriction scaling 513c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 514c79ef259SPeter Brune 515c79ef259SPeter Brune Level: advanced 516c79ef259SPeter Brune 517c79ef259SPeter Brune Notes: 518c79ef259SPeter Brune This is only used in the case that the injection is not set. 519c79ef259SPeter Brune 520c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 521c79ef259SPeter Brune 522c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 523c79ef259SPeter Brune @*/ 524421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 525421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 526d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 527bccf9bb3SJed Brown PetscErrorCode ierr; 528421d9b32SPeter Brune 529421d9b32SPeter Brune PetscFunctionBegin; 5302e8ce248SJed 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); 531421d9b32SPeter Brune /* get to the correct level */ 532d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 533421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 534421d9b32SPeter Brune } 5352e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 536bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 537bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 538421d9b32SPeter Brune fas->rscale = rscale; 539421d9b32SPeter Brune PetscFunctionReturn(0); 540421d9b32SPeter Brune } 541421d9b32SPeter Brune 542efe1f98aSPeter Brune 543efe1f98aSPeter Brune #undef __FUNCT__ 544efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 545c79ef259SPeter Brune /*@ 546c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 547c79ef259SPeter Brune from level l to l-1. 548c79ef259SPeter Brune 549c79ef259SPeter Brune Input Parameters: 550c79ef259SPeter Brune + snes - the multigrid context 551c79ef259SPeter Brune . mat - the restriction matrix 552c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 553c79ef259SPeter Brune 554c79ef259SPeter Brune Level: advanced 555c79ef259SPeter Brune 556c79ef259SPeter Brune Notes: 557c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 558c79ef259SPeter Brune project the solution instead. 559c79ef259SPeter Brune 560c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 561c79ef259SPeter Brune 562c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 563c79ef259SPeter Brune @*/ 564efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 565efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 566efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 567bccf9bb3SJed Brown PetscErrorCode ierr; 568efe1f98aSPeter Brune 569efe1f98aSPeter Brune PetscFunctionBegin; 5702e8ce248SJed 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); 571efe1f98aSPeter Brune /* get to the correct level */ 572efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 573efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 574efe1f98aSPeter Brune } 5752e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 576bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 577bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 578efe1f98aSPeter Brune fas->inject = mat; 579efe1f98aSPeter Brune PetscFunctionReturn(0); 580efe1f98aSPeter Brune } 581efe1f98aSPeter Brune 582421d9b32SPeter Brune #undef __FUNCT__ 583421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 584421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 585421d9b32SPeter Brune { 58677df8cc4SPeter Brune PetscErrorCode ierr = 0; 587421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 588421d9b32SPeter Brune 589421d9b32SPeter Brune PetscFunctionBegin; 590bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 591bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 5923dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 593bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 594bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 595bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 5963dccd265SPeter Brune 5973dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 5983dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 599742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 60022c1e704SPeter Brune 60122c1e704SPeter Brune ierr = LineSearchDestroy(&fas->linesearch_smooth);CHKERRQ(ierr); 60222c1e704SPeter Brune ierr = LineSearchDestroy(&fas->linesearch);CHKERRQ(ierr); 60322c1e704SPeter Brune 604421d9b32SPeter Brune PetscFunctionReturn(0); 605421d9b32SPeter Brune } 606421d9b32SPeter Brune 607421d9b32SPeter Brune #undef __FUNCT__ 608421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 609421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 610421d9b32SPeter Brune { 611421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 612742fe5e2SPeter Brune PetscErrorCode ierr = 0; 613421d9b32SPeter Brune 614421d9b32SPeter Brune PetscFunctionBegin; 615421d9b32SPeter Brune /* recursively resets and then destroys */ 61679d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 617421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 618421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 619ee1fd11aSPeter Brune 620421d9b32SPeter Brune PetscFunctionReturn(0); 621421d9b32SPeter Brune } 622421d9b32SPeter Brune 623421d9b32SPeter Brune #undef __FUNCT__ 624421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 625421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 626421d9b32SPeter Brune { 62748bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 628421d9b32SPeter Brune PetscErrorCode ierr; 62922c1e704SPeter Brune const char *optionsprefix; 630efe1f98aSPeter Brune VecScatter injscatter; 631d1adcc6fSPeter Brune PetscInt dm_levels; 6323dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 633d1adcc6fSPeter Brune 634421d9b32SPeter Brune PetscFunctionBegin; 635eff52c0eSPeter Brune 636cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 637d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 638d1adcc6fSPeter Brune dm_levels++; 639cc05f883SPeter Brune if (dm_levels > fas->levels) { 6403dccd265SPeter Brune 6412e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6423dccd265SPeter Brune vec_sol = snes->vec_sol; 6433dccd265SPeter Brune vec_func = snes->vec_func; 6443dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6453dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6463dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6473dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6483dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6493dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6503dccd265SPeter Brune 6513dccd265SPeter Brune /* reset the number of levels */ 652d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 653cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6543dccd265SPeter Brune 6553dccd265SPeter Brune snes->vec_sol = vec_sol; 6563dccd265SPeter Brune snes->vec_func = vec_func; 6573dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6583dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 659d1adcc6fSPeter Brune } 660d1adcc6fSPeter Brune } 661d1adcc6fSPeter Brune 6623dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6633dccd265SPeter Brune 66407144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 665fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 66607144faaSPeter Brune } else { 667ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 66807144faaSPeter Brune } 669cc05f883SPeter Brune 67079d9a41aSPeter Brune if (snes->dm) { 67179d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 67279d9a41aSPeter Brune if (fas->next) { 67379d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 67479d9a41aSPeter Brune if (!fas->next->dm) { 67579d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 67679d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 67779d9a41aSPeter Brune } 67879d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 67979d9a41aSPeter Brune if (!fas->interpolate) { 68079d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 681bccf9bb3SJed Brown if (!fas->restrct) { 682bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 68379d9a41aSPeter Brune fas->restrct = fas->interpolate; 68479d9a41aSPeter Brune } 685bccf9bb3SJed Brown } 68679d9a41aSPeter Brune /* set the injection from the DM */ 68779d9a41aSPeter Brune if (!fas->inject) { 68879d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 68979d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 69079d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 69179d9a41aSPeter Brune } 69279d9a41aSPeter Brune } 69379d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 69479d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 69579d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 69679d9a41aSPeter Brune } 69779d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 69879d9a41aSPeter Brune 69979d9a41aSPeter Brune if (fas->next) { 70079d9a41aSPeter Brune if (fas->galerkin) { 70179d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 70279d9a41aSPeter Brune } 70379d9a41aSPeter Brune } 70479d9a41aSPeter Brune 70579d9a41aSPeter Brune if (fas->next) { 70679d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 707938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 708938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 70979d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 71079d9a41aSPeter Brune } 71179d9a41aSPeter Brune 7126cab3a1bSJed Brown /* setup the pre and post smoothers */ 7136cab3a1bSJed Brown if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);} 7146cab3a1bSJed Brown if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);} 715cc05f883SPeter Brune 71622c1e704SPeter Brune /* if the pre and post smoothers don't exist, set up line searches in their place */ 71722c1e704SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 71822c1e704SPeter Brune if (!fas->upsmooth || !fas->downsmooth) { 71922c1e704SPeter Brune ierr = LineSearchCreate(((PetscObject)snes)->comm, &fas->linesearch_smooth);CHKERRQ(ierr); 72022c1e704SPeter Brune ierr = LineSearchSetSNES(fas->linesearch_smooth, snes);CHKERRQ(ierr); 72122c1e704SPeter Brune ierr = LineSearchSetType(fas->linesearch_smooth, LINESEARCHL2);CHKERRQ(ierr); 72222c1e704SPeter Brune ierr = LineSearchAppendOptionsPrefix(fas->linesearch_smooth, "fas_");CHKERRQ(ierr); 72322c1e704SPeter Brune ierr = LineSearchAppendOptionsPrefix(fas->linesearch_smooth, optionsprefix);CHKERRQ(ierr); 72422c1e704SPeter Brune ierr = LineSearchSetFromOptions(fas->linesearch_smooth);CHKERRQ(ierr); 72522c1e704SPeter Brune } 72622c1e704SPeter Brune 72722c1e704SPeter Brune /* set up the default line search for coarse grid corrections */ 72822c1e704SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 72922c1e704SPeter Brune ierr = LineSearchCreate(((PetscObject)snes)->comm, &fas->linesearch);CHKERRQ(ierr); 73022c1e704SPeter Brune ierr = LineSearchSetSNES(fas->linesearch, snes);CHKERRQ(ierr); 73122c1e704SPeter Brune ierr = LineSearchSetType(fas->linesearch, LINESEARCHL2);CHKERRQ(ierr); 73222c1e704SPeter Brune ierr = LineSearchAppendOptionsPrefix(fas->linesearch, optionsprefix);CHKERRQ(ierr); 73322c1e704SPeter Brune ierr = LineSearchSetFromOptions(fas->linesearch);CHKERRQ(ierr); 73422c1e704SPeter Brune } 73522c1e704SPeter Brune 7366273346dSPeter Brune /* setup FAS work vectors */ 7376273346dSPeter Brune if (fas->galerkin) { 7386273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7396273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7406273346dSPeter Brune } 741421d9b32SPeter Brune PetscFunctionReturn(0); 742421d9b32SPeter Brune } 743421d9b32SPeter Brune 744421d9b32SPeter Brune #undef __FUNCT__ 745421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 746421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 747421d9b32SPeter Brune { 748ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 749ee78dd50SPeter Brune PetscInt levels = 1; 750fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 751421d9b32SPeter Brune PetscErrorCode ierr; 752ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 75307144faaSPeter Brune SNESFASType fastype; 754fde0ff24SPeter Brune const char *optionsprefix; 755fde0ff24SPeter Brune const char *prefix; 756421d9b32SPeter Brune 757421d9b32SPeter Brune PetscFunctionBegin; 758c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 759ee78dd50SPeter Brune 760ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 761ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 762ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 763c732cbdbSBarry Smith if (!flg && snes->dm) { 764c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 765c732cbdbSBarry Smith levels++; 766d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 767c732cbdbSBarry Smith } 768ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 769ee78dd50SPeter Brune } 77007144faaSPeter Brune fastype = fas->fastype; 77107144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 77207144faaSPeter Brune if (flg) { 77307144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 77407144faaSPeter Brune } 775ee78dd50SPeter Brune 776fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 777fde0ff24SPeter Brune 778fde0ff24SPeter Brune /* smoother setup options */ 779fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 780fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 781fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 782fde0ff24SPeter Brune if (fas->level == 0) { 783fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 784fde0ff24SPeter Brune } 785ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 786fde0ff24SPeter Brune 787d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 788ee78dd50SPeter Brune 7896273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 7906273346dSPeter Brune 791646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 792ee78dd50SPeter Brune 793ee78dd50SPeter Brune 794421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7958cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 796eff52c0eSPeter Brune 797fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 798fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 799fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 800fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 801fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 802fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 803fde0ff24SPeter Brune } 804fde0ff24SPeter Brune 805fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8068cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8078cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8088cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 809fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8108cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8118cc86e31SPeter Brune } 8128cc86e31SPeter Brune 813fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 81467339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 815ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 81667339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 817fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 818293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 819ee78dd50SPeter Brune } 820fde0ff24SPeter Brune 821fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 822fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 823fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 824fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 825fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 826fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 827fde0ff24SPeter Brune } 828fde0ff24SPeter Brune 829fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 830fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 831fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 832fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 833fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 834fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 835fde0ff24SPeter Brune } 836fde0ff24SPeter Brune 8371a266240SBarry Smith if (fas->upsmooth) { 838fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8391a266240SBarry Smith } 8401a266240SBarry Smith 8411a266240SBarry Smith if (fas->downsmooth) { 842fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8431a266240SBarry Smith } 844ee78dd50SPeter Brune 845742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 846fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 847742fe5e2SPeter Brune } 848742fe5e2SPeter Brune 849ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 850ce11c761SPeter Brune if (!fas->downsmooth) { 85193dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 8525d115551SPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 853ce11c761SPeter Brune if (fas->level == 0) { 85479d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 855ce11c761SPeter Brune } 856ce11c761SPeter Brune } 857ce11c761SPeter Brune 858ce11c761SPeter Brune if (!fas->upsmooth) { 85993dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 8605d115551SPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 861ce11c761SPeter Brune } 862ce11c761SPeter Brune 863ee78dd50SPeter Brune if (monflg) { 864646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 865794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8662f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 867742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 868293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 869293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 870d28543b3SPeter Brune } else { 871d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 872d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 873d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 874d28543b3SPeter Brune } 875ee78dd50SPeter Brune } 876ee78dd50SPeter Brune 877ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 878d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 879421d9b32SPeter Brune PetscFunctionReturn(0); 880421d9b32SPeter Brune } 881421d9b32SPeter Brune 882421d9b32SPeter Brune #undef __FUNCT__ 883421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 884421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 885421d9b32SPeter Brune { 886421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 887421d9b32SPeter Brune PetscBool iascii; 888421d9b32SPeter Brune PetscErrorCode ierr; 889421d9b32SPeter Brune 890421d9b32SPeter Brune PetscFunctionBegin; 891421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 892421d9b32SPeter Brune if (iascii) { 8932e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 894421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 895bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 896ee78dd50SPeter Brune if (fas->upsmooth) { 89739a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 898421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 899ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 900421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 901421d9b32SPeter Brune } else { 902ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 903421d9b32SPeter Brune } 904ee78dd50SPeter Brune if (fas->downsmooth) { 90539a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 906421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 907ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 908421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 909421d9b32SPeter Brune } else { 910ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 911421d9b32SPeter Brune } 912421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 913421d9b32SPeter Brune } else { 914421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 915421d9b32SPeter Brune } 916421d9b32SPeter Brune PetscFunctionReturn(0); 917421d9b32SPeter Brune } 918421d9b32SPeter Brune 919421d9b32SPeter Brune #undef __FUNCT__ 92039bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92139bd7f45SPeter Brune /* 92239bd7f45SPeter Brune Defines the action of the downsmoother 92339bd7f45SPeter Brune */ 92439bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 92539bd7f45SPeter Brune PetscErrorCode ierr = 0; 92622c1e704SPeter Brune PetscReal fnorm; 927742fe5e2SPeter Brune SNESConvergedReason reason; 92839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 92922c1e704SPeter Brune Vec Y, FPC; 930fde0ff24SPeter Brune PetscBool lssuccess; 931fde0ff24SPeter Brune PetscInt k; 932421d9b32SPeter Brune PetscFunctionBegin; 933fde0ff24SPeter Brune Y = snes->work[3]; 934d1adcc6fSPeter Brune if (fas->downsmooth) { 935d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 936742fe5e2SPeter Brune /* check convergence reason for the smoother */ 937742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 938742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 939742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 940742fe5e2SPeter Brune PetscFunctionReturn(0); 941742fe5e2SPeter Brune } 9424b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9434b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 944fde0ff24SPeter Brune } else { 945fde0ff24SPeter Brune /* basic richardson smoother */ 946fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 947794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 948794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 949fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 95022c1e704SPeter Brune ierr = LineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 95122c1e704SPeter Brune ierr = LineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 952fde0ff24SPeter Brune if (!lssuccess) { 953fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 954fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9552f7ea302SPeter Brune PetscFunctionReturn(0); 9562f7ea302SPeter Brune } 957fe6f9142SPeter Brune } 958fde0ff24SPeter Brune } 959fde0ff24SPeter Brune } 96039bd7f45SPeter Brune PetscFunctionReturn(0); 96139bd7f45SPeter Brune } 96239bd7f45SPeter Brune 96339bd7f45SPeter Brune 96439bd7f45SPeter Brune #undef __FUNCT__ 96539bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 96639bd7f45SPeter Brune /* 96707144faaSPeter Brune Defines the action of the upsmoother 96839bd7f45SPeter Brune */ 96939bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 97039bd7f45SPeter Brune PetscErrorCode ierr = 0; 97122c1e704SPeter Brune PetscReal fnorm; 97239bd7f45SPeter Brune SNESConvergedReason reason; 97339bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 97422c1e704SPeter Brune Vec Y, FPC; 975fde0ff24SPeter Brune PetscBool lssuccess; 976fde0ff24SPeter Brune PetscInt k; 97739bd7f45SPeter Brune PetscFunctionBegin; 978fde0ff24SPeter Brune Y = snes->work[3]; 97939bd7f45SPeter Brune if (fas->upsmooth) { 980fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 98139bd7f45SPeter Brune /* check convergence reason for the smoother */ 982fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 98339bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 98439bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 98539bd7f45SPeter Brune PetscFunctionReturn(0); 98639bd7f45SPeter Brune } 9874b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9884b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 989fde0ff24SPeter Brune } else { 990fde0ff24SPeter Brune /* basic richardson smoother */ 991fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 99239bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 99339bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 994fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 99522c1e704SPeter Brune ierr = LineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 99622c1e704SPeter Brune ierr = LineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 997fde0ff24SPeter Brune if (!lssuccess) { 998fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 999fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 100039bd7f45SPeter Brune PetscFunctionReturn(0); 100139bd7f45SPeter Brune } 100239bd7f45SPeter Brune } 1003fde0ff24SPeter Brune } 1004fde0ff24SPeter Brune } 100539bd7f45SPeter Brune PetscFunctionReturn(0); 100639bd7f45SPeter Brune } 100739bd7f45SPeter Brune 100839bd7f45SPeter Brune #undef __FUNCT__ 1009938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1010938e4a01SJed Brown /*@ 1011938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1012938e4a01SJed Brown 1013938e4a01SJed Brown Collective 1014938e4a01SJed Brown 1015938e4a01SJed Brown Input Arguments: 1016938e4a01SJed Brown . snes - SNESFAS 1017938e4a01SJed Brown 1018938e4a01SJed Brown Output Arguments: 1019938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1020938e4a01SJed Brown 1021938e4a01SJed Brown Level: developer 1022938e4a01SJed Brown 1023938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1024938e4a01SJed Brown @*/ 1025938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1026938e4a01SJed Brown { 1027938e4a01SJed Brown PetscErrorCode ierr; 1028938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1029938e4a01SJed Brown 1030938e4a01SJed Brown PetscFunctionBegin; 1031938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1032938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1033938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1034938e4a01SJed Brown PetscFunctionReturn(0); 1035938e4a01SJed Brown } 1036938e4a01SJed Brown 1037e9923e8dSJed Brown #undef __FUNCT__ 1038e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1039e9923e8dSJed Brown /*@ 1040e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1041e9923e8dSJed Brown 1042e9923e8dSJed Brown Collective 1043e9923e8dSJed Brown 1044e9923e8dSJed Brown Input Arguments: 1045e9923e8dSJed Brown + fine - SNES from which to restrict 1046e9923e8dSJed Brown - Xfine - vector to restrict 1047e9923e8dSJed Brown 1048e9923e8dSJed Brown Output Arguments: 1049e9923e8dSJed Brown . Xcoarse - result of restriction 1050e9923e8dSJed Brown 1051e9923e8dSJed Brown Level: developer 1052e9923e8dSJed Brown 1053e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1054e9923e8dSJed Brown @*/ 1055e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1056e9923e8dSJed Brown { 1057e9923e8dSJed Brown PetscErrorCode ierr; 1058e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1059e9923e8dSJed Brown 1060e9923e8dSJed Brown PetscFunctionBegin; 1061e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1062e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1063e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1064e9923e8dSJed Brown if (fas->inject) { 1065e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1066e9923e8dSJed Brown } else { 1067e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1068e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1069e9923e8dSJed Brown } 1070e9923e8dSJed Brown PetscFunctionReturn(0); 1071e9923e8dSJed Brown } 1072e9923e8dSJed Brown 1073e9923e8dSJed Brown #undef __FUNCT__ 107439bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 107539bd7f45SPeter Brune /* 107639bd7f45SPeter Brune 107739bd7f45SPeter Brune Performs the FAS coarse correction as: 107839bd7f45SPeter Brune 107939bd7f45SPeter Brune fine problem: F(x) = 0 108039bd7f45SPeter Brune coarse problem: F^c(x) = b^c 108139bd7f45SPeter Brune 108239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 108339bd7f45SPeter Brune 108439bd7f45SPeter Brune with correction: 108539bd7f45SPeter Brune 108639bd7f45SPeter Brune 108739bd7f45SPeter Brune 108839bd7f45SPeter Brune */ 108939a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 109039bd7f45SPeter Brune PetscErrorCode ierr; 109139bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 109239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 109339bd7f45SPeter Brune SNESConvergedReason reason; 109439bd7f45SPeter Brune PetscFunctionBegin; 1095fa9694d7SPeter Brune if (fas->next) { 1096c90fad12SPeter Brune X_c = fas->next->vec_sol; 1097293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1098c90fad12SPeter Brune F_c = fas->next->vec_func; 1099742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1100efe1f98aSPeter Brune 1101938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1102293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1103293a7e31SPeter Brune 1104293a7e31SPeter Brune /* restrict the defect */ 1105293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1106293a7e31SPeter Brune 1107c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1108ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1109c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1110742fe5e2SPeter Brune 1111293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1112c90fad12SPeter Brune 1113ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1114ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1115c90fad12SPeter Brune 1116c90fad12SPeter Brune /* recurse to the next level */ 1117f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1118742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1119742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1120742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1121742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1122742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1123742fe5e2SPeter Brune PetscFunctionReturn(0); 1124742fe5e2SPeter Brune } 1125742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1126ee78dd50SPeter Brune 1127fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1128fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 112939bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1130293a7e31SPeter Brune } 113139bd7f45SPeter Brune PetscFunctionReturn(0); 113239bd7f45SPeter Brune } 113339bd7f45SPeter Brune 113439bd7f45SPeter Brune #undef __FUNCT__ 113539bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 113639bd7f45SPeter Brune /* 113739bd7f45SPeter Brune 113839bd7f45SPeter Brune The additive cycle looks like: 113939bd7f45SPeter Brune 114007144faaSPeter Brune xhat = x 114107144faaSPeter Brune xhat = dS(x, b) 114207144faaSPeter Brune x = coarsecorrection(xhat, b_d) 114307144faaSPeter Brune x = x + nu*(xhat - x); 114439bd7f45SPeter Brune (optional) x = uS(x, b) 114539bd7f45SPeter Brune 114639bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 114739bd7f45SPeter Brune 114839bd7f45SPeter Brune */ 114939bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 115007144faaSPeter Brune Vec F, B, Xhat; 115122c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 115239bd7f45SPeter Brune PetscErrorCode ierr; 115307144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 115407144faaSPeter Brune SNESConvergedReason reason; 115522c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 115622c1e704SPeter Brune PetscBool lssuccess; 115739bd7f45SPeter Brune PetscFunctionBegin; 115839bd7f45SPeter Brune 115939bd7f45SPeter Brune F = snes->vec_func; 116039bd7f45SPeter Brune B = snes->vec_rhs; 1161fde0ff24SPeter Brune Xhat = snes->work[3]; 116207144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 116307144faaSPeter Brune /* recurse first */ 116407144faaSPeter Brune if (fas->next) { 116507144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 116639bd7f45SPeter Brune 116707144faaSPeter Brune X_c = fas->next->vec_sol; 116807144faaSPeter Brune Xo_c = fas->next->work[0]; 116907144faaSPeter Brune F_c = fas->next->vec_func; 117007144faaSPeter Brune B_c = fas->next->vec_rhs; 117139bd7f45SPeter Brune 1172938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 117307144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 117407144faaSPeter Brune 117507144faaSPeter Brune /* restrict the defect */ 117607144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 117707144faaSPeter Brune 117807144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 117907144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 118007144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 118107144faaSPeter Brune 118207144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 118307144faaSPeter Brune 118407144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 118507144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 118607144faaSPeter Brune 118707144faaSPeter Brune /* recurse */ 118807144faaSPeter Brune fas->next->vec_rhs = B_c; 118907144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 119007144faaSPeter Brune 119107144faaSPeter Brune /* smooth on this level */ 119207144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 119307144faaSPeter Brune 119407144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 119507144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 119607144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 119707144faaSPeter Brune PetscFunctionReturn(0); 119807144faaSPeter Brune } 119907144faaSPeter Brune 120007144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1201*c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1202ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 120307144faaSPeter Brune 1204ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1205ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1206ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 120722c1e704SPeter Brune ierr = LineSearchApply(fas->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 120822c1e704SPeter Brune ierr = LineSearchGetSuccess(fas->linesearch, &lssuccess);CHKERRQ(ierr); 120922c1e704SPeter Brune ierr = LineSearchGetNorms(fas->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 121007144faaSPeter Brune } else { 121107144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 121207144faaSPeter Brune } 121339bd7f45SPeter Brune PetscFunctionReturn(0); 121439bd7f45SPeter Brune } 121539bd7f45SPeter Brune 121639bd7f45SPeter Brune #undef __FUNCT__ 121739bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 121839bd7f45SPeter Brune /* 121939bd7f45SPeter Brune 122039bd7f45SPeter Brune Defines the FAS cycle as: 122139bd7f45SPeter Brune 122239bd7f45SPeter Brune fine problem: F(x) = 0 122339bd7f45SPeter Brune coarse problem: F^c(x) = b^c 122439bd7f45SPeter Brune 122539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 122639bd7f45SPeter Brune 122739bd7f45SPeter Brune correction: 122839bd7f45SPeter Brune 122939bd7f45SPeter Brune x = x + I(x^c - Rx) 123039bd7f45SPeter Brune 123139bd7f45SPeter Brune */ 123239bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 123339bd7f45SPeter Brune 123439bd7f45SPeter Brune PetscErrorCode ierr; 123539bd7f45SPeter Brune Vec F,B; 123639bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 123739bd7f45SPeter Brune 123839bd7f45SPeter Brune PetscFunctionBegin; 123939bd7f45SPeter Brune F = snes->vec_func; 124039bd7f45SPeter Brune B = snes->vec_rhs; 124139bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 124239bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 124339bd7f45SPeter Brune 124439a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 124539bd7f45SPeter Brune 1246c90fad12SPeter Brune if (fas->level != 0) { 124739bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1248fe6f9142SPeter Brune } 1249fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1250fa9694d7SPeter Brune 1251fa9694d7SPeter Brune PetscFunctionReturn(0); 1252421d9b32SPeter Brune } 1253421d9b32SPeter Brune 1254421d9b32SPeter Brune #undef __FUNCT__ 1255421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1256421d9b32SPeter Brune 1257421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1258421d9b32SPeter Brune { 1259fa9694d7SPeter Brune PetscErrorCode ierr; 1260fe6f9142SPeter Brune PetscInt i, maxits; 1261ddb5aff1SPeter Brune Vec X, F; 1262fe6f9142SPeter Brune PetscReal fnorm; 1263b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1264b17ce1afSJed Brown DM dm; 1265b17ce1afSJed Brown 1266421d9b32SPeter Brune PetscFunctionBegin; 1267fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1268fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1269fa9694d7SPeter Brune X = snes->vec_sol; 1270f5a6d4f9SBarry Smith F = snes->vec_func; 1271293a7e31SPeter Brune 1272293a7e31SPeter Brune /*norm setup */ 1273fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1274fe6f9142SPeter Brune snes->iter = 0; 1275fe6f9142SPeter Brune snes->norm = 0.; 1276fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1277fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1278fe6f9142SPeter Brune if (snes->domainerror) { 1279fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1280fe6f9142SPeter Brune PetscFunctionReturn(0); 1281fe6f9142SPeter Brune } 1282fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1283fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1284fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1285fe6f9142SPeter Brune snes->norm = fnorm; 1286fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1287fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1288fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1289fe6f9142SPeter Brune 1290fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1291fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1292fe6f9142SPeter Brune /* test convergence */ 1293fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1294fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1295b17ce1afSJed Brown 1296b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1297b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1298b17ce1afSJed Brown DM dmcoarse; 1299b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1300b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1301b17ce1afSJed Brown dm = dmcoarse; 1302b17ce1afSJed Brown } 1303b17ce1afSJed Brown 1304fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1305fe6f9142SPeter Brune /* Call general purpose update function */ 1306646217ecSPeter Brune 1307fe6f9142SPeter Brune if (snes->ops->update) { 1308fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1309fe6f9142SPeter Brune } 131007144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 131139bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 131207144faaSPeter Brune } else { 131307144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 131407144faaSPeter Brune } 1315742fe5e2SPeter Brune 1316742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1317742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1318742fe5e2SPeter Brune PetscFunctionReturn(0); 1319742fe5e2SPeter Brune } 1320c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1321c90fad12SPeter Brune /* Monitor convergence */ 1322c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1323c90fad12SPeter Brune snes->iter = i+1; 1324c90fad12SPeter Brune snes->norm = fnorm; 1325c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1326c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1327c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1328c90fad12SPeter Brune /* Test for convergence */ 1329c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1330c90fad12SPeter Brune if (snes->reason) break; 1331fe6f9142SPeter Brune } 1332fe6f9142SPeter Brune if (i == maxits) { 1333fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1334fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1335fe6f9142SPeter Brune } 1336421d9b32SPeter Brune PetscFunctionReturn(0); 1337421d9b32SPeter Brune } 1338