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 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 13421d9b32SPeter Brune 141fbfccc6SJed Brown /*MC 151fbfccc6SJed Brown 161fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 171fbfccc6SJed Brown 181fbfccc6SJed Brown The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 191fbfccc6SJed Brown 201fbfccc6SJed Brown Level: advanced 211fbfccc6SJed Brown 221fbfccc6SJed Brown .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 231fbfccc6SJed Brown M*/ 24421d9b32SPeter Brune 25421d9b32SPeter Brune #undef __FUNCT__ 26421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 271fbfccc6SJed Brown PETSC_EXTERN_C 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 43*88976e71SPeter Brune if (!snes->tolerancesset) { 440e444f03SPeter Brune snes->max_funcs = 30000; 450e444f03SPeter Brune snes->max_its = 10000; 46*88976e71SPeter Brune } 470e444f03SPeter Brune 48421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 49421d9b32SPeter Brune snes->data = (void*) fas; 50421d9b32SPeter Brune fas->level = 0; 51293a7e31SPeter Brune fas->levels = 1; 52ee78dd50SPeter Brune fas->n_cycles = 1; 53ee78dd50SPeter Brune fas->max_up_it = 1; 54ee78dd50SPeter Brune fas->max_down_it = 1; 55ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 56ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 57421d9b32SPeter Brune fas->next = PETSC_NULL; 586273346dSPeter Brune fas->previous = PETSC_NULL; 59421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 60421d9b32SPeter Brune fas->restrct = PETSC_NULL; 61efe1f98aSPeter Brune fas->inject = PETSC_NULL; 62cc05f883SPeter Brune fas->monitor = PETSC_NULL; 63cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 64ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 6522c1e704SPeter Brune fas->linesearch_smooth = PETSC_NULL; 66ddebd997SPeter Brune 67421d9b32SPeter Brune PetscFunctionReturn(0); 68421d9b32SPeter Brune } 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 2781fbfccc6SJed Brown Logically Collective 27907144faaSPeter Brune 2801fbfccc6SJed Brown Input Arguments: 2811fbfccc6SJed Brown snes - nonlinear solver 2821fbfccc6SJed Brown fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 2831fbfccc6SJed Brown 2841fbfccc6SJed Brown .seealso: PCMGSetType() 28507144faaSPeter Brune @*/ 28607144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 28707144faaSPeter Brune { 28807144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 28907144faaSPeter Brune 29007144faaSPeter Brune PetscFunctionBegin; 29107144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 29207144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 29307144faaSPeter Brune fas->fastype = fastype; 29407144faaSPeter Brune PetscFunctionReturn(0); 29507144faaSPeter Brune } 29607144faaSPeter Brune 29707144faaSPeter Brune #undef __FUNCT__ 298421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 299aeed3662SMatthew G Knepley /*@C 300c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 301c79ef259SPeter Brune Must be called before any other FAS routine. 302c79ef259SPeter Brune 303c79ef259SPeter Brune Input Parameters: 304c79ef259SPeter Brune + snes - the snes context 305c79ef259SPeter Brune . levels - the number of levels 306c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 307c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 308c79ef259SPeter Brune Fortran. 309c79ef259SPeter Brune 310c79ef259SPeter Brune Level: intermediate 311c79ef259SPeter Brune 312c79ef259SPeter Brune Notes: 313c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 314c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 315c79ef259SPeter Brune 316c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 317c79ef259SPeter Brune 318c79ef259SPeter Brune .seealso: SNESFASGetLevels() 319c79ef259SPeter Brune @*/ 320421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 321421d9b32SPeter Brune PetscErrorCode ierr; 322421d9b32SPeter Brune PetscInt i; 323421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3246273346dSPeter Brune SNES prevsnes; 325421d9b32SPeter Brune MPI_Comm comm; 326421d9b32SPeter Brune PetscFunctionBegin; 327ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 328ee1fd11aSPeter Brune if (levels == fas->levels) { 329ee1fd11aSPeter Brune if (!comms) 330ee1fd11aSPeter Brune PetscFunctionReturn(0); 331ee1fd11aSPeter Brune } 332421d9b32SPeter Brune /* user has changed the number of levels; reset */ 333421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 334421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 335421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 336ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3376273346dSPeter Brune fas->previous = PETSC_NULL; 3386273346dSPeter Brune prevsnes = snes; 339421d9b32SPeter Brune /* setup the finest level */ 340421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 341421d9b32SPeter Brune if (comms) comm = comms[i]; 342421d9b32SPeter Brune fas->level = i; 343421d9b32SPeter Brune fas->levels = levels; 344ee1fd11aSPeter Brune fas->next = PETSC_NULL; 345421d9b32SPeter Brune if (i > 0) { 346421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3474a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 348421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 349794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3506273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3516273346dSPeter Brune prevsnes = fas->next; 3526273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 353421d9b32SPeter Brune } 354421d9b32SPeter Brune } 355421d9b32SPeter Brune PetscFunctionReturn(0); 356421d9b32SPeter Brune } 357421d9b32SPeter Brune 358421d9b32SPeter Brune #undef __FUNCT__ 359c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 360c79ef259SPeter Brune /*@ 361c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 362c79ef259SPeter Brune use on all levels. 363c79ef259SPeter Brune 364fde0ff24SPeter Brune Logically Collective on SNES 365c79ef259SPeter Brune 366c79ef259SPeter Brune Input Parameters: 367c79ef259SPeter Brune + snes - the multigrid context 368c79ef259SPeter Brune - n - the number of smoothing steps 369c79ef259SPeter Brune 370c79ef259SPeter Brune Options Database Key: 371d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 372c79ef259SPeter Brune 373c79ef259SPeter Brune Level: advanced 374c79ef259SPeter Brune 375c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 376c79ef259SPeter Brune 377c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 378c79ef259SPeter Brune @*/ 379c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 380c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 381c79ef259SPeter Brune PetscErrorCode ierr = 0; 382c79ef259SPeter Brune PetscFunctionBegin; 383c79ef259SPeter Brune fas->max_up_it = n; 384c79ef259SPeter Brune if (fas->next) { 385c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 386c79ef259SPeter Brune } 387c79ef259SPeter Brune PetscFunctionReturn(0); 388c79ef259SPeter Brune } 389c79ef259SPeter Brune 390c79ef259SPeter Brune #undef __FUNCT__ 391c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 392c79ef259SPeter Brune /*@ 393c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 394c79ef259SPeter Brune use on all levels. 395c79ef259SPeter Brune 396fde0ff24SPeter Brune Logically Collective on SNES 397c79ef259SPeter Brune 398c79ef259SPeter Brune Input Parameters: 399c79ef259SPeter Brune + snes - the multigrid context 400c79ef259SPeter Brune - n - the number of smoothing steps 401c79ef259SPeter Brune 402c79ef259SPeter Brune Options Database Key: 403d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 404c79ef259SPeter Brune 405c79ef259SPeter Brune Level: advanced 406c79ef259SPeter Brune 407c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 408c79ef259SPeter Brune 409c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 410c79ef259SPeter Brune @*/ 411c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 412c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 413c79ef259SPeter Brune PetscErrorCode ierr = 0; 414c79ef259SPeter Brune PetscFunctionBegin; 415c79ef259SPeter Brune fas->max_down_it = n; 416c79ef259SPeter Brune if (fas->next) { 417c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 418c79ef259SPeter Brune } 419c79ef259SPeter Brune PetscFunctionReturn(0); 420c79ef259SPeter Brune } 421c79ef259SPeter Brune 422c79ef259SPeter Brune #undef __FUNCT__ 423421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 424c79ef259SPeter Brune /*@ 425c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 426c79ef259SPeter Brune interpolation from l-1 to the lth level 427c79ef259SPeter Brune 428c79ef259SPeter Brune Input Parameters: 429c79ef259SPeter Brune + snes - the multigrid context 430c79ef259SPeter Brune . mat - the interpolation operator 431c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 432c79ef259SPeter Brune 433c79ef259SPeter Brune Level: advanced 434c79ef259SPeter Brune 435c79ef259SPeter Brune Notes: 436c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 437c79ef259SPeter Brune for the same level. 438c79ef259SPeter Brune 439c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 440c79ef259SPeter Brune out from the matrix size which one it is. 441c79ef259SPeter Brune 442c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 443c79ef259SPeter Brune 444bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 445c79ef259SPeter Brune @*/ 446421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 447421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 448d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 449bccf9bb3SJed Brown PetscErrorCode ierr; 450421d9b32SPeter Brune 451421d9b32SPeter Brune PetscFunctionBegin; 4522e8ce248SJed 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); 453421d9b32SPeter Brune /* get to the correct level */ 454d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 455421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 456421d9b32SPeter Brune } 4572e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 458bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 459bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 460421d9b32SPeter Brune fas->interpolate = mat; 461421d9b32SPeter Brune PetscFunctionReturn(0); 462421d9b32SPeter Brune } 463421d9b32SPeter Brune 464421d9b32SPeter Brune #undef __FUNCT__ 465421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 466c79ef259SPeter Brune /*@ 467c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 468c79ef259SPeter Brune from level l to l-1. 469c79ef259SPeter Brune 470c79ef259SPeter Brune Input Parameters: 471c79ef259SPeter Brune + snes - the multigrid context 472c79ef259SPeter Brune . mat - the restriction matrix 473c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 474c79ef259SPeter Brune 475c79ef259SPeter Brune Level: advanced 476c79ef259SPeter Brune 477c79ef259SPeter Brune Notes: 478c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 479c79ef259SPeter Brune for the same level. 480c79ef259SPeter Brune 481c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 482c79ef259SPeter Brune out from the matrix size which one it is. 483c79ef259SPeter Brune 484fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 485c79ef259SPeter Brune is used. 486c79ef259SPeter Brune 487c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 488c79ef259SPeter Brune 489c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 490c79ef259SPeter Brune @*/ 491421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 492421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 493d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 494bccf9bb3SJed Brown PetscErrorCode ierr; 495421d9b32SPeter Brune 496421d9b32SPeter Brune PetscFunctionBegin; 4972e8ce248SJed 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); 498421d9b32SPeter Brune /* get to the correct level */ 499d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 500421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 501421d9b32SPeter Brune } 5022e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 503bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 504bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 505421d9b32SPeter Brune fas->restrct = mat; 506421d9b32SPeter Brune PetscFunctionReturn(0); 507421d9b32SPeter Brune } 508421d9b32SPeter Brune 509421d9b32SPeter Brune #undef __FUNCT__ 510421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 511c79ef259SPeter Brune /*@ 512c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 513c79ef259SPeter Brune operator from level l to l-1. 514c79ef259SPeter Brune 515c79ef259SPeter Brune Input Parameters: 516c79ef259SPeter Brune + snes - the multigrid context 517c79ef259SPeter Brune . rscale - the restriction scaling 518c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 519c79ef259SPeter Brune 520c79ef259SPeter Brune Level: advanced 521c79ef259SPeter Brune 522c79ef259SPeter Brune Notes: 523c79ef259SPeter Brune This is only used in the case that the injection is not set. 524c79ef259SPeter Brune 525c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 526c79ef259SPeter Brune 527c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 528c79ef259SPeter Brune @*/ 529421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 530421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 531d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 532bccf9bb3SJed Brown PetscErrorCode ierr; 533421d9b32SPeter Brune 534421d9b32SPeter Brune PetscFunctionBegin; 5352e8ce248SJed 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); 536421d9b32SPeter Brune /* get to the correct level */ 537d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 538421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 539421d9b32SPeter Brune } 5402e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 541bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 542bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 543421d9b32SPeter Brune fas->rscale = rscale; 544421d9b32SPeter Brune PetscFunctionReturn(0); 545421d9b32SPeter Brune } 546421d9b32SPeter Brune 547efe1f98aSPeter Brune 548efe1f98aSPeter Brune #undef __FUNCT__ 549efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 550c79ef259SPeter Brune /*@ 551c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 552c79ef259SPeter Brune from level l to l-1. 553c79ef259SPeter Brune 554c79ef259SPeter Brune Input Parameters: 555c79ef259SPeter Brune + snes - the multigrid context 556c79ef259SPeter Brune . mat - the restriction matrix 557c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 558c79ef259SPeter Brune 559c79ef259SPeter Brune Level: advanced 560c79ef259SPeter Brune 561c79ef259SPeter Brune Notes: 562c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 563c79ef259SPeter Brune project the solution instead. 564c79ef259SPeter Brune 565c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 566c79ef259SPeter Brune 567c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 568c79ef259SPeter Brune @*/ 569efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 570efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 571efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 572bccf9bb3SJed Brown PetscErrorCode ierr; 573efe1f98aSPeter Brune 574efe1f98aSPeter Brune PetscFunctionBegin; 5752e8ce248SJed 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); 576efe1f98aSPeter Brune /* get to the correct level */ 577efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 578efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 579efe1f98aSPeter Brune } 5802e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 581bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 582bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 583efe1f98aSPeter Brune fas->inject = mat; 584efe1f98aSPeter Brune PetscFunctionReturn(0); 585efe1f98aSPeter Brune } 586efe1f98aSPeter Brune 587421d9b32SPeter Brune #undef __FUNCT__ 588421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 589421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 590421d9b32SPeter Brune { 59177df8cc4SPeter Brune PetscErrorCode ierr = 0; 592421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 593421d9b32SPeter Brune 594421d9b32SPeter Brune PetscFunctionBegin; 595bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 596bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 5973dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 598bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 599bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 600bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6013dccd265SPeter Brune 6023dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6033dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 604742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 60522c1e704SPeter Brune 606f1c6b773SPeter Brune ierr = SNESLineSearchDestroy(&fas->linesearch_smooth);CHKERRQ(ierr); 60722c1e704SPeter Brune 608421d9b32SPeter Brune PetscFunctionReturn(0); 609421d9b32SPeter Brune } 610421d9b32SPeter Brune 611421d9b32SPeter Brune #undef __FUNCT__ 612421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 613421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 614421d9b32SPeter Brune { 615421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 616742fe5e2SPeter Brune PetscErrorCode ierr = 0; 617421d9b32SPeter Brune 618421d9b32SPeter Brune PetscFunctionBegin; 619421d9b32SPeter Brune /* recursively resets and then destroys */ 62079d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 621421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 622421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 623ee1fd11aSPeter Brune 624421d9b32SPeter Brune PetscFunctionReturn(0); 625421d9b32SPeter Brune } 626421d9b32SPeter Brune 627421d9b32SPeter Brune #undef __FUNCT__ 628421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 629421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 630421d9b32SPeter Brune { 63148bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 632421d9b32SPeter Brune PetscErrorCode ierr; 63322c1e704SPeter Brune const char *optionsprefix; 634efe1f98aSPeter Brune VecScatter injscatter; 635d1adcc6fSPeter Brune PetscInt dm_levels; 6363dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 637d1adcc6fSPeter Brune 638421d9b32SPeter Brune PetscFunctionBegin; 639eff52c0eSPeter Brune 640cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 641d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 642d1adcc6fSPeter Brune dm_levels++; 643cc05f883SPeter Brune if (dm_levels > fas->levels) { 6443dccd265SPeter Brune 6452e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6463dccd265SPeter Brune vec_sol = snes->vec_sol; 6473dccd265SPeter Brune vec_func = snes->vec_func; 6483dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6493dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6503dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6513dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6523dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6533dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6543dccd265SPeter Brune 6553dccd265SPeter Brune /* reset the number of levels */ 656d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 657cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6583dccd265SPeter Brune 6593dccd265SPeter Brune snes->vec_sol = vec_sol; 6603dccd265SPeter Brune snes->vec_func = vec_func; 6613dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6623dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 663d1adcc6fSPeter Brune } 664d1adcc6fSPeter Brune } 665d1adcc6fSPeter Brune 6663dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6673dccd265SPeter Brune 66807144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 669fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 67007144faaSPeter Brune } else { 671ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 67207144faaSPeter Brune } 673cc05f883SPeter Brune 67479d9a41aSPeter Brune if (snes->dm) { 67579d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 67679d9a41aSPeter Brune if (fas->next) { 67779d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 67879d9a41aSPeter Brune if (!fas->next->dm) { 67979d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 68079d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 68179d9a41aSPeter Brune } 68279d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 68379d9a41aSPeter Brune if (!fas->interpolate) { 68479d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 685bccf9bb3SJed Brown if (!fas->restrct) { 686bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 68779d9a41aSPeter Brune fas->restrct = fas->interpolate; 68879d9a41aSPeter Brune } 689bccf9bb3SJed Brown } 69079d9a41aSPeter Brune /* set the injection from the DM */ 69179d9a41aSPeter Brune if (!fas->inject) { 69279d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 69379d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 69479d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 69579d9a41aSPeter Brune } 69679d9a41aSPeter Brune } 69779d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 69879d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 69979d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 70079d9a41aSPeter Brune } 70179d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 70279d9a41aSPeter Brune 70379d9a41aSPeter Brune if (fas->next) { 70479d9a41aSPeter Brune if (fas->galerkin) { 70579d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 70679d9a41aSPeter Brune } 70779d9a41aSPeter Brune } 70879d9a41aSPeter Brune 70979d9a41aSPeter Brune if (fas->next) { 71079d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 711938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 712938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 71379d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 71479d9a41aSPeter Brune } 71579d9a41aSPeter Brune 7166cab3a1bSJed Brown /* setup the pre and post smoothers */ 7176cab3a1bSJed Brown if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);} 7186cab3a1bSJed Brown if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);} 719cc05f883SPeter Brune 72022c1e704SPeter Brune /* if the pre and post smoothers don't exist, set up line searches in their place */ 72122c1e704SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 72222c1e704SPeter Brune if (!fas->upsmooth || !fas->downsmooth) { 723f1c6b773SPeter Brune ierr = SNESLineSearchCreate(((PetscObject)snes)->comm, &fas->linesearch_smooth);CHKERRQ(ierr); 724f1c6b773SPeter Brune ierr = SNESLineSearchSetSNES(fas->linesearch_smooth, snes);CHKERRQ(ierr); 725f1c6b773SPeter Brune ierr = SNESLineSearchSetType(fas->linesearch_smooth, SNES_LINESEARCH_L2);CHKERRQ(ierr); 726f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(fas->linesearch_smooth, "fas_");CHKERRQ(ierr); 727f1c6b773SPeter Brune ierr = SNESLineSearchAppendOptionsPrefix(fas->linesearch_smooth, optionsprefix);CHKERRQ(ierr); 728f1c6b773SPeter Brune ierr = SNESLineSearchSetFromOptions(fas->linesearch_smooth);CHKERRQ(ierr); 72922c1e704SPeter Brune } 73022c1e704SPeter Brune 7316273346dSPeter Brune /* setup FAS work vectors */ 7326273346dSPeter Brune if (fas->galerkin) { 7336273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7346273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7356273346dSPeter Brune } 736421d9b32SPeter Brune PetscFunctionReturn(0); 737421d9b32SPeter Brune } 738421d9b32SPeter Brune 739421d9b32SPeter Brune #undef __FUNCT__ 740421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 741421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 742421d9b32SPeter Brune { 743ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 744ee78dd50SPeter Brune PetscInt levels = 1; 745fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 746421d9b32SPeter Brune PetscErrorCode ierr; 747ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 74807144faaSPeter Brune SNESFASType fastype; 749fde0ff24SPeter Brune const char *optionsprefix; 750fde0ff24SPeter Brune const char *prefix; 751f1c6b773SPeter Brune SNESLineSearch linesearch; 752421d9b32SPeter Brune 753421d9b32SPeter Brune PetscFunctionBegin; 754c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 755ee78dd50SPeter Brune 756ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 757ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 758ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 759c732cbdbSBarry Smith if (!flg && snes->dm) { 760c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 761c732cbdbSBarry Smith levels++; 762d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 763c732cbdbSBarry Smith } 764ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 765ee78dd50SPeter Brune } 76607144faaSPeter Brune fastype = fas->fastype; 76707144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 76807144faaSPeter Brune if (flg) { 76907144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 77007144faaSPeter Brune } 771ee78dd50SPeter Brune 772fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 773fde0ff24SPeter Brune 774fde0ff24SPeter Brune /* smoother setup options */ 775fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 776fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 777fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 778fde0ff24SPeter Brune if (fas->level == 0) { 779fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 780fde0ff24SPeter Brune } 781ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 782fde0ff24SPeter Brune 783d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 784ee78dd50SPeter Brune 7856273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 7866273346dSPeter Brune 787646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 788ee78dd50SPeter Brune 789ee78dd50SPeter Brune 790421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7918cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 792eff52c0eSPeter Brune 793fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 794fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 795fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 796fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 797fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 798fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 799fde0ff24SPeter Brune } 800fde0ff24SPeter Brune 801fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8028cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8038cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8048cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 805fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8068cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8078cc86e31SPeter Brune } 8088cc86e31SPeter Brune 809fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 81067339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 811ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 81267339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 813fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 814293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 815ee78dd50SPeter Brune } 816fde0ff24SPeter Brune 817fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 818fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 819fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 820fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 821fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 822fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 823fde0ff24SPeter Brune } 824fde0ff24SPeter Brune 825fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 826fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 827fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 828fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 829fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 830fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 831fde0ff24SPeter Brune } 832fde0ff24SPeter Brune 8331a266240SBarry Smith if (fas->upsmooth) { 834c60f73f4SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->stol, 1, 1000);CHKERRQ(ierr); 8351a266240SBarry Smith } 8361a266240SBarry Smith 8371a266240SBarry Smith if (fas->downsmooth) { 838c60f73f4SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->stol, 1, 1000);CHKERRQ(ierr); 8391a266240SBarry Smith } 840ee78dd50SPeter Brune 841742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 842c60f73f4SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, fas->n_cycles, 1000);CHKERRQ(ierr); 843742fe5e2SPeter Brune } 844742fe5e2SPeter Brune 845ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 846ce11c761SPeter Brune if (!fas->downsmooth) { 84793dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 8485d115551SPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 849ce11c761SPeter Brune if (fas->level == 0) { 85079d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 851ce11c761SPeter Brune } 852ce11c761SPeter Brune } 853ce11c761SPeter Brune 854ce11c761SPeter Brune if (!fas->upsmooth) { 85593dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 8565d115551SPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 857ce11c761SPeter Brune } 858ce11c761SPeter Brune 859ee78dd50SPeter Brune if (monflg) { 860646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 861794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8622f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 863742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 864293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 865293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 866d28543b3SPeter Brune } else { 867d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 868d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 869d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 870d28543b3SPeter Brune } 871ee78dd50SPeter Brune } 872ee78dd50SPeter Brune 8739e764e56SPeter Brune 8749e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 8759e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 8769e764e56SPeter Brune if (!snes->linesearch) { 877f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 878f1c6b773SPeter Brune ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr); 8799e764e56SPeter Brune } 8809e764e56SPeter Brune } 8819e764e56SPeter Brune 882ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 883d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 884421d9b32SPeter Brune PetscFunctionReturn(0); 885421d9b32SPeter Brune } 886421d9b32SPeter Brune 887421d9b32SPeter Brune #undef __FUNCT__ 888421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 889421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 890421d9b32SPeter Brune { 891421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 892421d9b32SPeter Brune PetscBool iascii; 893421d9b32SPeter Brune PetscErrorCode ierr; 894421d9b32SPeter Brune 895421d9b32SPeter Brune PetscFunctionBegin; 896421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 897421d9b32SPeter Brune if (iascii) { 8982e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 899421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 900bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 901ee78dd50SPeter Brune if (fas->upsmooth) { 90239a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 903421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 904ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 905421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 906421d9b32SPeter Brune } else { 907ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 908421d9b32SPeter Brune } 909ee78dd50SPeter Brune if (fas->downsmooth) { 91039a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 911421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 912ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 913421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 914421d9b32SPeter Brune } else { 915ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 916421d9b32SPeter Brune } 917421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 918421d9b32SPeter Brune } else { 919421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 920421d9b32SPeter Brune } 921421d9b32SPeter Brune PetscFunctionReturn(0); 922421d9b32SPeter Brune } 923421d9b32SPeter Brune 924421d9b32SPeter Brune #undef __FUNCT__ 92539bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92639bd7f45SPeter Brune /* 92739bd7f45SPeter Brune Defines the action of the downsmoother 92839bd7f45SPeter Brune */ 92939bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 93039bd7f45SPeter Brune PetscErrorCode ierr = 0; 93122c1e704SPeter Brune PetscReal fnorm; 932742fe5e2SPeter Brune SNESConvergedReason reason; 93339bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 93422c1e704SPeter Brune Vec Y, FPC; 935fde0ff24SPeter Brune PetscBool lssuccess; 936fde0ff24SPeter Brune PetscInt k; 937421d9b32SPeter Brune PetscFunctionBegin; 938fde0ff24SPeter Brune Y = snes->work[3]; 939d1adcc6fSPeter Brune if (fas->downsmooth) { 940d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 941742fe5e2SPeter Brune /* check convergence reason for the smoother */ 942742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 943742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 944742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 945742fe5e2SPeter Brune PetscFunctionReturn(0); 946742fe5e2SPeter Brune } 9474b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9484b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 949fde0ff24SPeter Brune } else { 950fde0ff24SPeter Brune /* basic richardson smoother */ 951fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 952794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 953794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 954fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 955f1c6b773SPeter Brune ierr = SNESLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 956f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 957fde0ff24SPeter Brune if (!lssuccess) { 958fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 959fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9602f7ea302SPeter Brune PetscFunctionReturn(0); 9612f7ea302SPeter Brune } 962fe6f9142SPeter Brune } 963fde0ff24SPeter Brune } 964fde0ff24SPeter Brune } 96539bd7f45SPeter Brune PetscFunctionReturn(0); 96639bd7f45SPeter Brune } 96739bd7f45SPeter Brune 96839bd7f45SPeter Brune 96939bd7f45SPeter Brune #undef __FUNCT__ 97039bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 97139bd7f45SPeter Brune /* 97207144faaSPeter Brune Defines the action of the upsmoother 97339bd7f45SPeter Brune */ 97439bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 97539bd7f45SPeter Brune PetscErrorCode ierr = 0; 97622c1e704SPeter Brune PetscReal fnorm; 97739bd7f45SPeter Brune SNESConvergedReason reason; 97839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 97922c1e704SPeter Brune Vec Y, FPC; 980fde0ff24SPeter Brune PetscBool lssuccess; 981fde0ff24SPeter Brune PetscInt k; 98239bd7f45SPeter Brune PetscFunctionBegin; 983fde0ff24SPeter Brune Y = snes->work[3]; 98439bd7f45SPeter Brune if (fas->upsmooth) { 985fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 98639bd7f45SPeter Brune /* check convergence reason for the smoother */ 987fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 98839bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 98939bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 99039bd7f45SPeter Brune PetscFunctionReturn(0); 99139bd7f45SPeter Brune } 9924b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9934b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 994fde0ff24SPeter Brune } else { 995fde0ff24SPeter Brune /* basic richardson smoother */ 996fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 99739bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 99839bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 999fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1000f1c6b773SPeter Brune ierr = SNESLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr); 1001f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr); 1002fde0ff24SPeter Brune if (!lssuccess) { 1003fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1004fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 100539bd7f45SPeter Brune PetscFunctionReturn(0); 100639bd7f45SPeter Brune } 100739bd7f45SPeter Brune } 1008fde0ff24SPeter Brune } 1009fde0ff24SPeter Brune } 101039bd7f45SPeter Brune PetscFunctionReturn(0); 101139bd7f45SPeter Brune } 101239bd7f45SPeter Brune 101339bd7f45SPeter Brune #undef __FUNCT__ 1014938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1015938e4a01SJed Brown /*@ 1016938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1017938e4a01SJed Brown 1018938e4a01SJed Brown Collective 1019938e4a01SJed Brown 1020938e4a01SJed Brown Input Arguments: 1021938e4a01SJed Brown . snes - SNESFAS 1022938e4a01SJed Brown 1023938e4a01SJed Brown Output Arguments: 1024938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1025938e4a01SJed Brown 1026938e4a01SJed Brown Level: developer 1027938e4a01SJed Brown 1028938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1029938e4a01SJed Brown @*/ 1030938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1031938e4a01SJed Brown { 1032938e4a01SJed Brown PetscErrorCode ierr; 1033938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1034938e4a01SJed Brown 1035938e4a01SJed Brown PetscFunctionBegin; 1036938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1037938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1038938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1039938e4a01SJed Brown PetscFunctionReturn(0); 1040938e4a01SJed Brown } 1041938e4a01SJed Brown 1042e9923e8dSJed Brown #undef __FUNCT__ 1043e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1044e9923e8dSJed Brown /*@ 1045e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1046e9923e8dSJed Brown 1047e9923e8dSJed Brown Collective 1048e9923e8dSJed Brown 1049e9923e8dSJed Brown Input Arguments: 1050e9923e8dSJed Brown + fine - SNES from which to restrict 1051e9923e8dSJed Brown - Xfine - vector to restrict 1052e9923e8dSJed Brown 1053e9923e8dSJed Brown Output Arguments: 1054e9923e8dSJed Brown . Xcoarse - result of restriction 1055e9923e8dSJed Brown 1056e9923e8dSJed Brown Level: developer 1057e9923e8dSJed Brown 1058e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1059e9923e8dSJed Brown @*/ 1060e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1061e9923e8dSJed Brown { 1062e9923e8dSJed Brown PetscErrorCode ierr; 1063e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1064e9923e8dSJed Brown 1065e9923e8dSJed Brown PetscFunctionBegin; 1066e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1067e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1068e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1069e9923e8dSJed Brown if (fas->inject) { 1070e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1071e9923e8dSJed Brown } else { 1072e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1073e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1074e9923e8dSJed Brown } 1075e9923e8dSJed Brown PetscFunctionReturn(0); 1076e9923e8dSJed Brown } 1077e9923e8dSJed Brown 1078e9923e8dSJed Brown #undef __FUNCT__ 107939bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 108039bd7f45SPeter Brune /* 108139bd7f45SPeter Brune 108239bd7f45SPeter Brune Performs the FAS coarse correction as: 108339bd7f45SPeter Brune 108439bd7f45SPeter Brune fine problem: F(x) = 0 108539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 108639bd7f45SPeter Brune 108739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 108839bd7f45SPeter Brune 108939bd7f45SPeter Brune with correction: 109039bd7f45SPeter Brune 109139bd7f45SPeter Brune 109239bd7f45SPeter Brune 109339bd7f45SPeter Brune */ 109439a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 109539bd7f45SPeter Brune PetscErrorCode ierr; 109639bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 109739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 109839bd7f45SPeter Brune SNESConvergedReason reason; 109939bd7f45SPeter Brune PetscFunctionBegin; 1100fa9694d7SPeter Brune if (fas->next) { 1101c90fad12SPeter Brune X_c = fas->next->vec_sol; 1102293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1103c90fad12SPeter Brune F_c = fas->next->vec_func; 1104742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1105efe1f98aSPeter Brune 1106938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1107293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1108293a7e31SPeter Brune 1109293a7e31SPeter Brune /* restrict the defect */ 1110293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1111293a7e31SPeter Brune 1112c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1113ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1114c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1115742fe5e2SPeter Brune 1116293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1117c90fad12SPeter Brune 1118ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1119ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1120c90fad12SPeter Brune 1121c90fad12SPeter Brune /* recurse to the next level */ 1122f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1123742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1124742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1125742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1126742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1127742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1128742fe5e2SPeter Brune PetscFunctionReturn(0); 1129742fe5e2SPeter Brune } 1130742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1131ee78dd50SPeter Brune 1132fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1133fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 113439bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1135293a7e31SPeter Brune } 113639bd7f45SPeter Brune PetscFunctionReturn(0); 113739bd7f45SPeter Brune } 113839bd7f45SPeter Brune 113939bd7f45SPeter Brune #undef __FUNCT__ 114039bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 114139bd7f45SPeter Brune /* 114239bd7f45SPeter Brune 114339bd7f45SPeter Brune The additive cycle looks like: 114439bd7f45SPeter Brune 114507144faaSPeter Brune xhat = x 114607144faaSPeter Brune xhat = dS(x, b) 114707144faaSPeter Brune x = coarsecorrection(xhat, b_d) 114807144faaSPeter Brune x = x + nu*(xhat - x); 114939bd7f45SPeter Brune (optional) x = uS(x, b) 115039bd7f45SPeter Brune 115139bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 115239bd7f45SPeter Brune 115339bd7f45SPeter Brune */ 115439bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 115507144faaSPeter Brune Vec F, B, Xhat; 115622c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 115739bd7f45SPeter Brune PetscErrorCode ierr; 115807144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 115907144faaSPeter Brune SNESConvergedReason reason; 116022c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 116122c1e704SPeter Brune PetscBool lssuccess; 116239bd7f45SPeter Brune PetscFunctionBegin; 116339bd7f45SPeter Brune 116439bd7f45SPeter Brune F = snes->vec_func; 116539bd7f45SPeter Brune B = snes->vec_rhs; 1166fde0ff24SPeter Brune Xhat = snes->work[3]; 116707144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 116807144faaSPeter Brune /* recurse first */ 116907144faaSPeter Brune if (fas->next) { 117007144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 117139bd7f45SPeter Brune 117207144faaSPeter Brune X_c = fas->next->vec_sol; 117307144faaSPeter Brune Xo_c = fas->next->work[0]; 117407144faaSPeter Brune F_c = fas->next->vec_func; 117507144faaSPeter Brune B_c = fas->next->vec_rhs; 117639bd7f45SPeter Brune 1177938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 117807144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 117907144faaSPeter Brune 118007144faaSPeter Brune /* restrict the defect */ 118107144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 118207144faaSPeter Brune 118307144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 118407144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 118507144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 118607144faaSPeter Brune 118707144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 118807144faaSPeter Brune 118907144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 119007144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 119107144faaSPeter Brune 119207144faaSPeter Brune /* recurse */ 119307144faaSPeter Brune fas->next->vec_rhs = B_c; 119407144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 119507144faaSPeter Brune 119607144faaSPeter Brune /* smooth on this level */ 119707144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 119807144faaSPeter Brune 119907144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 120007144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 120107144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 120207144faaSPeter Brune PetscFunctionReturn(0); 120307144faaSPeter Brune } 120407144faaSPeter Brune 120507144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1206c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1207ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 120807144faaSPeter Brune 1209ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1210ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1211ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1212f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 1213f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 12149e764e56SPeter Brune if (!lssuccess) { 12159e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 12169e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 12179e764e56SPeter Brune PetscFunctionReturn(0); 12189e764e56SPeter Brune } 12199e764e56SPeter Brune } 1220f1c6b773SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); 122107144faaSPeter Brune } else { 122207144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 122307144faaSPeter Brune } 122439bd7f45SPeter Brune PetscFunctionReturn(0); 122539bd7f45SPeter Brune } 122639bd7f45SPeter Brune 122739bd7f45SPeter Brune #undef __FUNCT__ 122839bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 122939bd7f45SPeter Brune /* 123039bd7f45SPeter Brune 123139bd7f45SPeter Brune Defines the FAS cycle as: 123239bd7f45SPeter Brune 123339bd7f45SPeter Brune fine problem: F(x) = 0 123439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 123539bd7f45SPeter Brune 123639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 123739bd7f45SPeter Brune 123839bd7f45SPeter Brune correction: 123939bd7f45SPeter Brune 124039bd7f45SPeter Brune x = x + I(x^c - Rx) 124139bd7f45SPeter Brune 124239bd7f45SPeter Brune */ 124339bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 124439bd7f45SPeter Brune 124539bd7f45SPeter Brune PetscErrorCode ierr; 124639bd7f45SPeter Brune Vec F,B; 124739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 124839bd7f45SPeter Brune 124939bd7f45SPeter Brune PetscFunctionBegin; 125039bd7f45SPeter Brune F = snes->vec_func; 125139bd7f45SPeter Brune B = snes->vec_rhs; 125239bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 125339bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 125439bd7f45SPeter Brune 125539a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 125639bd7f45SPeter Brune 1257c90fad12SPeter Brune if (fas->level != 0) { 125839bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1259fe6f9142SPeter Brune } 1260fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1261fa9694d7SPeter Brune 1262fa9694d7SPeter Brune PetscFunctionReturn(0); 1263421d9b32SPeter Brune } 1264421d9b32SPeter Brune 1265421d9b32SPeter Brune #undef __FUNCT__ 1266421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1267421d9b32SPeter Brune 1268421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1269421d9b32SPeter Brune { 1270fa9694d7SPeter Brune PetscErrorCode ierr; 1271fe6f9142SPeter Brune PetscInt i, maxits; 1272ddb5aff1SPeter Brune Vec X, F; 1273fe6f9142SPeter Brune PetscReal fnorm; 1274b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1275b17ce1afSJed Brown DM dm; 1276b17ce1afSJed Brown 1277421d9b32SPeter Brune PetscFunctionBegin; 1278fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1279fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1280fa9694d7SPeter Brune X = snes->vec_sol; 1281f5a6d4f9SBarry Smith F = snes->vec_func; 1282293a7e31SPeter Brune 1283293a7e31SPeter Brune /*norm setup */ 1284fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1285fe6f9142SPeter Brune snes->iter = 0; 1286fe6f9142SPeter Brune snes->norm = 0.; 1287fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1288fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1289fe6f9142SPeter Brune if (snes->domainerror) { 1290fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1291fe6f9142SPeter Brune PetscFunctionReturn(0); 1292fe6f9142SPeter Brune } 1293fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1294fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1295fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1296fe6f9142SPeter Brune snes->norm = fnorm; 1297fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1298fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1299fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1300fe6f9142SPeter Brune 1301fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1302fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1303fe6f9142SPeter Brune /* test convergence */ 1304fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1305fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1306b17ce1afSJed Brown 1307b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1308b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1309b17ce1afSJed Brown DM dmcoarse; 1310b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1311b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1312b17ce1afSJed Brown dm = dmcoarse; 1313b17ce1afSJed Brown } 1314b17ce1afSJed Brown 1315fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1316fe6f9142SPeter Brune /* Call general purpose update function */ 1317646217ecSPeter Brune 1318fe6f9142SPeter Brune if (snes->ops->update) { 1319fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1320fe6f9142SPeter Brune } 132107144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 132239bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 132307144faaSPeter Brune } else { 132407144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 132507144faaSPeter Brune } 1326742fe5e2SPeter Brune 1327742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1328742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1329742fe5e2SPeter Brune PetscFunctionReturn(0); 1330742fe5e2SPeter Brune } 1331c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1332c90fad12SPeter Brune /* Monitor convergence */ 1333c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1334c90fad12SPeter Brune snes->iter = i+1; 1335c90fad12SPeter Brune snes->norm = fnorm; 1336c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1337c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1338c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1339c90fad12SPeter Brune /* Test for convergence */ 1340c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1341c90fad12SPeter Brune if (snes->reason) break; 1342fe6f9142SPeter Brune } 1343fe6f9142SPeter Brune if (i == maxits) { 1344fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1345fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1346fe6f9142SPeter Brune } 1347421d9b32SPeter Brune PetscFunctionReturn(0); 1348421d9b32SPeter Brune } 1349