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 24ddebd997SPeter Brune #undef __FUNCT__ 25ddebd997SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_FAS" 26ddebd997SPeter Brune PetscErrorCode SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type) 27ddebd997SPeter Brune { 28ddebd997SPeter Brune PetscErrorCode ierr; 29ddebd997SPeter Brune PetscFunctionBegin; 30ddebd997SPeter Brune 31ddebd997SPeter Brune switch (type) { 32ddebd997SPeter Brune case SNES_LS_BASIC: 33ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 34ddebd997SPeter Brune break; 35ddebd997SPeter Brune case SNES_LS_BASIC_NONORMS: 36ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 37ddebd997SPeter Brune break; 38ddebd997SPeter Brune case SNES_LS_QUADRATIC: 39ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 40ddebd997SPeter Brune break; 41ddebd997SPeter Brune default: 42ddebd997SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 43ddebd997SPeter Brune break; 44ddebd997SPeter Brune } 45ddebd997SPeter Brune snes->ls_type = type; 46ddebd997SPeter Brune PetscFunctionReturn(0); 47ddebd997SPeter Brune } 48ddebd997SPeter Brune EXTERN_C_END 49ddebd997SPeter Brune 50ddebd997SPeter Brune EXTERN_C_BEGIN 51421d9b32SPeter Brune 52421d9b32SPeter Brune #undef __FUNCT__ 53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 54421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 55421d9b32SPeter Brune { 56421d9b32SPeter Brune SNES_FAS * fas; 57421d9b32SPeter Brune PetscErrorCode ierr; 58421d9b32SPeter Brune 59421d9b32SPeter Brune PetscFunctionBegin; 60421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 61421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 62421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 63421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 64421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 65421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 66421d9b32SPeter Brune 67ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 68ed020824SBarry Smith snes->usespc = PETSC_FALSE; 69ed020824SBarry Smith 70421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 71421d9b32SPeter Brune snes->data = (void*) fas; 72421d9b32SPeter Brune fas->level = 0; 73293a7e31SPeter Brune fas->levels = 1; 74ee78dd50SPeter Brune fas->n_cycles = 1; 75ee78dd50SPeter Brune fas->max_up_it = 1; 76ee78dd50SPeter Brune fas->max_down_it = 1; 77ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 78ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 79421d9b32SPeter Brune fas->next = PETSC_NULL; 806273346dSPeter Brune fas->previous = PETSC_NULL; 81421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 82421d9b32SPeter Brune fas->restrct = PETSC_NULL; 83efe1f98aSPeter Brune fas->inject = PETSC_NULL; 84cc05f883SPeter Brune fas->monitor = PETSC_NULL; 85cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 86ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 87ddebd997SPeter Brune 88ddebd997SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 89ddebd997SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 90ddebd997SPeter Brune 91421d9b32SPeter Brune PetscFunctionReturn(0); 92421d9b32SPeter Brune } 93421d9b32SPeter Brune EXTERN_C_END 94421d9b32SPeter Brune 95421d9b32SPeter Brune #undef __FUNCT__ 96421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 97c79ef259SPeter Brune /*@ 98722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 99c79ef259SPeter Brune 100c79ef259SPeter Brune Input Parameter: 101c79ef259SPeter Brune . snes - the preconditioner context 102c79ef259SPeter Brune 103c79ef259SPeter Brune Output parameter: 104c79ef259SPeter Brune . levels - the number of levels 105c79ef259SPeter Brune 106c79ef259SPeter Brune Level: advanced 107c79ef259SPeter Brune 108c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 109c79ef259SPeter Brune 110c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 111c79ef259SPeter Brune @*/ 112421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 113421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 114421d9b32SPeter Brune PetscFunctionBegin; 115ee1fd11aSPeter Brune *levels = fas->levels; 116421d9b32SPeter Brune PetscFunctionReturn(0); 117421d9b32SPeter Brune } 118421d9b32SPeter Brune 119421d9b32SPeter Brune #undef __FUNCT__ 120646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 121c79ef259SPeter Brune /*@ 122c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 123c79ef259SPeter Brune complicated cycling. 124c79ef259SPeter Brune 125c79ef259SPeter Brune Logically Collective on SNES 126c79ef259SPeter Brune 127c79ef259SPeter Brune Input Parameters: 128c79ef259SPeter Brune + snes - the multigrid context 129c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 130c79ef259SPeter Brune 131c79ef259SPeter Brune Options Database Key: 132c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 133c79ef259SPeter Brune 134c79ef259SPeter Brune Level: advanced 135c79ef259SPeter Brune 136c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 137c79ef259SPeter Brune 138c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 139c79ef259SPeter Brune @*/ 140646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 141646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 142646217ecSPeter Brune PetscErrorCode ierr; 143646217ecSPeter Brune PetscFunctionBegin; 144646217ecSPeter Brune fas->n_cycles = cycles; 145646217ecSPeter Brune if (fas->next) { 146646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 147646217ecSPeter Brune } 148646217ecSPeter Brune PetscFunctionReturn(0); 149646217ecSPeter Brune } 150646217ecSPeter Brune 151eff52c0eSPeter Brune #undef __FUNCT__ 152c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 153c79ef259SPeter Brune /*@ 154722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 155c79ef259SPeter Brune 156c79ef259SPeter Brune Logically Collective on SNES 157c79ef259SPeter Brune 158c79ef259SPeter Brune Input Parameters: 159c79ef259SPeter Brune + snes - the multigrid context 160c79ef259SPeter Brune . level - the level to set the number of cycles on 161c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 162c79ef259SPeter Brune 163c79ef259SPeter Brune Level: advanced 164c79ef259SPeter Brune 165c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 166c79ef259SPeter Brune 167c79ef259SPeter Brune .seealso: SNESFASSetCycles() 168c79ef259SPeter Brune @*/ 169c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 170c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 171c79ef259SPeter Brune PetscInt top_level = fas->level,i; 172c79ef259SPeter Brune 173c79ef259SPeter Brune PetscFunctionBegin; 174c79ef259SPeter Brune if (level > top_level) 175c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 176c79ef259SPeter Brune /* get to the correct level */ 177c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 178c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 179c79ef259SPeter Brune } 180c79ef259SPeter Brune if (fas->level != level) 181c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 182c79ef259SPeter Brune fas->n_cycles = cycles; 183c79ef259SPeter Brune PetscFunctionReturn(0); 184c79ef259SPeter Brune } 185c79ef259SPeter Brune 186c79ef259SPeter Brune #undef __FUNCT__ 187eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 188aeed3662SMatthew G Knepley /*@C 189c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 190c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 191c79ef259SPeter Brune and nonlinear preconditioners. 192c79ef259SPeter Brune 193c79ef259SPeter Brune Logically Collective on SNES 194c79ef259SPeter Brune 195c79ef259SPeter Brune Input Parameters: 196c79ef259SPeter Brune + snes - the multigrid context 197c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 198c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 199c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 200c79ef259SPeter Brune 201c79ef259SPeter Brune Level: advanced 202c79ef259SPeter Brune 203c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 204c79ef259SPeter Brune 205c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 206c79ef259SPeter Brune @*/ 207eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 208eff52c0eSPeter Brune PetscErrorCode ierr = 0; 209eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 210eff52c0eSPeter Brune PetscFunctionBegin; 211eff52c0eSPeter Brune 212eff52c0eSPeter Brune if (gsfunc) { 213eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 214eff52c0eSPeter Brune /* push the provided GS up the tree */ 215eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 216eff52c0eSPeter Brune } else if (snes->ops->computegs) { 217eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 218eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 219eff52c0eSPeter Brune } else if (use_gs) { 220eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 221eff52c0eSPeter Brune } 222eff52c0eSPeter Brune PetscFunctionReturn(0); 223eff52c0eSPeter Brune } 224eff52c0eSPeter Brune 225eff52c0eSPeter Brune #undef __FUNCT__ 226eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 227aeed3662SMatthew G Knepley /*@C 228c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 229c79ef259SPeter Brune 230c79ef259SPeter Brune Logically Collective on SNES 231c79ef259SPeter Brune 232c79ef259SPeter Brune Input Parameters: 233c79ef259SPeter Brune + snes - the multigrid context 234c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 235c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 236c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 237c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 238c79ef259SPeter Brune 239c79ef259SPeter Brune Level: advanced 240c79ef259SPeter Brune 241c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 242c79ef259SPeter Brune 243c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 244c79ef259SPeter Brune @*/ 245eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 246eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 247eff52c0eSPeter Brune PetscErrorCode ierr; 248eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 249eff52c0eSPeter Brune SNES cur_snes = snes; 250eff52c0eSPeter Brune PetscFunctionBegin; 251eff52c0eSPeter Brune if (level > top_level) 252eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 253eff52c0eSPeter Brune /* get to the correct level */ 254eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 255eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 256eff52c0eSPeter Brune cur_snes = fas->next; 257eff52c0eSPeter Brune } 258eff52c0eSPeter Brune if (fas->level != level) 259eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 260eff52c0eSPeter Brune if (gsfunc) { 2616273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 262eff52c0eSPeter Brune } 263eff52c0eSPeter Brune PetscFunctionReturn(0); 264eff52c0eSPeter Brune } 265eff52c0eSPeter Brune 266646217ecSPeter Brune #undef __FUNCT__ 267421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 268c79ef259SPeter Brune /*@ 269c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 270c79ef259SPeter Brune level of the FAS hierarchy. 271c79ef259SPeter Brune 272c79ef259SPeter Brune Input Parameters: 273c79ef259SPeter Brune + snes - the multigrid context 274c79ef259SPeter Brune level - the level to get 275c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 276c79ef259SPeter Brune 277c79ef259SPeter Brune Level: advanced 278c79ef259SPeter Brune 279c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 280c79ef259SPeter Brune 281c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 282c79ef259SPeter Brune @*/ 283421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 284421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 285421d9b32SPeter Brune PetscInt levels = fas->level; 286421d9b32SPeter Brune PetscInt i; 287421d9b32SPeter Brune PetscFunctionBegin; 288421d9b32SPeter Brune *lsnes = snes; 289421d9b32SPeter Brune if (fas->level < level) { 290421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 291421d9b32SPeter Brune } 292421d9b32SPeter Brune if (level > levels - 1) { 293421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 294421d9b32SPeter Brune } 295421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 296421d9b32SPeter Brune *lsnes = fas->next; 297421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 298421d9b32SPeter Brune } 299421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 300421d9b32SPeter Brune PetscFunctionReturn(0); 301421d9b32SPeter Brune } 302421d9b32SPeter Brune 303421d9b32SPeter Brune #undef __FUNCT__ 30407144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 30507144faaSPeter Brune /*@ 30607144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 30707144faaSPeter Brune e 30807144faaSPeter Brune 30907144faaSPeter Brune 31007144faaSPeter Brune @*/ 31107144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 31207144faaSPeter Brune { 31307144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31407144faaSPeter Brune 31507144faaSPeter Brune PetscFunctionBegin; 31607144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31707144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 31807144faaSPeter Brune fas->fastype = fastype; 31907144faaSPeter Brune PetscFunctionReturn(0); 32007144faaSPeter Brune } 32107144faaSPeter Brune 32207144faaSPeter Brune 32307144faaSPeter Brune 32407144faaSPeter Brune #undef __FUNCT__ 325421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 326aeed3662SMatthew G Knepley /*@C 327c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 328c79ef259SPeter Brune Must be called before any other FAS routine. 329c79ef259SPeter Brune 330c79ef259SPeter Brune Input Parameters: 331c79ef259SPeter Brune + snes - the snes context 332c79ef259SPeter Brune . levels - the number of levels 333c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 334c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 335c79ef259SPeter Brune Fortran. 336c79ef259SPeter Brune 337c79ef259SPeter Brune Level: intermediate 338c79ef259SPeter Brune 339c79ef259SPeter Brune Notes: 340c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 341c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 342c79ef259SPeter Brune 343c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 344c79ef259SPeter Brune 345c79ef259SPeter Brune .seealso: SNESFASGetLevels() 346c79ef259SPeter Brune @*/ 347421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 348421d9b32SPeter Brune PetscErrorCode ierr; 349421d9b32SPeter Brune PetscInt i; 350421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3516273346dSPeter Brune SNES prevsnes; 352421d9b32SPeter Brune MPI_Comm comm; 353421d9b32SPeter Brune PetscFunctionBegin; 354ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 355ee1fd11aSPeter Brune if (levels == fas->levels) { 356ee1fd11aSPeter Brune if (!comms) 357ee1fd11aSPeter Brune PetscFunctionReturn(0); 358ee1fd11aSPeter Brune } 359421d9b32SPeter Brune /* user has changed the number of levels; reset */ 360421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 361421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 362421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 363ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3646273346dSPeter Brune fas->previous = PETSC_NULL; 3656273346dSPeter Brune prevsnes = snes; 366421d9b32SPeter Brune /* setup the finest level */ 367421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 368421d9b32SPeter Brune if (comms) comm = comms[i]; 369421d9b32SPeter Brune fas->level = i; 370421d9b32SPeter Brune fas->levels = levels; 371ee1fd11aSPeter Brune fas->next = PETSC_NULL; 372421d9b32SPeter Brune if (i > 0) { 373421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3744a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 375421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 376794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3776273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3786273346dSPeter Brune prevsnes = fas->next; 3796273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 380421d9b32SPeter Brune } 381421d9b32SPeter Brune } 382421d9b32SPeter Brune PetscFunctionReturn(0); 383421d9b32SPeter Brune } 384421d9b32SPeter Brune 385421d9b32SPeter Brune #undef __FUNCT__ 386c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 387c79ef259SPeter Brune /*@ 388c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-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_smoothup <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: SNESFASSetNumberSmoothDown() 405c79ef259SPeter Brune @*/ 406c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(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_up_it = n; 411c79ef259SPeter Brune if (fas->next) { 412c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 413c79ef259SPeter Brune } 414c79ef259SPeter Brune PetscFunctionReturn(0); 415c79ef259SPeter Brune } 416c79ef259SPeter Brune 417c79ef259SPeter Brune #undef __FUNCT__ 418c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 419c79ef259SPeter Brune /*@ 420c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 421c79ef259SPeter Brune use on all levels. 422c79ef259SPeter Brune 423fde0ff24SPeter Brune Logically Collective on SNES 424c79ef259SPeter Brune 425c79ef259SPeter Brune Input Parameters: 426c79ef259SPeter Brune + snes - the multigrid context 427c79ef259SPeter Brune - n - the number of smoothing steps 428c79ef259SPeter Brune 429c79ef259SPeter Brune Options Database Key: 430d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 431c79ef259SPeter Brune 432c79ef259SPeter Brune Level: advanced 433c79ef259SPeter Brune 434c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 435c79ef259SPeter Brune 436c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 437c79ef259SPeter Brune @*/ 438c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 439c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 440c79ef259SPeter Brune PetscErrorCode ierr = 0; 441c79ef259SPeter Brune PetscFunctionBegin; 442c79ef259SPeter Brune fas->max_down_it = n; 443c79ef259SPeter Brune if (fas->next) { 444c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 445c79ef259SPeter Brune } 446c79ef259SPeter Brune PetscFunctionReturn(0); 447c79ef259SPeter Brune } 448c79ef259SPeter Brune 449c79ef259SPeter Brune #undef __FUNCT__ 450421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 451c79ef259SPeter Brune /*@ 452c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 453c79ef259SPeter Brune interpolation from l-1 to the lth level 454c79ef259SPeter Brune 455c79ef259SPeter Brune Input Parameters: 456c79ef259SPeter Brune + snes - the multigrid context 457c79ef259SPeter Brune . mat - the interpolation operator 458c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 459c79ef259SPeter Brune 460c79ef259SPeter Brune Level: advanced 461c79ef259SPeter Brune 462c79ef259SPeter Brune Notes: 463c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 464c79ef259SPeter Brune for the same level. 465c79ef259SPeter Brune 466c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 467c79ef259SPeter Brune out from the matrix size which one it is. 468c79ef259SPeter Brune 469c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 470c79ef259SPeter Brune 471*bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 472c79ef259SPeter Brune @*/ 473421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 474421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 475d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 476*bccf9bb3SJed Brown PetscErrorCode ierr; 477421d9b32SPeter Brune 478421d9b32SPeter Brune PetscFunctionBegin; 479421d9b32SPeter Brune if (level > top_level) 480421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 481421d9b32SPeter Brune /* get to the correct level */ 482d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 483421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 484421d9b32SPeter Brune } 485421d9b32SPeter Brune if (fas->level != level) 486421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 487*bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 488421d9b32SPeter Brune fas->interpolate = mat; 489421d9b32SPeter Brune PetscFunctionReturn(0); 490421d9b32SPeter Brune } 491421d9b32SPeter Brune 492421d9b32SPeter Brune #undef __FUNCT__ 493421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 494c79ef259SPeter Brune /*@ 495c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 496c79ef259SPeter Brune from level l to l-1. 497c79ef259SPeter Brune 498c79ef259SPeter Brune Input Parameters: 499c79ef259SPeter Brune + snes - the multigrid context 500c79ef259SPeter Brune . mat - the restriction matrix 501c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 502c79ef259SPeter Brune 503c79ef259SPeter Brune Level: advanced 504c79ef259SPeter Brune 505c79ef259SPeter Brune Notes: 506c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 507c79ef259SPeter Brune for the same level. 508c79ef259SPeter Brune 509c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 510c79ef259SPeter Brune out from the matrix size which one it is. 511c79ef259SPeter Brune 512fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 513c79ef259SPeter Brune is used. 514c79ef259SPeter Brune 515c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 516c79ef259SPeter Brune 517c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 518c79ef259SPeter Brune @*/ 519421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 520421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 521d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 522*bccf9bb3SJed Brown PetscErrorCode ierr; 523421d9b32SPeter Brune 524421d9b32SPeter Brune PetscFunctionBegin; 525421d9b32SPeter Brune if (level > top_level) 526421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 527421d9b32SPeter Brune /* get to the correct level */ 528d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 529421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 530421d9b32SPeter Brune } 531421d9b32SPeter Brune if (fas->level != level) 532421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 533*bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 534421d9b32SPeter Brune fas->restrct = mat; 535421d9b32SPeter Brune PetscFunctionReturn(0); 536421d9b32SPeter Brune } 537421d9b32SPeter Brune 538efe1f98aSPeter Brune 539421d9b32SPeter Brune #undef __FUNCT__ 540421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 541c79ef259SPeter Brune /*@ 542c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 543c79ef259SPeter Brune operator from level l to l-1. 544c79ef259SPeter Brune 545c79ef259SPeter Brune Input Parameters: 546c79ef259SPeter Brune + snes - the multigrid context 547c79ef259SPeter Brune . rscale - the restriction scaling 548c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 549c79ef259SPeter Brune 550c79ef259SPeter Brune Level: advanced 551c79ef259SPeter Brune 552c79ef259SPeter Brune Notes: 553c79ef259SPeter Brune This is only used in the case that the injection is not set. 554c79ef259SPeter Brune 555c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 556c79ef259SPeter Brune 557c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 558c79ef259SPeter Brune @*/ 559421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 560421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 561d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 562*bccf9bb3SJed Brown PetscErrorCode ierr; 563421d9b32SPeter Brune 564421d9b32SPeter Brune PetscFunctionBegin; 565421d9b32SPeter Brune if (level > top_level) 566421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 567421d9b32SPeter Brune /* get to the correct level */ 568d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 569421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 570421d9b32SPeter Brune } 571421d9b32SPeter Brune if (fas->level != level) 572421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 573*bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 574421d9b32SPeter Brune fas->rscale = rscale; 575421d9b32SPeter Brune PetscFunctionReturn(0); 576421d9b32SPeter Brune } 577421d9b32SPeter Brune 578efe1f98aSPeter Brune 579efe1f98aSPeter Brune #undef __FUNCT__ 580efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 581c79ef259SPeter Brune /*@ 582c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 583c79ef259SPeter Brune from level l to l-1. 584c79ef259SPeter Brune 585c79ef259SPeter Brune Input Parameters: 586c79ef259SPeter Brune + snes - the multigrid context 587c79ef259SPeter Brune . mat - the restriction matrix 588c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 589c79ef259SPeter Brune 590c79ef259SPeter Brune Level: advanced 591c79ef259SPeter Brune 592c79ef259SPeter Brune Notes: 593c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 594c79ef259SPeter Brune project the solution instead. 595c79ef259SPeter Brune 596c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 597c79ef259SPeter Brune 598c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 599c79ef259SPeter Brune @*/ 600efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 601efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 602efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 603*bccf9bb3SJed Brown PetscErrorCode ierr; 604efe1f98aSPeter Brune 605efe1f98aSPeter Brune PetscFunctionBegin; 606efe1f98aSPeter Brune if (level > top_level) 607efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 608efe1f98aSPeter Brune /* get to the correct level */ 609efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 610efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 611efe1f98aSPeter Brune } 612efe1f98aSPeter Brune if (fas->level != level) 613efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 614*bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 615efe1f98aSPeter Brune fas->inject = mat; 616efe1f98aSPeter Brune PetscFunctionReturn(0); 617efe1f98aSPeter Brune } 618efe1f98aSPeter Brune 619421d9b32SPeter Brune #undef __FUNCT__ 620421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 621421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 622421d9b32SPeter Brune { 62377df8cc4SPeter Brune PetscErrorCode ierr = 0; 624421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 625421d9b32SPeter Brune 626421d9b32SPeter Brune PetscFunctionBegin; 627*bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 628*bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 6293dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 630*bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 631*bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 632*bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6333dccd265SPeter Brune 6343dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6353dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 636742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 637421d9b32SPeter Brune PetscFunctionReturn(0); 638421d9b32SPeter Brune } 639421d9b32SPeter Brune 640421d9b32SPeter Brune #undef __FUNCT__ 641421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 642421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 643421d9b32SPeter Brune { 644421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 645742fe5e2SPeter Brune PetscErrorCode ierr = 0; 646421d9b32SPeter Brune 647421d9b32SPeter Brune PetscFunctionBegin; 648421d9b32SPeter Brune /* recursively resets and then destroys */ 64979d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 650421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 651421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 652ee1fd11aSPeter Brune 653421d9b32SPeter Brune PetscFunctionReturn(0); 654421d9b32SPeter Brune } 655421d9b32SPeter Brune 656421d9b32SPeter Brune #undef __FUNCT__ 657421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 658421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 659421d9b32SPeter Brune { 66048bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 661421d9b32SPeter Brune PetscErrorCode ierr; 662efe1f98aSPeter Brune VecScatter injscatter; 663d1adcc6fSPeter Brune PetscInt dm_levels; 6643dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 665d1adcc6fSPeter Brune 666421d9b32SPeter Brune PetscFunctionBegin; 667eff52c0eSPeter Brune 668cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 669d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 670d1adcc6fSPeter Brune dm_levels++; 671cc05f883SPeter Brune if (dm_levels > fas->levels) { 6723dccd265SPeter Brune 6733dccd265SPeter Brune /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/ 6743dccd265SPeter Brune vec_sol = snes->vec_sol; 6753dccd265SPeter Brune vec_func = snes->vec_func; 6763dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6773dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6783dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6793dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6803dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6813dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6823dccd265SPeter Brune 6833dccd265SPeter Brune /* reset the number of levels */ 684d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 685cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6863dccd265SPeter Brune 6873dccd265SPeter Brune snes->vec_sol = vec_sol; 6883dccd265SPeter Brune snes->vec_func = vec_func; 6893dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6903dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 691d1adcc6fSPeter Brune } 692d1adcc6fSPeter Brune } 693d1adcc6fSPeter Brune 6943dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6953dccd265SPeter Brune 69607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 697fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 69807144faaSPeter Brune } else { 699ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 70007144faaSPeter Brune } 701cc05f883SPeter Brune 70279d9a41aSPeter Brune if (snes->dm) { 70379d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 70479d9a41aSPeter Brune if (fas->next) { 70579d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 70679d9a41aSPeter Brune if (!fas->next->dm) { 70779d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 70879d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 70979d9a41aSPeter Brune } 71079d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 71179d9a41aSPeter Brune if (!fas->interpolate) { 71279d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 713*bccf9bb3SJed Brown if (!fas->restrct) { 714*bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 71579d9a41aSPeter Brune fas->restrct = fas->interpolate; 71679d9a41aSPeter Brune } 717*bccf9bb3SJed Brown } 71879d9a41aSPeter Brune /* set the injection from the DM */ 71979d9a41aSPeter Brune if (!fas->inject) { 72079d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 72179d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 72279d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 72379d9a41aSPeter Brune } 72479d9a41aSPeter Brune } 72579d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 72679d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 72779d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 72879d9a41aSPeter Brune } 72979d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 73079d9a41aSPeter Brune 73179d9a41aSPeter Brune if (fas->next) { 73279d9a41aSPeter Brune if (fas->galerkin) { 73379d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 73479d9a41aSPeter Brune } else { 73579d9a41aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 73679d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 73779d9a41aSPeter Brune } 73879d9a41aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 73979d9a41aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 74079d9a41aSPeter Brune } 74179d9a41aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 74279d9a41aSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 74379d9a41aSPeter Brune } 74479d9a41aSPeter Brune } 74579d9a41aSPeter Brune } 74679d9a41aSPeter Brune 74779d9a41aSPeter Brune if (fas->next) { 74879d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 74979d9a41aSPeter Brune if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 75079d9a41aSPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 75179d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 75279d9a41aSPeter Brune } 75379d9a41aSPeter Brune 75448bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 75548bfdf8aSPeter Brune if (fas->upsmooth) { 75648bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 75748bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 75848bfdf8aSPeter Brune } 75948bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 76048bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 76148bfdf8aSPeter Brune } 76248bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 76348bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 76448bfdf8aSPeter Brune } 765fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 76648bfdf8aSPeter Brune } 76748bfdf8aSPeter Brune if (fas->downsmooth) { 76848bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 76948bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 77048bfdf8aSPeter Brune } 77148bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 77248bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 77348bfdf8aSPeter Brune } 77448bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 77548bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 77648bfdf8aSPeter Brune } 777fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 7781a266240SBarry Smith } 779cc05f883SPeter Brune 7806273346dSPeter Brune /* setup FAS work vectors */ 7816273346dSPeter Brune if (fas->galerkin) { 7826273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7836273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7846273346dSPeter Brune } 7856273346dSPeter Brune 78679d9a41aSPeter Brune 787fa9694d7SPeter Brune /* got to set them all up at once */ 788421d9b32SPeter Brune PetscFunctionReturn(0); 789421d9b32SPeter Brune } 790421d9b32SPeter Brune 791421d9b32SPeter Brune #undef __FUNCT__ 792421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 793421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 794421d9b32SPeter Brune { 795ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 796ee78dd50SPeter Brune PetscInt levels = 1; 797fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 798421d9b32SPeter Brune PetscErrorCode ierr; 799ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 80007144faaSPeter Brune SNESFASType fastype; 801fde0ff24SPeter Brune const char *optionsprefix; 802fde0ff24SPeter Brune const char *prefix; 803421d9b32SPeter Brune 804421d9b32SPeter Brune PetscFunctionBegin; 805c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 806ee78dd50SPeter Brune 807ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 808ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 809ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 810c732cbdbSBarry Smith if (!flg && snes->dm) { 811c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 812c732cbdbSBarry Smith levels++; 813d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 814c732cbdbSBarry Smith } 815ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 816ee78dd50SPeter Brune } 81707144faaSPeter Brune fastype = fas->fastype; 81807144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 81907144faaSPeter Brune if (flg) { 82007144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 82107144faaSPeter Brune } 822ee78dd50SPeter Brune 823fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 824fde0ff24SPeter Brune 825fde0ff24SPeter Brune /* smoother setup options */ 826fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 827fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 828fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 829fde0ff24SPeter Brune if (fas->level == 0) { 830fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 831fde0ff24SPeter Brune } 832ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 833fde0ff24SPeter Brune 834d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 835ee78dd50SPeter Brune 8366273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8376273346dSPeter Brune 838646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 839ee78dd50SPeter Brune 840ee78dd50SPeter Brune 841421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8428cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 843eff52c0eSPeter Brune 844fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 845fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 846fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 847fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 848fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 849fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 850fde0ff24SPeter Brune } 851fde0ff24SPeter Brune 852fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8538cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8548cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8558cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 856fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8578cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8588cc86e31SPeter Brune } 8598cc86e31SPeter Brune 860fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 86167339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 862ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 86367339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 864fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 865293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 866ee78dd50SPeter Brune } 867fde0ff24SPeter Brune 868fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 869fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 870fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 871fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 872fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 873fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 874fde0ff24SPeter Brune } 875fde0ff24SPeter Brune 876fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 877fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 878fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 879fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 880fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 881fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 882fde0ff24SPeter Brune } 883fde0ff24SPeter Brune 8841a266240SBarry Smith if (fas->upsmooth) { 885fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8861a266240SBarry Smith } 8871a266240SBarry Smith 8881a266240SBarry Smith if (fas->downsmooth) { 889fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8901a266240SBarry Smith } 891ee78dd50SPeter Brune 892742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 893fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 894742fe5e2SPeter Brune } 895742fe5e2SPeter Brune 896ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 897ce11c761SPeter Brune if (!fas->downsmooth) { 89879d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 89993dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 900ce11c761SPeter Brune if (fas->level == 0) { 90179d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 902ce11c761SPeter Brune } 903ce11c761SPeter Brune } 904ce11c761SPeter Brune 905ce11c761SPeter Brune if (!fas->upsmooth) { 90679d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 90793dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 908ce11c761SPeter Brune } 909ce11c761SPeter Brune 910ee78dd50SPeter Brune if (monflg) { 911646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 912794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 9132f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 914742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 915293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 916293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 917d28543b3SPeter Brune } else { 918d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 919d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 920d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 921d28543b3SPeter Brune } 922ee78dd50SPeter Brune } 923ee78dd50SPeter Brune 924ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 925d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 926421d9b32SPeter Brune PetscFunctionReturn(0); 927421d9b32SPeter Brune } 928421d9b32SPeter Brune 929421d9b32SPeter Brune #undef __FUNCT__ 930421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 931421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 932421d9b32SPeter Brune { 933421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 934421d9b32SPeter Brune PetscBool iascii; 935421d9b32SPeter Brune PetscErrorCode ierr; 936421d9b32SPeter Brune 937421d9b32SPeter Brune PetscFunctionBegin; 938421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 939421d9b32SPeter Brune if (iascii) { 940421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 941421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 942421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 943ee78dd50SPeter Brune if (fas->upsmooth) { 94439a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 945421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 946ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 947421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 948421d9b32SPeter Brune } else { 949ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 950421d9b32SPeter Brune } 951ee78dd50SPeter Brune if (fas->downsmooth) { 95239a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 953421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 954ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 955421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 956421d9b32SPeter Brune } else { 957ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 958421d9b32SPeter Brune } 959421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 960421d9b32SPeter Brune } else { 961421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 962421d9b32SPeter Brune } 963421d9b32SPeter Brune PetscFunctionReturn(0); 964421d9b32SPeter Brune } 965421d9b32SPeter Brune 966421d9b32SPeter Brune #undef __FUNCT__ 96739bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 96839bd7f45SPeter Brune /* 96939bd7f45SPeter Brune Defines the action of the downsmoother 97039bd7f45SPeter Brune */ 97139bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 97239bd7f45SPeter Brune PetscErrorCode ierr = 0; 973fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 974742fe5e2SPeter Brune SNESConvergedReason reason; 97539bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 976fde0ff24SPeter Brune Vec G, W, Y; 977fde0ff24SPeter Brune PetscBool lssuccess; 978fde0ff24SPeter Brune PetscInt k; 979421d9b32SPeter Brune PetscFunctionBegin; 980fde0ff24SPeter Brune G = snes->work[1]; 981fde0ff24SPeter Brune W = snes->work[2]; 982fde0ff24SPeter Brune Y = snes->work[3]; 983d1adcc6fSPeter Brune if (fas->downsmooth) { 984d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 985742fe5e2SPeter Brune /* check convergence reason for the smoother */ 986742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 987742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 988742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 989742fe5e2SPeter Brune PetscFunctionReturn(0); 990742fe5e2SPeter Brune } 991fde0ff24SPeter Brune } else { 992fde0ff24SPeter Brune /* basic richardson smoother */ 993fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 994794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 995794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 996fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 997fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 998fde0ff24SPeter Brune if (!lssuccess) { 999fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1000fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 10012f7ea302SPeter Brune PetscFunctionReturn(0); 10022f7ea302SPeter Brune } 1003fe6f9142SPeter Brune } 1004fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1005fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1006fde0ff24SPeter Brune fnorm = gnorm; 1007fde0ff24SPeter Brune } 1008fde0ff24SPeter Brune } 100939bd7f45SPeter Brune PetscFunctionReturn(0); 101039bd7f45SPeter Brune } 101139bd7f45SPeter Brune 101239bd7f45SPeter Brune 101339bd7f45SPeter Brune #undef __FUNCT__ 101439bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 101539bd7f45SPeter Brune /* 101607144faaSPeter Brune Defines the action of the upsmoother 101739bd7f45SPeter Brune */ 101839bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 101939bd7f45SPeter Brune PetscErrorCode ierr = 0; 1020fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 102139bd7f45SPeter Brune SNESConvergedReason reason; 102239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 1023fde0ff24SPeter Brune Vec G, W, Y; 1024fde0ff24SPeter Brune PetscBool lssuccess; 1025fde0ff24SPeter Brune PetscInt k; 102639bd7f45SPeter Brune PetscFunctionBegin; 1027fde0ff24SPeter Brune G = snes->work[1]; 1028fde0ff24SPeter Brune W = snes->work[2]; 1029fde0ff24SPeter Brune Y = snes->work[3]; 103039bd7f45SPeter Brune if (fas->upsmooth) { 1031fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 103239bd7f45SPeter Brune /* check convergence reason for the smoother */ 1033fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 103439bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 103539bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 103639bd7f45SPeter Brune PetscFunctionReturn(0); 103739bd7f45SPeter Brune } 1038fde0ff24SPeter Brune } else { 1039fde0ff24SPeter Brune /* basic richardson smoother */ 1040fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 104139bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 104239bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1043fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1044fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1045fde0ff24SPeter Brune if (!lssuccess) { 1046fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1047fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 104839bd7f45SPeter Brune PetscFunctionReturn(0); 104939bd7f45SPeter Brune } 105039bd7f45SPeter Brune } 1051fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1052fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1053fde0ff24SPeter Brune fnorm = gnorm; 1054fde0ff24SPeter Brune } 1055fde0ff24SPeter Brune } 105639bd7f45SPeter Brune PetscFunctionReturn(0); 105739bd7f45SPeter Brune } 105839bd7f45SPeter Brune 105939bd7f45SPeter Brune #undef __FUNCT__ 106039bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 106139bd7f45SPeter Brune /* 106239bd7f45SPeter Brune 106339bd7f45SPeter Brune Performs the FAS coarse correction as: 106439bd7f45SPeter Brune 106539bd7f45SPeter Brune fine problem: F(x) = 0 106639bd7f45SPeter Brune coarse problem: F^c(x) = b^c 106739bd7f45SPeter Brune 106839bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 106939bd7f45SPeter Brune 107039bd7f45SPeter Brune with correction: 107139bd7f45SPeter Brune 107239bd7f45SPeter Brune 107339bd7f45SPeter Brune 107439bd7f45SPeter Brune */ 107539a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 107639bd7f45SPeter Brune PetscErrorCode ierr; 107739bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 107839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 107939bd7f45SPeter Brune SNESConvergedReason reason; 108039bd7f45SPeter Brune PetscFunctionBegin; 1081fa9694d7SPeter Brune if (fas->next) { 1082ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1083794bee33SPeter Brune 1084c90fad12SPeter Brune X_c = fas->next->vec_sol; 1085293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1086c90fad12SPeter Brune F_c = fas->next->vec_func; 1087742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1088efe1f98aSPeter Brune 1089efe1f98aSPeter Brune /* inject the solution */ 1090efe1f98aSPeter Brune if (fas->inject) { 1091a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 1092efe1f98aSPeter Brune } else { 1093a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 1094a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1095efe1f98aSPeter Brune } 1096293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1097293a7e31SPeter Brune 1098293a7e31SPeter Brune /* restrict the defect */ 1099293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1100293a7e31SPeter Brune 1101c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1102ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1103c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1104742fe5e2SPeter Brune 1105293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1106c90fad12SPeter Brune 1107ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1108ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1109c90fad12SPeter Brune 1110c90fad12SPeter Brune /* recurse to the next level */ 1111f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1112742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1113742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1114742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1115742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1116742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1117742fe5e2SPeter Brune PetscFunctionReturn(0); 1118742fe5e2SPeter Brune } 1119742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1120ee78dd50SPeter Brune 1121fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1122fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 112339bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1124293a7e31SPeter Brune } 112539bd7f45SPeter Brune PetscFunctionReturn(0); 112639bd7f45SPeter Brune } 112739bd7f45SPeter Brune 112839bd7f45SPeter Brune #undef __FUNCT__ 112939bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 113039bd7f45SPeter Brune /* 113139bd7f45SPeter Brune 113239bd7f45SPeter Brune The additive cycle looks like: 113339bd7f45SPeter Brune 113407144faaSPeter Brune xhat = x 113507144faaSPeter Brune xhat = dS(x, b) 113607144faaSPeter Brune x = coarsecorrection(xhat, b_d) 113707144faaSPeter Brune x = x + nu*(xhat - x); 113839bd7f45SPeter Brune (optional) x = uS(x, b) 113939bd7f45SPeter Brune 114039bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 114139bd7f45SPeter Brune 114239bd7f45SPeter Brune */ 114339bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 114407144faaSPeter Brune Vec F, B, Xhat; 1145ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 114639bd7f45SPeter Brune PetscErrorCode ierr; 114707144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 114807144faaSPeter Brune SNESConvergedReason reason; 1149ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1150ddebd997SPeter Brune PetscBool lssucceed; 115139bd7f45SPeter Brune PetscFunctionBegin; 115239bd7f45SPeter Brune 115339bd7f45SPeter Brune F = snes->vec_func; 115439bd7f45SPeter Brune B = snes->vec_rhs; 1155fde0ff24SPeter Brune Xhat = snes->work[3]; 1156fde0ff24SPeter Brune G = snes->work[1]; 1157fde0ff24SPeter Brune W = snes->work[2]; 115807144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 115907144faaSPeter Brune /* recurse first */ 116007144faaSPeter Brune if (fas->next) { 116107144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 116239bd7f45SPeter Brune 116307144faaSPeter Brune X_c = fas->next->vec_sol; 116407144faaSPeter Brune Xo_c = fas->next->work[0]; 116507144faaSPeter Brune F_c = fas->next->vec_func; 116607144faaSPeter Brune B_c = fas->next->vec_rhs; 116739bd7f45SPeter Brune 116807144faaSPeter Brune /* inject the solution */ 116907144faaSPeter Brune if (fas->inject) { 117007144faaSPeter Brune ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr); 117107144faaSPeter Brune } else { 117207144faaSPeter Brune ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr); 117307144faaSPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 117407144faaSPeter Brune } 117507144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 117607144faaSPeter Brune 117707144faaSPeter Brune /* restrict the defect */ 117807144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 117907144faaSPeter Brune 118007144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 118107144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 118207144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 118307144faaSPeter Brune 118407144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 118507144faaSPeter Brune 118607144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 118707144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 118807144faaSPeter Brune 118907144faaSPeter Brune /* recurse */ 119007144faaSPeter Brune fas->next->vec_rhs = B_c; 119107144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 119207144faaSPeter Brune 119307144faaSPeter Brune /* smooth on this level */ 119407144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 119507144faaSPeter Brune 119607144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 119707144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 119807144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 119907144faaSPeter Brune PetscFunctionReturn(0); 120007144faaSPeter Brune } 120107144faaSPeter Brune 120207144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 120307144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1204ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 120507144faaSPeter Brune 1206ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1207ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1208ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1209eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 1210ddebd997SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1211ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1212ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1213ddebd997SPeter Brune fnorm = gnorm; 121407144faaSPeter Brune } else { 121507144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 121607144faaSPeter Brune } 121739bd7f45SPeter Brune PetscFunctionReturn(0); 121839bd7f45SPeter Brune } 121939bd7f45SPeter Brune 122039bd7f45SPeter Brune #undef __FUNCT__ 122139bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 122239bd7f45SPeter Brune /* 122339bd7f45SPeter Brune 122439bd7f45SPeter Brune Defines the FAS cycle as: 122539bd7f45SPeter Brune 122639bd7f45SPeter Brune fine problem: F(x) = 0 122739bd7f45SPeter Brune coarse problem: F^c(x) = b^c 122839bd7f45SPeter Brune 122939bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 123039bd7f45SPeter Brune 123139bd7f45SPeter Brune correction: 123239bd7f45SPeter Brune 123339bd7f45SPeter Brune x = x + I(x^c - Rx) 123439bd7f45SPeter Brune 123539bd7f45SPeter Brune */ 123639bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 123739bd7f45SPeter Brune 123839bd7f45SPeter Brune PetscErrorCode ierr; 123939bd7f45SPeter Brune Vec F,B; 124039bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 124139bd7f45SPeter Brune 124239bd7f45SPeter Brune PetscFunctionBegin; 124339bd7f45SPeter Brune F = snes->vec_func; 124439bd7f45SPeter Brune B = snes->vec_rhs; 124539bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 124639bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 124739bd7f45SPeter Brune 124839a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 124939bd7f45SPeter Brune 1250c90fad12SPeter Brune if (fas->level != 0) { 125139bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1252fe6f9142SPeter Brune } 1253fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1254fa9694d7SPeter Brune 1255fa9694d7SPeter Brune PetscFunctionReturn(0); 1256421d9b32SPeter Brune } 1257421d9b32SPeter Brune 1258421d9b32SPeter Brune #undef __FUNCT__ 1259421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1260421d9b32SPeter Brune 1261421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1262421d9b32SPeter Brune { 1263fa9694d7SPeter Brune PetscErrorCode ierr; 1264fe6f9142SPeter Brune PetscInt i, maxits; 1265ddb5aff1SPeter Brune Vec X, F; 1266fe6f9142SPeter Brune PetscReal fnorm; 126707144faaSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 1268421d9b32SPeter Brune PetscFunctionBegin; 1269fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1270fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1271fa9694d7SPeter Brune X = snes->vec_sol; 1272f5a6d4f9SBarry Smith F = snes->vec_func; 1273293a7e31SPeter Brune 1274293a7e31SPeter Brune /*norm setup */ 1275fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1276fe6f9142SPeter Brune snes->iter = 0; 1277fe6f9142SPeter Brune snes->norm = 0.; 1278fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1279fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1280fe6f9142SPeter Brune if (snes->domainerror) { 1281fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1282fe6f9142SPeter Brune PetscFunctionReturn(0); 1283fe6f9142SPeter Brune } 1284fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1285fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1286fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1287fe6f9142SPeter Brune snes->norm = fnorm; 1288fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1289fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1290fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1291fe6f9142SPeter Brune 1292fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1293fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1294fe6f9142SPeter Brune /* test convergence */ 1295fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1296fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1297fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1298fe6f9142SPeter Brune /* Call general purpose update function */ 1299646217ecSPeter Brune 1300fe6f9142SPeter Brune if (snes->ops->update) { 1301fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1302fe6f9142SPeter Brune } 130307144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 130439bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 130507144faaSPeter Brune } else { 130607144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 130707144faaSPeter Brune } 1308742fe5e2SPeter Brune 1309742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1310742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1311742fe5e2SPeter Brune PetscFunctionReturn(0); 1312742fe5e2SPeter Brune } 1313c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1314c90fad12SPeter Brune /* Monitor convergence */ 1315c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1316c90fad12SPeter Brune snes->iter = i+1; 1317c90fad12SPeter Brune snes->norm = fnorm; 1318c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1319c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1320c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1321c90fad12SPeter Brune /* Test for convergence */ 1322c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1323c90fad12SPeter Brune if (snes->reason) break; 1324fe6f9142SPeter Brune } 1325fe6f9142SPeter Brune if (i == maxits) { 1326fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1327fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1328fe6f9142SPeter Brune } 1329421d9b32SPeter Brune PetscFunctionReturn(0); 1330421d9b32SPeter Brune } 1331