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 700e444f03SPeter Brune snes->max_funcs = 30000; 710e444f03SPeter Brune snes->max_its = 10000; 720e444f03SPeter Brune 73421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 74421d9b32SPeter Brune snes->data = (void*) fas; 75421d9b32SPeter Brune fas->level = 0; 76293a7e31SPeter Brune fas->levels = 1; 77ee78dd50SPeter Brune fas->n_cycles = 1; 78ee78dd50SPeter Brune fas->max_up_it = 1; 79ee78dd50SPeter Brune fas->max_down_it = 1; 80ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 81ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 82421d9b32SPeter Brune fas->next = PETSC_NULL; 836273346dSPeter Brune fas->previous = PETSC_NULL; 84421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 85421d9b32SPeter Brune fas->restrct = PETSC_NULL; 86efe1f98aSPeter Brune fas->inject = PETSC_NULL; 87cc05f883SPeter Brune fas->monitor = PETSC_NULL; 88cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 89ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 90ddebd997SPeter Brune 91ddebd997SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 92ddebd997SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 93ddebd997SPeter Brune 94421d9b32SPeter Brune PetscFunctionReturn(0); 95421d9b32SPeter Brune } 96421d9b32SPeter Brune EXTERN_C_END 97421d9b32SPeter Brune 98421d9b32SPeter Brune #undef __FUNCT__ 99421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 100c79ef259SPeter Brune /*@ 1012e8ce248SJed Brown SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 102c79ef259SPeter Brune 103c79ef259SPeter Brune Input Parameter: 1042e8ce248SJed Brown . snes - the nonlinear solver context 105c79ef259SPeter Brune 106c79ef259SPeter Brune Output parameter: 107c79ef259SPeter Brune . levels - the number of levels 108c79ef259SPeter Brune 109c79ef259SPeter Brune Level: advanced 110c79ef259SPeter Brune 111c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 112c79ef259SPeter Brune 113c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 114c79ef259SPeter Brune @*/ 115421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 116421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 117421d9b32SPeter Brune PetscFunctionBegin; 118ee1fd11aSPeter Brune *levels = fas->levels; 119421d9b32SPeter Brune PetscFunctionReturn(0); 120421d9b32SPeter Brune } 121421d9b32SPeter Brune 122421d9b32SPeter Brune #undef __FUNCT__ 123646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 124c79ef259SPeter Brune /*@ 1252e8ce248SJed Brown SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 126c79ef259SPeter Brune complicated cycling. 127c79ef259SPeter Brune 128c79ef259SPeter Brune Logically Collective on SNES 129c79ef259SPeter Brune 130c79ef259SPeter Brune Input Parameters: 131c79ef259SPeter Brune + snes - the multigrid context 132c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 133c79ef259SPeter Brune 134c79ef259SPeter Brune Options Database Key: 135c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 136c79ef259SPeter Brune 137c79ef259SPeter Brune Level: advanced 138c79ef259SPeter Brune 139c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 140c79ef259SPeter Brune 141c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 142c79ef259SPeter Brune @*/ 143646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 144646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 145646217ecSPeter Brune PetscErrorCode ierr; 146646217ecSPeter Brune PetscFunctionBegin; 147646217ecSPeter Brune fas->n_cycles = cycles; 148646217ecSPeter Brune if (fas->next) { 149646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 150646217ecSPeter Brune } 151646217ecSPeter Brune PetscFunctionReturn(0); 152646217ecSPeter Brune } 153646217ecSPeter Brune 154eff52c0eSPeter Brune #undef __FUNCT__ 155c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 156c79ef259SPeter Brune /*@ 157722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 158c79ef259SPeter Brune 159c79ef259SPeter Brune Logically Collective on SNES 160c79ef259SPeter Brune 161c79ef259SPeter Brune Input Parameters: 162c79ef259SPeter Brune + snes - the multigrid context 163c79ef259SPeter Brune . level - the level to set the number of cycles on 164c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 165c79ef259SPeter Brune 166c79ef259SPeter Brune Level: advanced 167c79ef259SPeter Brune 168c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 169c79ef259SPeter Brune 170c79ef259SPeter Brune .seealso: SNESFASSetCycles() 171c79ef259SPeter Brune @*/ 172c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 173c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 174c79ef259SPeter Brune PetscInt top_level = fas->level,i; 175c79ef259SPeter Brune 176c79ef259SPeter Brune PetscFunctionBegin; 1772e8ce248SJed 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); 178c79ef259SPeter Brune /* get to the correct level */ 179c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 180c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 181c79ef259SPeter Brune } 1822e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 183c79ef259SPeter Brune fas->n_cycles = cycles; 184c79ef259SPeter Brune PetscFunctionReturn(0); 185c79ef259SPeter Brune } 186c79ef259SPeter Brune 187c79ef259SPeter Brune #undef __FUNCT__ 188eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 189aeed3662SMatthew G Knepley /*@C 190c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 191c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 192c79ef259SPeter Brune and nonlinear preconditioners. 193c79ef259SPeter Brune 194c79ef259SPeter Brune Logically Collective on SNES 195c79ef259SPeter Brune 196c79ef259SPeter Brune Input Parameters: 197c79ef259SPeter Brune + snes - the multigrid context 198c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 199c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 200c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 201c79ef259SPeter Brune 202c79ef259SPeter Brune Level: advanced 203c79ef259SPeter Brune 204c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 205c79ef259SPeter Brune 206c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 207c79ef259SPeter Brune @*/ 208eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 209eff52c0eSPeter Brune PetscErrorCode ierr = 0; 210eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 211eff52c0eSPeter Brune PetscFunctionBegin; 212eff52c0eSPeter Brune 213eff52c0eSPeter Brune if (gsfunc) { 214eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 215eff52c0eSPeter Brune /* push the provided GS up the tree */ 216eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 217*6cab3a1bSJed Brown } else { 218*6cab3a1bSJed Brown ierr = SNESGetGS(snes,&gsfunc,&ctx);CHKERRQ(ierr); 219*6cab3a1bSJed Brown if (gsfunc) { 220eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 221eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 222eff52c0eSPeter Brune } else if (use_gs) { 2232e8ce248SJed Brown SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level); 224eff52c0eSPeter Brune } 225*6cab3a1bSJed Brown } 226eff52c0eSPeter Brune PetscFunctionReturn(0); 227eff52c0eSPeter Brune } 228eff52c0eSPeter Brune 229eff52c0eSPeter Brune #undef __FUNCT__ 230eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 231aeed3662SMatthew G Knepley /*@C 232c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 233c79ef259SPeter Brune 234c79ef259SPeter Brune Logically Collective on SNES 235c79ef259SPeter Brune 236c79ef259SPeter Brune Input Parameters: 237c79ef259SPeter Brune + snes - the multigrid context 238c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 239c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 240c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 241c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 242c79ef259SPeter Brune 243c79ef259SPeter Brune Level: advanced 244c79ef259SPeter Brune 245c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 246c79ef259SPeter Brune 247c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 248c79ef259SPeter Brune @*/ 249eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 250eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 251eff52c0eSPeter Brune PetscErrorCode ierr; 252eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 253eff52c0eSPeter Brune SNES cur_snes = snes; 254eff52c0eSPeter Brune PetscFunctionBegin; 2552e8ce248SJed 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); 256eff52c0eSPeter Brune /* get to the correct level */ 257eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 258eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 259eff52c0eSPeter Brune cur_snes = fas->next; 260eff52c0eSPeter Brune } 2612e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 262eff52c0eSPeter Brune if (gsfunc) { 2636273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 264eff52c0eSPeter Brune } 265eff52c0eSPeter Brune PetscFunctionReturn(0); 266eff52c0eSPeter Brune } 267eff52c0eSPeter Brune 268646217ecSPeter Brune #undef __FUNCT__ 269421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 270c79ef259SPeter Brune /*@ 271c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 272c79ef259SPeter Brune level of the FAS hierarchy. 273c79ef259SPeter Brune 274c79ef259SPeter Brune Input Parameters: 275c79ef259SPeter Brune + snes - the multigrid context 276c79ef259SPeter Brune level - the level to get 277c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 278c79ef259SPeter Brune 279c79ef259SPeter Brune Level: advanced 280c79ef259SPeter Brune 281c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 282c79ef259SPeter Brune 283c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 284c79ef259SPeter Brune @*/ 285421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) { 286421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 287421d9b32SPeter Brune PetscInt i; 2882e8ce248SJed Brown 289421d9b32SPeter Brune PetscFunctionBegin; 2902e8ce248SJed 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); 2912e8ce248SJed 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); 292421d9b32SPeter Brune *lsnes = snes; 293421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 294421d9b32SPeter Brune *lsnes = fas->next; 295421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 296421d9b32SPeter Brune } 2972e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 298421d9b32SPeter Brune PetscFunctionReturn(0); 299421d9b32SPeter Brune } 300421d9b32SPeter Brune 301421d9b32SPeter Brune #undef __FUNCT__ 30207144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 30307144faaSPeter Brune /*@ 30407144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 30507144faaSPeter Brune e 30607144faaSPeter Brune 30707144faaSPeter Brune 30807144faaSPeter Brune @*/ 30907144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 31007144faaSPeter Brune { 31107144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31207144faaSPeter Brune 31307144faaSPeter Brune PetscFunctionBegin; 31407144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31507144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 31607144faaSPeter Brune fas->fastype = fastype; 31707144faaSPeter Brune PetscFunctionReturn(0); 31807144faaSPeter Brune } 31907144faaSPeter Brune 32007144faaSPeter Brune #undef __FUNCT__ 321421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 322aeed3662SMatthew G Knepley /*@C 323c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 324c79ef259SPeter Brune Must be called before any other FAS routine. 325c79ef259SPeter Brune 326c79ef259SPeter Brune Input Parameters: 327c79ef259SPeter Brune + snes - the snes context 328c79ef259SPeter Brune . levels - the number of levels 329c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 330c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 331c79ef259SPeter Brune Fortran. 332c79ef259SPeter Brune 333c79ef259SPeter Brune Level: intermediate 334c79ef259SPeter Brune 335c79ef259SPeter Brune Notes: 336c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 337c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 338c79ef259SPeter Brune 339c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 340c79ef259SPeter Brune 341c79ef259SPeter Brune .seealso: SNESFASGetLevels() 342c79ef259SPeter Brune @*/ 343421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 344421d9b32SPeter Brune PetscErrorCode ierr; 345421d9b32SPeter Brune PetscInt i; 346421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3476273346dSPeter Brune SNES prevsnes; 348421d9b32SPeter Brune MPI_Comm comm; 349421d9b32SPeter Brune PetscFunctionBegin; 350ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 351ee1fd11aSPeter Brune if (levels == fas->levels) { 352ee1fd11aSPeter Brune if (!comms) 353ee1fd11aSPeter Brune PetscFunctionReturn(0); 354ee1fd11aSPeter Brune } 355421d9b32SPeter Brune /* user has changed the number of levels; reset */ 356421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 357421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 358421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 359ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3606273346dSPeter Brune fas->previous = PETSC_NULL; 3616273346dSPeter Brune prevsnes = snes; 362421d9b32SPeter Brune /* setup the finest level */ 363421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 364421d9b32SPeter Brune if (comms) comm = comms[i]; 365421d9b32SPeter Brune fas->level = i; 366421d9b32SPeter Brune fas->levels = levels; 367ee1fd11aSPeter Brune fas->next = PETSC_NULL; 368421d9b32SPeter Brune if (i > 0) { 369421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3704a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 371421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 372794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3736273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3746273346dSPeter Brune prevsnes = fas->next; 3756273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 376421d9b32SPeter Brune } 377421d9b32SPeter Brune } 378421d9b32SPeter Brune PetscFunctionReturn(0); 379421d9b32SPeter Brune } 380421d9b32SPeter Brune 381421d9b32SPeter Brune #undef __FUNCT__ 382c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 383c79ef259SPeter Brune /*@ 384c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 385c79ef259SPeter Brune use on all levels. 386c79ef259SPeter Brune 387fde0ff24SPeter Brune Logically Collective on SNES 388c79ef259SPeter Brune 389c79ef259SPeter Brune Input Parameters: 390c79ef259SPeter Brune + snes - the multigrid context 391c79ef259SPeter Brune - n - the number of smoothing steps 392c79ef259SPeter Brune 393c79ef259SPeter Brune Options Database Key: 394d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 395c79ef259SPeter Brune 396c79ef259SPeter Brune Level: advanced 397c79ef259SPeter Brune 398c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 399c79ef259SPeter Brune 400c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 401c79ef259SPeter Brune @*/ 402c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 403c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 404c79ef259SPeter Brune PetscErrorCode ierr = 0; 405c79ef259SPeter Brune PetscFunctionBegin; 406c79ef259SPeter Brune fas->max_up_it = n; 407c79ef259SPeter Brune if (fas->next) { 408c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 409c79ef259SPeter Brune } 410c79ef259SPeter Brune PetscFunctionReturn(0); 411c79ef259SPeter Brune } 412c79ef259SPeter Brune 413c79ef259SPeter Brune #undef __FUNCT__ 414c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 415c79ef259SPeter Brune /*@ 416c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 417c79ef259SPeter Brune use on all levels. 418c79ef259SPeter Brune 419fde0ff24SPeter Brune Logically Collective on SNES 420c79ef259SPeter Brune 421c79ef259SPeter Brune Input Parameters: 422c79ef259SPeter Brune + snes - the multigrid context 423c79ef259SPeter Brune - n - the number of smoothing steps 424c79ef259SPeter Brune 425c79ef259SPeter Brune Options Database Key: 426d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 427c79ef259SPeter Brune 428c79ef259SPeter Brune Level: advanced 429c79ef259SPeter Brune 430c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 431c79ef259SPeter Brune 432c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 433c79ef259SPeter Brune @*/ 434c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 435c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 436c79ef259SPeter Brune PetscErrorCode ierr = 0; 437c79ef259SPeter Brune PetscFunctionBegin; 438c79ef259SPeter Brune fas->max_down_it = n; 439c79ef259SPeter Brune if (fas->next) { 440c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 441c79ef259SPeter Brune } 442c79ef259SPeter Brune PetscFunctionReturn(0); 443c79ef259SPeter Brune } 444c79ef259SPeter Brune 445c79ef259SPeter Brune #undef __FUNCT__ 446421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 447c79ef259SPeter Brune /*@ 448c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 449c79ef259SPeter Brune interpolation from l-1 to the lth level 450c79ef259SPeter Brune 451c79ef259SPeter Brune Input Parameters: 452c79ef259SPeter Brune + snes - the multigrid context 453c79ef259SPeter Brune . mat - the interpolation operator 454c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 455c79ef259SPeter Brune 456c79ef259SPeter Brune Level: advanced 457c79ef259SPeter Brune 458c79ef259SPeter Brune Notes: 459c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 460c79ef259SPeter Brune for the same level. 461c79ef259SPeter Brune 462c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 463c79ef259SPeter Brune out from the matrix size which one it is. 464c79ef259SPeter Brune 465c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 466c79ef259SPeter Brune 467bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 468c79ef259SPeter Brune @*/ 469421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 470421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 471d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 472bccf9bb3SJed Brown PetscErrorCode ierr; 473421d9b32SPeter Brune 474421d9b32SPeter Brune PetscFunctionBegin; 4752e8ce248SJed 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); 476421d9b32SPeter Brune /* get to the correct level */ 477d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 478421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 479421d9b32SPeter Brune } 4802e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 481bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 482bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 483421d9b32SPeter Brune fas->interpolate = mat; 484421d9b32SPeter Brune PetscFunctionReturn(0); 485421d9b32SPeter Brune } 486421d9b32SPeter Brune 487421d9b32SPeter Brune #undef __FUNCT__ 488421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 489c79ef259SPeter Brune /*@ 490c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 491c79ef259SPeter Brune from level l to l-1. 492c79ef259SPeter Brune 493c79ef259SPeter Brune Input Parameters: 494c79ef259SPeter Brune + snes - the multigrid context 495c79ef259SPeter Brune . mat - the restriction matrix 496c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 497c79ef259SPeter Brune 498c79ef259SPeter Brune Level: advanced 499c79ef259SPeter Brune 500c79ef259SPeter Brune Notes: 501c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 502c79ef259SPeter Brune for the same level. 503c79ef259SPeter Brune 504c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 505c79ef259SPeter Brune out from the matrix size which one it is. 506c79ef259SPeter Brune 507fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 508c79ef259SPeter Brune is used. 509c79ef259SPeter Brune 510c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 511c79ef259SPeter Brune 512c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 513c79ef259SPeter Brune @*/ 514421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 515421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 516d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 517bccf9bb3SJed Brown PetscErrorCode ierr; 518421d9b32SPeter Brune 519421d9b32SPeter Brune PetscFunctionBegin; 5202e8ce248SJed 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); 521421d9b32SPeter Brune /* get to the correct level */ 522d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 523421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 524421d9b32SPeter Brune } 5252e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 526bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 527bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 528421d9b32SPeter Brune fas->restrct = mat; 529421d9b32SPeter Brune PetscFunctionReturn(0); 530421d9b32SPeter Brune } 531421d9b32SPeter Brune 532421d9b32SPeter Brune #undef __FUNCT__ 533421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 534c79ef259SPeter Brune /*@ 535c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 536c79ef259SPeter Brune operator from level l to l-1. 537c79ef259SPeter Brune 538c79ef259SPeter Brune Input Parameters: 539c79ef259SPeter Brune + snes - the multigrid context 540c79ef259SPeter Brune . rscale - the restriction scaling 541c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 542c79ef259SPeter Brune 543c79ef259SPeter Brune Level: advanced 544c79ef259SPeter Brune 545c79ef259SPeter Brune Notes: 546c79ef259SPeter Brune This is only used in the case that the injection is not set. 547c79ef259SPeter Brune 548c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 549c79ef259SPeter Brune 550c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 551c79ef259SPeter Brune @*/ 552421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 553421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 554d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 555bccf9bb3SJed Brown PetscErrorCode ierr; 556421d9b32SPeter Brune 557421d9b32SPeter Brune PetscFunctionBegin; 5582e8ce248SJed 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); 559421d9b32SPeter Brune /* get to the correct level */ 560d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 561421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 562421d9b32SPeter Brune } 5632e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 564bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 565bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 566421d9b32SPeter Brune fas->rscale = rscale; 567421d9b32SPeter Brune PetscFunctionReturn(0); 568421d9b32SPeter Brune } 569421d9b32SPeter Brune 570efe1f98aSPeter Brune 571efe1f98aSPeter Brune #undef __FUNCT__ 572efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 573c79ef259SPeter Brune /*@ 574c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 575c79ef259SPeter Brune from level l to l-1. 576c79ef259SPeter Brune 577c79ef259SPeter Brune Input Parameters: 578c79ef259SPeter Brune + snes - the multigrid context 579c79ef259SPeter Brune . mat - the restriction matrix 580c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 581c79ef259SPeter Brune 582c79ef259SPeter Brune Level: advanced 583c79ef259SPeter Brune 584c79ef259SPeter Brune Notes: 585c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 586c79ef259SPeter Brune project the solution instead. 587c79ef259SPeter Brune 588c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 589c79ef259SPeter Brune 590c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 591c79ef259SPeter Brune @*/ 592efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 593efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 594efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 595bccf9bb3SJed Brown PetscErrorCode ierr; 596efe1f98aSPeter Brune 597efe1f98aSPeter Brune PetscFunctionBegin; 5982e8ce248SJed 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); 599efe1f98aSPeter Brune /* get to the correct level */ 600efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 601efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 602efe1f98aSPeter Brune } 6032e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 604bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 605bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 606efe1f98aSPeter Brune fas->inject = mat; 607efe1f98aSPeter Brune PetscFunctionReturn(0); 608efe1f98aSPeter Brune } 609efe1f98aSPeter Brune 610421d9b32SPeter Brune #undef __FUNCT__ 611421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 612421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 613421d9b32SPeter Brune { 61477df8cc4SPeter Brune PetscErrorCode ierr = 0; 615421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 616421d9b32SPeter Brune 617421d9b32SPeter Brune PetscFunctionBegin; 618bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 619bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 6203dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 621bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 622bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 623bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6243dccd265SPeter Brune 6253dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6263dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 627742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 62809dc8347SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","",PETSC_NULL);CHKERRQ(ierr); 629421d9b32SPeter Brune PetscFunctionReturn(0); 630421d9b32SPeter Brune } 631421d9b32SPeter Brune 632421d9b32SPeter Brune #undef __FUNCT__ 633421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 634421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 635421d9b32SPeter Brune { 636421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 637742fe5e2SPeter Brune PetscErrorCode ierr = 0; 638421d9b32SPeter Brune 639421d9b32SPeter Brune PetscFunctionBegin; 640421d9b32SPeter Brune /* recursively resets and then destroys */ 64179d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 642421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 643421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 644ee1fd11aSPeter Brune 645421d9b32SPeter Brune PetscFunctionReturn(0); 646421d9b32SPeter Brune } 647421d9b32SPeter Brune 648421d9b32SPeter Brune #undef __FUNCT__ 649421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 650421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 651421d9b32SPeter Brune { 65248bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 653421d9b32SPeter Brune PetscErrorCode ierr; 654efe1f98aSPeter Brune VecScatter injscatter; 655d1adcc6fSPeter Brune PetscInt dm_levels; 6563dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 657d1adcc6fSPeter Brune 658421d9b32SPeter Brune PetscFunctionBegin; 659eff52c0eSPeter Brune 660cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 661d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 662d1adcc6fSPeter Brune dm_levels++; 663cc05f883SPeter Brune if (dm_levels > fas->levels) { 6643dccd265SPeter Brune 6652e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6663dccd265SPeter Brune vec_sol = snes->vec_sol; 6673dccd265SPeter Brune vec_func = snes->vec_func; 6683dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6693dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6703dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6713dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6723dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6733dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6743dccd265SPeter Brune 6753dccd265SPeter Brune /* reset the number of levels */ 676d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 677cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6783dccd265SPeter Brune 6793dccd265SPeter Brune snes->vec_sol = vec_sol; 6803dccd265SPeter Brune snes->vec_func = vec_func; 6813dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6823dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 683d1adcc6fSPeter Brune } 684d1adcc6fSPeter Brune } 685d1adcc6fSPeter Brune 6863dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6873dccd265SPeter Brune 68807144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 689fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 69007144faaSPeter Brune } else { 691ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 69207144faaSPeter Brune } 693cc05f883SPeter Brune 69479d9a41aSPeter Brune if (snes->dm) { 69579d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 69679d9a41aSPeter Brune if (fas->next) { 69779d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 69879d9a41aSPeter Brune if (!fas->next->dm) { 69979d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 70079d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 70179d9a41aSPeter Brune } 70279d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 70379d9a41aSPeter Brune if (!fas->interpolate) { 70479d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 705bccf9bb3SJed Brown if (!fas->restrct) { 706bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 70779d9a41aSPeter Brune fas->restrct = fas->interpolate; 70879d9a41aSPeter Brune } 709bccf9bb3SJed Brown } 71079d9a41aSPeter Brune /* set the injection from the DM */ 71179d9a41aSPeter Brune if (!fas->inject) { 71279d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 71379d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 71479d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 71579d9a41aSPeter Brune } 71679d9a41aSPeter Brune } 71779d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 71879d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 71979d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 72079d9a41aSPeter Brune } 72179d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 72279d9a41aSPeter Brune 72379d9a41aSPeter Brune if (fas->next) { 72479d9a41aSPeter Brune if (fas->galerkin) { 72579d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 72679d9a41aSPeter Brune } 72779d9a41aSPeter Brune } 72879d9a41aSPeter Brune 72979d9a41aSPeter Brune if (fas->next) { 73079d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 731938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 732938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 73379d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 73479d9a41aSPeter Brune } 73579d9a41aSPeter Brune 736*6cab3a1bSJed Brown /* setup the pre and post smoothers */ 737*6cab3a1bSJed Brown if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);} 738*6cab3a1bSJed Brown if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);} 739cc05f883SPeter Brune 7406273346dSPeter Brune /* setup FAS work vectors */ 7416273346dSPeter Brune if (fas->galerkin) { 7426273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7436273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7446273346dSPeter Brune } 745421d9b32SPeter Brune PetscFunctionReturn(0); 746421d9b32SPeter Brune } 747421d9b32SPeter Brune 748421d9b32SPeter Brune #undef __FUNCT__ 749421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 750421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 751421d9b32SPeter Brune { 752ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 753ee78dd50SPeter Brune PetscInt levels = 1; 754fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 755421d9b32SPeter Brune PetscErrorCode ierr; 756ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 75707144faaSPeter Brune SNESFASType fastype; 758fde0ff24SPeter Brune const char *optionsprefix; 759fde0ff24SPeter Brune const char *prefix; 760421d9b32SPeter Brune 761421d9b32SPeter Brune PetscFunctionBegin; 762c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 763ee78dd50SPeter Brune 764ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 765ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 766ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 767c732cbdbSBarry Smith if (!flg && snes->dm) { 768c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 769c732cbdbSBarry Smith levels++; 770d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 771c732cbdbSBarry Smith } 772ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 773ee78dd50SPeter Brune } 77407144faaSPeter Brune fastype = fas->fastype; 77507144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 77607144faaSPeter Brune if (flg) { 77707144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 77807144faaSPeter Brune } 779ee78dd50SPeter Brune 780fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 781fde0ff24SPeter Brune 782fde0ff24SPeter Brune /* smoother setup options */ 783fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 784fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 785fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 786fde0ff24SPeter Brune if (fas->level == 0) { 787fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 788fde0ff24SPeter Brune } 789ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 790fde0ff24SPeter Brune 791d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 792ee78dd50SPeter Brune 7936273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 7946273346dSPeter Brune 795646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 796ee78dd50SPeter Brune 797ee78dd50SPeter Brune 798421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7998cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 800eff52c0eSPeter Brune 801fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 802fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 803fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 804fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 805fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 806fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 807fde0ff24SPeter Brune } 808fde0ff24SPeter Brune 809fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8108cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8118cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8128cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 813fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8148cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8158cc86e31SPeter Brune } 8168cc86e31SPeter Brune 817fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 81867339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 819ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 82067339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 821fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 822293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 823ee78dd50SPeter Brune } 824fde0ff24SPeter Brune 825fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 826fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 827fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 828fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 829fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 830fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 831fde0ff24SPeter Brune } 832fde0ff24SPeter Brune 833fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 834fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 835fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 836fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 837fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 838fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 839fde0ff24SPeter Brune } 840fde0ff24SPeter Brune 8411a266240SBarry Smith if (fas->upsmooth) { 842fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8431a266240SBarry Smith } 8441a266240SBarry Smith 8451a266240SBarry Smith if (fas->downsmooth) { 846fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8471a266240SBarry Smith } 848ee78dd50SPeter Brune 849742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 850fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 851742fe5e2SPeter Brune } 852742fe5e2SPeter Brune 853ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 854ce11c761SPeter Brune if (!fas->downsmooth) { 85593dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 8565d115551SPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 857ce11c761SPeter Brune if (fas->level == 0) { 85879d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 859ce11c761SPeter Brune } 860ce11c761SPeter Brune } 861ce11c761SPeter Brune 862ce11c761SPeter Brune if (!fas->upsmooth) { 86393dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 8645d115551SPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 865ce11c761SPeter Brune } 866ce11c761SPeter Brune 867ee78dd50SPeter Brune if (monflg) { 868646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 869794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8702f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 871742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 872293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 873293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 874d28543b3SPeter Brune } else { 875d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 876d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 877d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 878d28543b3SPeter Brune } 879ee78dd50SPeter Brune } 880ee78dd50SPeter Brune 881ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 882d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 883421d9b32SPeter Brune PetscFunctionReturn(0); 884421d9b32SPeter Brune } 885421d9b32SPeter Brune 886421d9b32SPeter Brune #undef __FUNCT__ 887421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 888421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 889421d9b32SPeter Brune { 890421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 891421d9b32SPeter Brune PetscBool iascii; 892421d9b32SPeter Brune PetscErrorCode ierr; 893421d9b32SPeter Brune 894421d9b32SPeter Brune PetscFunctionBegin; 895421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 896421d9b32SPeter Brune if (iascii) { 8972e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 898421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 899bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 900ee78dd50SPeter Brune if (fas->upsmooth) { 90139a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 902421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 903ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 904421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 905421d9b32SPeter Brune } else { 906ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 907421d9b32SPeter Brune } 908ee78dd50SPeter Brune if (fas->downsmooth) { 90939a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 910421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 911ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 912421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 913421d9b32SPeter Brune } else { 914ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 915421d9b32SPeter Brune } 916421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 917421d9b32SPeter Brune } else { 918421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 919421d9b32SPeter Brune } 920421d9b32SPeter Brune PetscFunctionReturn(0); 921421d9b32SPeter Brune } 922421d9b32SPeter Brune 923421d9b32SPeter Brune #undef __FUNCT__ 92439bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92539bd7f45SPeter Brune /* 92639bd7f45SPeter Brune Defines the action of the downsmoother 92739bd7f45SPeter Brune */ 92839bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 92939bd7f45SPeter Brune PetscErrorCode ierr = 0; 930fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 931742fe5e2SPeter Brune SNESConvergedReason reason; 93239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 9334b32a720SPeter Brune Vec G, W, Y, FPC; 934fde0ff24SPeter Brune PetscBool lssuccess; 935fde0ff24SPeter Brune PetscInt k; 936421d9b32SPeter Brune PetscFunctionBegin; 937fde0ff24SPeter Brune G = snes->work[1]; 938fde0ff24SPeter Brune W = snes->work[2]; 939fde0ff24SPeter Brune Y = snes->work[3]; 940d1adcc6fSPeter Brune if (fas->downsmooth) { 941d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 942742fe5e2SPeter Brune /* check convergence reason for the smoother */ 943742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 944742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 945742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 946742fe5e2SPeter Brune PetscFunctionReturn(0); 947742fe5e2SPeter Brune } 9484b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9494b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 950fde0ff24SPeter Brune } else { 951fde0ff24SPeter Brune /* basic richardson smoother */ 952fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 953794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 954794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 955fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 95605b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr); 95705b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 958fde0ff24SPeter Brune if (!lssuccess) { 959fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 960fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9612f7ea302SPeter Brune PetscFunctionReturn(0); 9622f7ea302SPeter Brune } 963fe6f9142SPeter Brune } 964fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 965fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 966fde0ff24SPeter Brune fnorm = gnorm; 967fde0ff24SPeter Brune } 968fde0ff24SPeter Brune } 96939bd7f45SPeter Brune PetscFunctionReturn(0); 97039bd7f45SPeter Brune } 97139bd7f45SPeter Brune 97239bd7f45SPeter Brune 97339bd7f45SPeter Brune #undef __FUNCT__ 97439bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 97539bd7f45SPeter Brune /* 97607144faaSPeter Brune Defines the action of the upsmoother 97739bd7f45SPeter Brune */ 97839bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 97939bd7f45SPeter Brune PetscErrorCode ierr = 0; 980fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 98139bd7f45SPeter Brune SNESConvergedReason reason; 98239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 9834b32a720SPeter Brune Vec G, W, Y, FPC; 984fde0ff24SPeter Brune PetscBool lssuccess; 985fde0ff24SPeter Brune PetscInt k; 98639bd7f45SPeter Brune PetscFunctionBegin; 987fde0ff24SPeter Brune G = snes->work[1]; 988fde0ff24SPeter Brune W = snes->work[2]; 989fde0ff24SPeter Brune Y = snes->work[3]; 99039bd7f45SPeter Brune if (fas->upsmooth) { 991fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 99239bd7f45SPeter Brune /* check convergence reason for the smoother */ 993fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 99439bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 99539bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 99639bd7f45SPeter Brune PetscFunctionReturn(0); 99739bd7f45SPeter Brune } 9984b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9994b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 1000fde0ff24SPeter Brune } else { 1001fde0ff24SPeter Brune /* basic richardson smoother */ 1002fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 100339bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 100439bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1005fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 100605b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr); 100705b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1008fde0ff24SPeter Brune if (!lssuccess) { 1009fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1010fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 101139bd7f45SPeter Brune PetscFunctionReturn(0); 101239bd7f45SPeter Brune } 101339bd7f45SPeter Brune } 1014fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1015fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1016fde0ff24SPeter Brune fnorm = gnorm; 1017fde0ff24SPeter Brune } 1018fde0ff24SPeter Brune } 101939bd7f45SPeter Brune PetscFunctionReturn(0); 102039bd7f45SPeter Brune } 102139bd7f45SPeter Brune 102239bd7f45SPeter Brune #undef __FUNCT__ 1023938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1024938e4a01SJed Brown /*@ 1025938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1026938e4a01SJed Brown 1027938e4a01SJed Brown Collective 1028938e4a01SJed Brown 1029938e4a01SJed Brown Input Arguments: 1030938e4a01SJed Brown . snes - SNESFAS 1031938e4a01SJed Brown 1032938e4a01SJed Brown Output Arguments: 1033938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1034938e4a01SJed Brown 1035938e4a01SJed Brown Level: developer 1036938e4a01SJed Brown 1037938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1038938e4a01SJed Brown @*/ 1039938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1040938e4a01SJed Brown { 1041938e4a01SJed Brown PetscErrorCode ierr; 1042938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1043938e4a01SJed Brown 1044938e4a01SJed Brown PetscFunctionBegin; 1045938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1046938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1047938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1048938e4a01SJed Brown PetscFunctionReturn(0); 1049938e4a01SJed Brown } 1050938e4a01SJed Brown 1051e9923e8dSJed Brown #undef __FUNCT__ 1052e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1053e9923e8dSJed Brown /*@ 1054e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1055e9923e8dSJed Brown 1056e9923e8dSJed Brown Collective 1057e9923e8dSJed Brown 1058e9923e8dSJed Brown Input Arguments: 1059e9923e8dSJed Brown + fine - SNES from which to restrict 1060e9923e8dSJed Brown - Xfine - vector to restrict 1061e9923e8dSJed Brown 1062e9923e8dSJed Brown Output Arguments: 1063e9923e8dSJed Brown . Xcoarse - result of restriction 1064e9923e8dSJed Brown 1065e9923e8dSJed Brown Level: developer 1066e9923e8dSJed Brown 1067e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1068e9923e8dSJed Brown @*/ 1069e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1070e9923e8dSJed Brown { 1071e9923e8dSJed Brown PetscErrorCode ierr; 1072e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1073e9923e8dSJed Brown 1074e9923e8dSJed Brown PetscFunctionBegin; 1075e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1076e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1077e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1078e9923e8dSJed Brown if (fas->inject) { 1079e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1080e9923e8dSJed Brown } else { 1081e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1082e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1083e9923e8dSJed Brown } 1084e9923e8dSJed Brown PetscFunctionReturn(0); 1085e9923e8dSJed Brown } 1086e9923e8dSJed Brown 1087e9923e8dSJed Brown #undef __FUNCT__ 108839bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 108939bd7f45SPeter Brune /* 109039bd7f45SPeter Brune 109139bd7f45SPeter Brune Performs the FAS coarse correction as: 109239bd7f45SPeter Brune 109339bd7f45SPeter Brune fine problem: F(x) = 0 109439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 109539bd7f45SPeter Brune 109639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 109739bd7f45SPeter Brune 109839bd7f45SPeter Brune with correction: 109939bd7f45SPeter Brune 110039bd7f45SPeter Brune 110139bd7f45SPeter Brune 110239bd7f45SPeter Brune */ 110339a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 110439bd7f45SPeter Brune PetscErrorCode ierr; 110539bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 110639bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 110739bd7f45SPeter Brune SNESConvergedReason reason; 110839bd7f45SPeter Brune PetscFunctionBegin; 1109fa9694d7SPeter Brune if (fas->next) { 1110c90fad12SPeter Brune X_c = fas->next->vec_sol; 1111293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1112c90fad12SPeter Brune F_c = fas->next->vec_func; 1113742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1114efe1f98aSPeter Brune 1115938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1116293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1117293a7e31SPeter Brune 1118293a7e31SPeter Brune /* restrict the defect */ 1119293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1120293a7e31SPeter Brune 1121c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1122ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1123c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1124742fe5e2SPeter Brune 1125293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1126c90fad12SPeter Brune 1127ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1128ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1129c90fad12SPeter Brune 1130c90fad12SPeter Brune /* recurse to the next level */ 1131f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1132742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1133742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1134742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1135742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1136742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1137742fe5e2SPeter Brune PetscFunctionReturn(0); 1138742fe5e2SPeter Brune } 1139742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1140ee78dd50SPeter Brune 1141fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1142fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 114339bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1144293a7e31SPeter Brune } 114539bd7f45SPeter Brune PetscFunctionReturn(0); 114639bd7f45SPeter Brune } 114739bd7f45SPeter Brune 114839bd7f45SPeter Brune #undef __FUNCT__ 114939bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 115039bd7f45SPeter Brune /* 115139bd7f45SPeter Brune 115239bd7f45SPeter Brune The additive cycle looks like: 115339bd7f45SPeter Brune 115407144faaSPeter Brune xhat = x 115507144faaSPeter Brune xhat = dS(x, b) 115607144faaSPeter Brune x = coarsecorrection(xhat, b_d) 115707144faaSPeter Brune x = x + nu*(xhat - x); 115839bd7f45SPeter Brune (optional) x = uS(x, b) 115939bd7f45SPeter Brune 116039bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 116139bd7f45SPeter Brune 116239bd7f45SPeter Brune */ 116339bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 116407144faaSPeter Brune Vec F, B, Xhat; 1165ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 116639bd7f45SPeter Brune PetscErrorCode ierr; 116707144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 116807144faaSPeter Brune SNESConvergedReason reason; 1169ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1170ddebd997SPeter Brune PetscBool lssucceed; 117139bd7f45SPeter Brune PetscFunctionBegin; 117239bd7f45SPeter Brune 117339bd7f45SPeter Brune F = snes->vec_func; 117439bd7f45SPeter Brune B = snes->vec_rhs; 1175fde0ff24SPeter Brune Xhat = snes->work[3]; 1176fde0ff24SPeter Brune G = snes->work[1]; 1177fde0ff24SPeter Brune W = snes->work[2]; 117807144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 117907144faaSPeter Brune /* recurse first */ 118007144faaSPeter Brune if (fas->next) { 118107144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 118239bd7f45SPeter Brune 118307144faaSPeter Brune X_c = fas->next->vec_sol; 118407144faaSPeter Brune Xo_c = fas->next->work[0]; 118507144faaSPeter Brune F_c = fas->next->vec_func; 118607144faaSPeter Brune B_c = fas->next->vec_rhs; 118739bd7f45SPeter Brune 1188938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 118907144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 119007144faaSPeter Brune 119107144faaSPeter Brune /* restrict the defect */ 119207144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 119307144faaSPeter Brune 119407144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 119507144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 119607144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 119707144faaSPeter Brune 119807144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 119907144faaSPeter Brune 120007144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 120107144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 120207144faaSPeter Brune 120307144faaSPeter Brune /* recurse */ 120407144faaSPeter Brune fas->next->vec_rhs = B_c; 120507144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 120607144faaSPeter Brune 120707144faaSPeter Brune /* smooth on this level */ 120807144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 120907144faaSPeter Brune 121007144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 121107144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 121207144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 121307144faaSPeter Brune PetscFunctionReturn(0); 121407144faaSPeter Brune } 121507144faaSPeter Brune 121607144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 121707144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1218ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 121907144faaSPeter Brune 1220ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1221ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1222ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1223eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 122405b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes, X,Xhat,PETSC_NULL);CHKERRQ(ierr); 122505b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Xhat,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1226ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1227ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1228ddebd997SPeter Brune fnorm = gnorm; 122907144faaSPeter Brune } else { 123007144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 123107144faaSPeter Brune } 123239bd7f45SPeter Brune PetscFunctionReturn(0); 123339bd7f45SPeter Brune } 123439bd7f45SPeter Brune 123539bd7f45SPeter Brune #undef __FUNCT__ 123639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 123739bd7f45SPeter Brune /* 123839bd7f45SPeter Brune 123939bd7f45SPeter Brune Defines the FAS cycle as: 124039bd7f45SPeter Brune 124139bd7f45SPeter Brune fine problem: F(x) = 0 124239bd7f45SPeter Brune coarse problem: F^c(x) = b^c 124339bd7f45SPeter Brune 124439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 124539bd7f45SPeter Brune 124639bd7f45SPeter Brune correction: 124739bd7f45SPeter Brune 124839bd7f45SPeter Brune x = x + I(x^c - Rx) 124939bd7f45SPeter Brune 125039bd7f45SPeter Brune */ 125139bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 125239bd7f45SPeter Brune 125339bd7f45SPeter Brune PetscErrorCode ierr; 125439bd7f45SPeter Brune Vec F,B; 125539bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 125639bd7f45SPeter Brune 125739bd7f45SPeter Brune PetscFunctionBegin; 125839bd7f45SPeter Brune F = snes->vec_func; 125939bd7f45SPeter Brune B = snes->vec_rhs; 126039bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 126139bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 126239bd7f45SPeter Brune 126339a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 126439bd7f45SPeter Brune 1265c90fad12SPeter Brune if (fas->level != 0) { 126639bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1267fe6f9142SPeter Brune } 1268fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1269fa9694d7SPeter Brune 1270fa9694d7SPeter Brune PetscFunctionReturn(0); 1271421d9b32SPeter Brune } 1272421d9b32SPeter Brune 1273421d9b32SPeter Brune #undef __FUNCT__ 1274421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1275421d9b32SPeter Brune 1276421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1277421d9b32SPeter Brune { 1278fa9694d7SPeter Brune PetscErrorCode ierr; 1279fe6f9142SPeter Brune PetscInt i, maxits; 1280ddb5aff1SPeter Brune Vec X, F; 1281fe6f9142SPeter Brune PetscReal fnorm; 1282b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1283b17ce1afSJed Brown DM dm; 1284b17ce1afSJed Brown 1285421d9b32SPeter Brune PetscFunctionBegin; 1286fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1287fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1288fa9694d7SPeter Brune X = snes->vec_sol; 1289f5a6d4f9SBarry Smith F = snes->vec_func; 1290293a7e31SPeter Brune 1291293a7e31SPeter Brune /*norm setup */ 1292fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1293fe6f9142SPeter Brune snes->iter = 0; 1294fe6f9142SPeter Brune snes->norm = 0.; 1295fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1296fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1297fe6f9142SPeter Brune if (snes->domainerror) { 1298fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1299fe6f9142SPeter Brune PetscFunctionReturn(0); 1300fe6f9142SPeter Brune } 1301fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1302fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1303fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1304fe6f9142SPeter Brune snes->norm = fnorm; 1305fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1306fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1307fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1308fe6f9142SPeter Brune 1309fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1310fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1311fe6f9142SPeter Brune /* test convergence */ 1312fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1313fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1314b17ce1afSJed Brown 1315b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1316b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1317b17ce1afSJed Brown DM dmcoarse; 1318b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1319b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1320b17ce1afSJed Brown dm = dmcoarse; 1321b17ce1afSJed Brown } 1322b17ce1afSJed Brown 1323fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1324fe6f9142SPeter Brune /* Call general purpose update function */ 1325646217ecSPeter Brune 1326fe6f9142SPeter Brune if (snes->ops->update) { 1327fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1328fe6f9142SPeter Brune } 132907144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 133039bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 133107144faaSPeter Brune } else { 133207144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 133307144faaSPeter Brune } 1334742fe5e2SPeter Brune 1335742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1336742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1337742fe5e2SPeter Brune PetscFunctionReturn(0); 1338742fe5e2SPeter Brune } 1339c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1340c90fad12SPeter Brune /* Monitor convergence */ 1341c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1342c90fad12SPeter Brune snes->iter = i+1; 1343c90fad12SPeter Brune snes->norm = fnorm; 1344c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1345c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1346c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1347c90fad12SPeter Brune /* Test for convergence */ 1348c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1349c90fad12SPeter Brune if (snes->reason) break; 1350fe6f9142SPeter Brune } 1351fe6f9142SPeter Brune if (i == maxits) { 1352fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1353fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1354fe6f9142SPeter Brune } 1355421d9b32SPeter Brune PetscFunctionReturn(0); 1356421d9b32SPeter Brune } 1357