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 /*@ 982e8ce248SJed Brown SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 99c79ef259SPeter Brune 100c79ef259SPeter Brune Input Parameter: 1012e8ce248SJed Brown . snes - the nonlinear solver 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 /*@ 1222e8ce248SJed Brown SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. 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; 1742e8ce248SJed 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); 175c79ef259SPeter Brune /* get to the correct level */ 176c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 177c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 178c79ef259SPeter Brune } 1792e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 180c79ef259SPeter Brune fas->n_cycles = cycles; 181c79ef259SPeter Brune PetscFunctionReturn(0); 182c79ef259SPeter Brune } 183c79ef259SPeter Brune 184c79ef259SPeter Brune #undef __FUNCT__ 185eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 186aeed3662SMatthew G Knepley /*@C 187c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 188c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 189c79ef259SPeter Brune and nonlinear preconditioners. 190c79ef259SPeter Brune 191c79ef259SPeter Brune Logically Collective on SNES 192c79ef259SPeter Brune 193c79ef259SPeter Brune Input Parameters: 194c79ef259SPeter Brune + snes - the multigrid context 195c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 196c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 197c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 198c79ef259SPeter Brune 199c79ef259SPeter Brune Level: advanced 200c79ef259SPeter Brune 201c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 202c79ef259SPeter Brune 203c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 204c79ef259SPeter Brune @*/ 205eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 206eff52c0eSPeter Brune PetscErrorCode ierr = 0; 207eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 208eff52c0eSPeter Brune PetscFunctionBegin; 209eff52c0eSPeter Brune 210eff52c0eSPeter Brune if (gsfunc) { 211eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 212eff52c0eSPeter Brune /* push the provided GS up the tree */ 213eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 214eff52c0eSPeter Brune } else if (snes->ops->computegs) { 215eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 216eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 217eff52c0eSPeter Brune } else if (use_gs) { 2182e8ce248SJed Brown SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level); 219eff52c0eSPeter Brune } 220eff52c0eSPeter Brune PetscFunctionReturn(0); 221eff52c0eSPeter Brune } 222eff52c0eSPeter Brune 223eff52c0eSPeter Brune #undef __FUNCT__ 224eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 225aeed3662SMatthew G Knepley /*@C 226c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 227c79ef259SPeter Brune 228c79ef259SPeter Brune Logically Collective on SNES 229c79ef259SPeter Brune 230c79ef259SPeter Brune Input Parameters: 231c79ef259SPeter Brune + snes - the multigrid context 232c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 233c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 234c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 235c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 236c79ef259SPeter Brune 237c79ef259SPeter Brune Level: advanced 238c79ef259SPeter Brune 239c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 240c79ef259SPeter Brune 241c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 242c79ef259SPeter Brune @*/ 243eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 244eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 245eff52c0eSPeter Brune PetscErrorCode ierr; 246eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 247eff52c0eSPeter Brune SNES cur_snes = snes; 248eff52c0eSPeter Brune PetscFunctionBegin; 2492e8ce248SJed 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); 250eff52c0eSPeter Brune /* get to the correct level */ 251eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 252eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 253eff52c0eSPeter Brune cur_snes = fas->next; 254eff52c0eSPeter Brune } 2552e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 256eff52c0eSPeter Brune if (gsfunc) { 2576273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 258eff52c0eSPeter Brune } 259eff52c0eSPeter Brune PetscFunctionReturn(0); 260eff52c0eSPeter Brune } 261eff52c0eSPeter Brune 262646217ecSPeter Brune #undef __FUNCT__ 263421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 264c79ef259SPeter Brune /*@ 265c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 266c79ef259SPeter Brune level of the FAS hierarchy. 267c79ef259SPeter Brune 268c79ef259SPeter Brune Input Parameters: 269c79ef259SPeter Brune + snes - the multigrid context 270c79ef259SPeter Brune level - the level to get 271c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 272c79ef259SPeter Brune 273c79ef259SPeter Brune Level: advanced 274c79ef259SPeter Brune 275c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 276c79ef259SPeter Brune 277c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 278c79ef259SPeter Brune @*/ 279421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) { 280421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 281421d9b32SPeter Brune PetscInt i; 2822e8ce248SJed Brown 283421d9b32SPeter Brune PetscFunctionBegin; 2842e8ce248SJed 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); 2852e8ce248SJed 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); 286421d9b32SPeter Brune *lsnes = snes; 287421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 288421d9b32SPeter Brune *lsnes = fas->next; 289421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 290421d9b32SPeter Brune } 2912e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 292421d9b32SPeter Brune PetscFunctionReturn(0); 293421d9b32SPeter Brune } 294421d9b32SPeter Brune 295421d9b32SPeter Brune #undef __FUNCT__ 29607144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 29707144faaSPeter Brune /*@ 29807144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 29907144faaSPeter Brune e 30007144faaSPeter Brune 30107144faaSPeter Brune 30207144faaSPeter Brune @*/ 30307144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 30407144faaSPeter Brune { 30507144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 30607144faaSPeter Brune 30707144faaSPeter Brune PetscFunctionBegin; 30807144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 30907144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 31007144faaSPeter Brune fas->fastype = fastype; 31107144faaSPeter Brune PetscFunctionReturn(0); 31207144faaSPeter Brune } 31307144faaSPeter Brune 31407144faaSPeter Brune #undef __FUNCT__ 315421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 316aeed3662SMatthew G Knepley /*@C 317c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 318c79ef259SPeter Brune Must be called before any other FAS routine. 319c79ef259SPeter Brune 320c79ef259SPeter Brune Input Parameters: 321c79ef259SPeter Brune + snes - the snes context 322c79ef259SPeter Brune . levels - the number of levels 323c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 324c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 325c79ef259SPeter Brune Fortran. 326c79ef259SPeter Brune 327c79ef259SPeter Brune Level: intermediate 328c79ef259SPeter Brune 329c79ef259SPeter Brune Notes: 330c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 331c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 332c79ef259SPeter Brune 333c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 334c79ef259SPeter Brune 335c79ef259SPeter Brune .seealso: SNESFASGetLevels() 336c79ef259SPeter Brune @*/ 337421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 338421d9b32SPeter Brune PetscErrorCode ierr; 339421d9b32SPeter Brune PetscInt i; 340421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3416273346dSPeter Brune SNES prevsnes; 342421d9b32SPeter Brune MPI_Comm comm; 343421d9b32SPeter Brune PetscFunctionBegin; 344ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 345ee1fd11aSPeter Brune if (levels == fas->levels) { 346ee1fd11aSPeter Brune if (!comms) 347ee1fd11aSPeter Brune PetscFunctionReturn(0); 348ee1fd11aSPeter Brune } 349421d9b32SPeter Brune /* user has changed the number of levels; reset */ 350421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 351421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 352421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 353ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3546273346dSPeter Brune fas->previous = PETSC_NULL; 3556273346dSPeter Brune prevsnes = snes; 356421d9b32SPeter Brune /* setup the finest level */ 357421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 358421d9b32SPeter Brune if (comms) comm = comms[i]; 359421d9b32SPeter Brune fas->level = i; 360421d9b32SPeter Brune fas->levels = levels; 361ee1fd11aSPeter Brune fas->next = PETSC_NULL; 362421d9b32SPeter Brune if (i > 0) { 363421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3644a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 365421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 366794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3676273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3686273346dSPeter Brune prevsnes = fas->next; 3696273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 370421d9b32SPeter Brune } 371421d9b32SPeter Brune } 372421d9b32SPeter Brune PetscFunctionReturn(0); 373421d9b32SPeter Brune } 374421d9b32SPeter Brune 375421d9b32SPeter Brune #undef __FUNCT__ 376c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 377c79ef259SPeter Brune /*@ 378c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 379c79ef259SPeter Brune use on all levels. 380c79ef259SPeter Brune 381fde0ff24SPeter Brune Logically Collective on SNES 382c79ef259SPeter Brune 383c79ef259SPeter Brune Input Parameters: 384c79ef259SPeter Brune + snes - the multigrid context 385c79ef259SPeter Brune - n - the number of smoothing steps 386c79ef259SPeter Brune 387c79ef259SPeter Brune Options Database Key: 388d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 389c79ef259SPeter Brune 390c79ef259SPeter Brune Level: advanced 391c79ef259SPeter Brune 392c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 393c79ef259SPeter Brune 394c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 395c79ef259SPeter Brune @*/ 396c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 397c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 398c79ef259SPeter Brune PetscErrorCode ierr = 0; 399c79ef259SPeter Brune PetscFunctionBegin; 400c79ef259SPeter Brune fas->max_up_it = n; 401c79ef259SPeter Brune if (fas->next) { 402c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 403c79ef259SPeter Brune } 404c79ef259SPeter Brune PetscFunctionReturn(0); 405c79ef259SPeter Brune } 406c79ef259SPeter Brune 407c79ef259SPeter Brune #undef __FUNCT__ 408c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 409c79ef259SPeter Brune /*@ 410c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 411c79ef259SPeter Brune use on all levels. 412c79ef259SPeter Brune 413fde0ff24SPeter Brune Logically Collective on SNES 414c79ef259SPeter Brune 415c79ef259SPeter Brune Input Parameters: 416c79ef259SPeter Brune + snes - the multigrid context 417c79ef259SPeter Brune - n - the number of smoothing steps 418c79ef259SPeter Brune 419c79ef259SPeter Brune Options Database Key: 420d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 421c79ef259SPeter Brune 422c79ef259SPeter Brune Level: advanced 423c79ef259SPeter Brune 424c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 425c79ef259SPeter Brune 426c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 427c79ef259SPeter Brune @*/ 428c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 429c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 430c79ef259SPeter Brune PetscErrorCode ierr = 0; 431c79ef259SPeter Brune PetscFunctionBegin; 432c79ef259SPeter Brune fas->max_down_it = n; 433c79ef259SPeter Brune if (fas->next) { 434c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 435c79ef259SPeter Brune } 436c79ef259SPeter Brune PetscFunctionReturn(0); 437c79ef259SPeter Brune } 438c79ef259SPeter Brune 439c79ef259SPeter Brune #undef __FUNCT__ 440421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 441c79ef259SPeter Brune /*@ 442c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 443c79ef259SPeter Brune interpolation from l-1 to the lth level 444c79ef259SPeter Brune 445c79ef259SPeter Brune Input Parameters: 446c79ef259SPeter Brune + snes - the multigrid context 447c79ef259SPeter Brune . mat - the interpolation operator 448c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 449c79ef259SPeter Brune 450c79ef259SPeter Brune Level: advanced 451c79ef259SPeter Brune 452c79ef259SPeter Brune Notes: 453c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 454c79ef259SPeter Brune for the same level. 455c79ef259SPeter Brune 456c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 457c79ef259SPeter Brune out from the matrix size which one it is. 458c79ef259SPeter Brune 459c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 460c79ef259SPeter Brune 461bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 462c79ef259SPeter Brune @*/ 463421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 464421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 465d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 466bccf9bb3SJed Brown PetscErrorCode ierr; 467421d9b32SPeter Brune 468421d9b32SPeter Brune PetscFunctionBegin; 4692e8ce248SJed 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); 470421d9b32SPeter Brune /* get to the correct level */ 471d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 472421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 473421d9b32SPeter Brune } 4742e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 475bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 476bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 477421d9b32SPeter Brune fas->interpolate = mat; 478421d9b32SPeter Brune PetscFunctionReturn(0); 479421d9b32SPeter Brune } 480421d9b32SPeter Brune 481421d9b32SPeter Brune #undef __FUNCT__ 482421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 483c79ef259SPeter Brune /*@ 484c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 485c79ef259SPeter Brune from level l to l-1. 486c79ef259SPeter Brune 487c79ef259SPeter Brune Input Parameters: 488c79ef259SPeter Brune + snes - the multigrid context 489c79ef259SPeter Brune . mat - the restriction matrix 490c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 491c79ef259SPeter Brune 492c79ef259SPeter Brune Level: advanced 493c79ef259SPeter Brune 494c79ef259SPeter Brune Notes: 495c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 496c79ef259SPeter Brune for the same level. 497c79ef259SPeter Brune 498c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 499c79ef259SPeter Brune out from the matrix size which one it is. 500c79ef259SPeter Brune 501fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 502c79ef259SPeter Brune is used. 503c79ef259SPeter Brune 504c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 505c79ef259SPeter Brune 506c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 507c79ef259SPeter Brune @*/ 508421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 509421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 510d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 511bccf9bb3SJed Brown PetscErrorCode ierr; 512421d9b32SPeter Brune 513421d9b32SPeter Brune PetscFunctionBegin; 5142e8ce248SJed 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); 515421d9b32SPeter Brune /* get to the correct level */ 516d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 517421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 518421d9b32SPeter Brune } 5192e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 520bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 521bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 522421d9b32SPeter Brune fas->restrct = mat; 523421d9b32SPeter Brune PetscFunctionReturn(0); 524421d9b32SPeter Brune } 525421d9b32SPeter Brune 526421d9b32SPeter Brune #undef __FUNCT__ 527421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 528c79ef259SPeter Brune /*@ 529c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 530c79ef259SPeter Brune operator from level l to l-1. 531c79ef259SPeter Brune 532c79ef259SPeter Brune Input Parameters: 533c79ef259SPeter Brune + snes - the multigrid context 534c79ef259SPeter Brune . rscale - the restriction scaling 535c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 536c79ef259SPeter Brune 537c79ef259SPeter Brune Level: advanced 538c79ef259SPeter Brune 539c79ef259SPeter Brune Notes: 540c79ef259SPeter Brune This is only used in the case that the injection is not set. 541c79ef259SPeter Brune 542c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 543c79ef259SPeter Brune 544c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 545c79ef259SPeter Brune @*/ 546421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 547421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 548d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 549bccf9bb3SJed Brown PetscErrorCode ierr; 550421d9b32SPeter Brune 551421d9b32SPeter Brune PetscFunctionBegin; 5522e8ce248SJed 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); 553421d9b32SPeter Brune /* get to the correct level */ 554d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 555421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 556421d9b32SPeter Brune } 5572e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 558bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 559bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 560421d9b32SPeter Brune fas->rscale = rscale; 561421d9b32SPeter Brune PetscFunctionReturn(0); 562421d9b32SPeter Brune } 563421d9b32SPeter Brune 564efe1f98aSPeter Brune 565efe1f98aSPeter Brune #undef __FUNCT__ 566efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 567c79ef259SPeter Brune /*@ 568c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 569c79ef259SPeter Brune from level l to l-1. 570c79ef259SPeter Brune 571c79ef259SPeter Brune Input Parameters: 572c79ef259SPeter Brune + snes - the multigrid context 573c79ef259SPeter Brune . mat - the restriction matrix 574c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 575c79ef259SPeter Brune 576c79ef259SPeter Brune Level: advanced 577c79ef259SPeter Brune 578c79ef259SPeter Brune Notes: 579c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 580c79ef259SPeter Brune project the solution instead. 581c79ef259SPeter Brune 582c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 583c79ef259SPeter Brune 584c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 585c79ef259SPeter Brune @*/ 586efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 587efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 588efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 589bccf9bb3SJed Brown PetscErrorCode ierr; 590efe1f98aSPeter Brune 591efe1f98aSPeter Brune PetscFunctionBegin; 5922e8ce248SJed 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); 593efe1f98aSPeter Brune /* get to the correct level */ 594efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 595efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 596efe1f98aSPeter Brune } 5972e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 598bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 599bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 600efe1f98aSPeter Brune fas->inject = mat; 601efe1f98aSPeter Brune PetscFunctionReturn(0); 602efe1f98aSPeter Brune } 603efe1f98aSPeter Brune 604421d9b32SPeter Brune #undef __FUNCT__ 605421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 606421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 607421d9b32SPeter Brune { 60877df8cc4SPeter Brune PetscErrorCode ierr = 0; 609421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 610421d9b32SPeter Brune 611421d9b32SPeter Brune PetscFunctionBegin; 612bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 613bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 6143dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 615bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 616bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 617bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6183dccd265SPeter Brune 6193dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6203dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 621742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 62209dc8347SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","",PETSC_NULL);CHKERRQ(ierr); 623421d9b32SPeter Brune PetscFunctionReturn(0); 624421d9b32SPeter Brune } 625421d9b32SPeter Brune 626421d9b32SPeter Brune #undef __FUNCT__ 627421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 628421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 629421d9b32SPeter Brune { 630421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 631742fe5e2SPeter Brune PetscErrorCode ierr = 0; 632421d9b32SPeter Brune 633421d9b32SPeter Brune PetscFunctionBegin; 634421d9b32SPeter Brune /* recursively resets and then destroys */ 63579d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 636421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 637421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 638ee1fd11aSPeter Brune 639421d9b32SPeter Brune PetscFunctionReturn(0); 640421d9b32SPeter Brune } 641421d9b32SPeter Brune 642421d9b32SPeter Brune #undef __FUNCT__ 643421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 644421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 645421d9b32SPeter Brune { 64648bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 647421d9b32SPeter Brune PetscErrorCode ierr; 648efe1f98aSPeter Brune VecScatter injscatter; 649d1adcc6fSPeter Brune PetscInt dm_levels; 6503dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 651d1adcc6fSPeter Brune 652421d9b32SPeter Brune PetscFunctionBegin; 653eff52c0eSPeter Brune 654cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 655d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 656d1adcc6fSPeter Brune dm_levels++; 657cc05f883SPeter Brune if (dm_levels > fas->levels) { 6583dccd265SPeter Brune 6592e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6603dccd265SPeter Brune vec_sol = snes->vec_sol; 6613dccd265SPeter Brune vec_func = snes->vec_func; 6623dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6633dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6643dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6653dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6663dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6673dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6683dccd265SPeter Brune 6693dccd265SPeter Brune /* reset the number of levels */ 670d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 671cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6723dccd265SPeter Brune 6733dccd265SPeter Brune snes->vec_sol = vec_sol; 6743dccd265SPeter Brune snes->vec_func = vec_func; 6753dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6763dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 677d1adcc6fSPeter Brune } 678d1adcc6fSPeter Brune } 679d1adcc6fSPeter Brune 6803dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6813dccd265SPeter Brune 68207144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 683fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 68407144faaSPeter Brune } else { 685ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 68607144faaSPeter Brune } 687cc05f883SPeter Brune 68879d9a41aSPeter Brune if (snes->dm) { 68979d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 69079d9a41aSPeter Brune if (fas->next) { 69179d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 69279d9a41aSPeter Brune if (!fas->next->dm) { 69379d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 69479d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 69579d9a41aSPeter Brune } 69679d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 69779d9a41aSPeter Brune if (!fas->interpolate) { 69879d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 699bccf9bb3SJed Brown if (!fas->restrct) { 700bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 70179d9a41aSPeter Brune fas->restrct = fas->interpolate; 70279d9a41aSPeter Brune } 703bccf9bb3SJed Brown } 70479d9a41aSPeter Brune /* set the injection from the DM */ 70579d9a41aSPeter Brune if (!fas->inject) { 70679d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 70779d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 70879d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 70979d9a41aSPeter Brune } 71079d9a41aSPeter Brune } 71179d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 71279d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 71379d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 71479d9a41aSPeter Brune } 71579d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 71679d9a41aSPeter Brune 71779d9a41aSPeter Brune if (fas->next) { 71879d9a41aSPeter Brune if (fas->galerkin) { 71979d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 72079d9a41aSPeter Brune } else { 72179d9a41aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 72279d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 72379d9a41aSPeter Brune } 72479d9a41aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 72579d9a41aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 72679d9a41aSPeter Brune } 72779d9a41aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 72879d9a41aSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 72979d9a41aSPeter Brune } 73079d9a41aSPeter Brune } 73179d9a41aSPeter Brune } 73279d9a41aSPeter Brune 73379d9a41aSPeter Brune if (fas->next) { 73479d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 735938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 736938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 73779d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 73879d9a41aSPeter Brune } 73979d9a41aSPeter Brune 74048bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 74148bfdf8aSPeter Brune if (fas->upsmooth) { 74248bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 74348bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 74448bfdf8aSPeter Brune } 74548bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 74648bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 74748bfdf8aSPeter Brune } 74848bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 74948bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 75048bfdf8aSPeter Brune } 751fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 75248bfdf8aSPeter Brune } 75348bfdf8aSPeter Brune if (fas->downsmooth) { 75448bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 75548bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 75648bfdf8aSPeter Brune } 75748bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 75848bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 75948bfdf8aSPeter Brune } 76048bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 76148bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 76248bfdf8aSPeter Brune } 763fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 7641a266240SBarry Smith } 765cc05f883SPeter Brune 7666273346dSPeter Brune /* setup FAS work vectors */ 7676273346dSPeter Brune if (fas->galerkin) { 7686273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7696273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7706273346dSPeter Brune } 7716273346dSPeter Brune 77279d9a41aSPeter Brune 773fa9694d7SPeter Brune /* got to set them all up at once */ 774421d9b32SPeter Brune PetscFunctionReturn(0); 775421d9b32SPeter Brune } 776421d9b32SPeter Brune 777421d9b32SPeter Brune #undef __FUNCT__ 778421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 779421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 780421d9b32SPeter Brune { 781ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 782ee78dd50SPeter Brune PetscInt levels = 1; 783fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 784421d9b32SPeter Brune PetscErrorCode ierr; 785ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 78607144faaSPeter Brune SNESFASType fastype; 787fde0ff24SPeter Brune const char *optionsprefix; 788fde0ff24SPeter Brune const char *prefix; 789421d9b32SPeter Brune 790421d9b32SPeter Brune PetscFunctionBegin; 791c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 792ee78dd50SPeter Brune 793ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 794ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 795ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 796c732cbdbSBarry Smith if (!flg && snes->dm) { 797c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 798c732cbdbSBarry Smith levels++; 799d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 800c732cbdbSBarry Smith } 801ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 802ee78dd50SPeter Brune } 80307144faaSPeter Brune fastype = fas->fastype; 80407144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 80507144faaSPeter Brune if (flg) { 80607144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 80707144faaSPeter Brune } 808ee78dd50SPeter Brune 809fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 810fde0ff24SPeter Brune 811fde0ff24SPeter Brune /* smoother setup options */ 812fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 813fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 814fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 815fde0ff24SPeter Brune if (fas->level == 0) { 816fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 817fde0ff24SPeter Brune } 818ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 819fde0ff24SPeter Brune 820d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 821ee78dd50SPeter Brune 8226273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8236273346dSPeter Brune 824646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 825ee78dd50SPeter Brune 826ee78dd50SPeter Brune 827421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8288cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 829eff52c0eSPeter Brune 830fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 831fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 832fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 833fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 834fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 835fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 836fde0ff24SPeter Brune } 837fde0ff24SPeter Brune 838fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8398cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8408cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8418cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 842fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8438cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8448cc86e31SPeter Brune } 8458cc86e31SPeter Brune 846fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 84767339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 848ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 84967339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 850fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 851293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 852ee78dd50SPeter Brune } 853fde0ff24SPeter Brune 854fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 855fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 856fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 857fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 858fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 859fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 860fde0ff24SPeter Brune } 861fde0ff24SPeter Brune 862fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 863fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 864fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 865fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 866fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 867fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 868fde0ff24SPeter Brune } 869fde0ff24SPeter Brune 8701a266240SBarry Smith if (fas->upsmooth) { 871fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8721a266240SBarry Smith } 8731a266240SBarry Smith 8741a266240SBarry Smith if (fas->downsmooth) { 875fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8761a266240SBarry Smith } 877ee78dd50SPeter Brune 878742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 879fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 880742fe5e2SPeter Brune } 881742fe5e2SPeter Brune 882ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 883ce11c761SPeter Brune if (!fas->downsmooth) { 88479d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 88593dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 886ce11c761SPeter Brune if (fas->level == 0) { 88779d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 888ce11c761SPeter Brune } 889ce11c761SPeter Brune } 890ce11c761SPeter Brune 891ce11c761SPeter Brune if (!fas->upsmooth) { 89279d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 89393dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 894ce11c761SPeter Brune } 895ce11c761SPeter Brune 896ee78dd50SPeter Brune if (monflg) { 897646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 898794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8992f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 900742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 901293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 902293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 903d28543b3SPeter Brune } else { 904d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 905d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 906d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 907d28543b3SPeter Brune } 908ee78dd50SPeter Brune } 909ee78dd50SPeter Brune 910ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 911d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 912421d9b32SPeter Brune PetscFunctionReturn(0); 913421d9b32SPeter Brune } 914421d9b32SPeter Brune 915421d9b32SPeter Brune #undef __FUNCT__ 916421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 917421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 918421d9b32SPeter Brune { 919421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 920421d9b32SPeter Brune PetscBool iascii; 921421d9b32SPeter Brune PetscErrorCode ierr; 922421d9b32SPeter Brune 923421d9b32SPeter Brune PetscFunctionBegin; 924421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 925421d9b32SPeter Brune if (iascii) { 9262e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 927421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 928bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 929ee78dd50SPeter Brune if (fas->upsmooth) { 93039a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 931421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 932ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 933421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 934421d9b32SPeter Brune } else { 935ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 936421d9b32SPeter Brune } 937ee78dd50SPeter Brune if (fas->downsmooth) { 93839a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 939421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 940ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 941421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 942421d9b32SPeter Brune } else { 943ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 944421d9b32SPeter Brune } 945421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 946421d9b32SPeter Brune } else { 947421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 948421d9b32SPeter Brune } 949421d9b32SPeter Brune PetscFunctionReturn(0); 950421d9b32SPeter Brune } 951421d9b32SPeter Brune 952421d9b32SPeter Brune #undef __FUNCT__ 95339bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 95439bd7f45SPeter Brune /* 95539bd7f45SPeter Brune Defines the action of the downsmoother 95639bd7f45SPeter Brune */ 95739bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 95839bd7f45SPeter Brune PetscErrorCode ierr = 0; 959fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 960742fe5e2SPeter Brune SNESConvergedReason reason; 96139bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 9624b32a720SPeter Brune Vec G, W, Y, FPC; 963fde0ff24SPeter Brune PetscBool lssuccess; 964fde0ff24SPeter Brune PetscInt k; 965421d9b32SPeter Brune PetscFunctionBegin; 966fde0ff24SPeter Brune G = snes->work[1]; 967fde0ff24SPeter Brune W = snes->work[2]; 968fde0ff24SPeter Brune Y = snes->work[3]; 969d1adcc6fSPeter Brune if (fas->downsmooth) { 970d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 971742fe5e2SPeter Brune /* check convergence reason for the smoother */ 972742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 973742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 974742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 975742fe5e2SPeter Brune PetscFunctionReturn(0); 976742fe5e2SPeter Brune } 9774b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9784b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 979fde0ff24SPeter Brune } else { 980fde0ff24SPeter Brune /* basic richardson smoother */ 981fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 982794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 983794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 984fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 985fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 986fde0ff24SPeter Brune if (!lssuccess) { 987fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 988fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9892f7ea302SPeter Brune PetscFunctionReturn(0); 9902f7ea302SPeter Brune } 991fe6f9142SPeter Brune } 992fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 993fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 994fde0ff24SPeter Brune fnorm = gnorm; 995fde0ff24SPeter Brune } 996fde0ff24SPeter Brune } 99739bd7f45SPeter Brune PetscFunctionReturn(0); 99839bd7f45SPeter Brune } 99939bd7f45SPeter Brune 100039bd7f45SPeter Brune 100139bd7f45SPeter Brune #undef __FUNCT__ 100239bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 100339bd7f45SPeter Brune /* 100407144faaSPeter Brune Defines the action of the upsmoother 100539bd7f45SPeter Brune */ 100639bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 100739bd7f45SPeter Brune PetscErrorCode ierr = 0; 1008fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 100939bd7f45SPeter Brune SNESConvergedReason reason; 101039bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 10114b32a720SPeter Brune Vec G, W, Y, FPC; 1012fde0ff24SPeter Brune PetscBool lssuccess; 1013fde0ff24SPeter Brune PetscInt k; 101439bd7f45SPeter Brune PetscFunctionBegin; 1015fde0ff24SPeter Brune G = snes->work[1]; 1016fde0ff24SPeter Brune W = snes->work[2]; 1017fde0ff24SPeter Brune Y = snes->work[3]; 101839bd7f45SPeter Brune if (fas->upsmooth) { 1019fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 102039bd7f45SPeter Brune /* check convergence reason for the smoother */ 1021fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 102239bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 102339bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 102439bd7f45SPeter Brune PetscFunctionReturn(0); 102539bd7f45SPeter Brune } 10264b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 10274b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 1028fde0ff24SPeter Brune } else { 1029fde0ff24SPeter Brune /* basic richardson smoother */ 1030fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 103139bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 103239bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1033fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1034fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1035fde0ff24SPeter Brune if (!lssuccess) { 1036fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1037fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 103839bd7f45SPeter Brune PetscFunctionReturn(0); 103939bd7f45SPeter Brune } 104039bd7f45SPeter Brune } 1041fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1042fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1043fde0ff24SPeter Brune fnorm = gnorm; 1044fde0ff24SPeter Brune } 1045fde0ff24SPeter Brune } 104639bd7f45SPeter Brune PetscFunctionReturn(0); 104739bd7f45SPeter Brune } 104839bd7f45SPeter Brune 104939bd7f45SPeter Brune #undef __FUNCT__ 1050938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1051938e4a01SJed Brown /*@ 1052938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1053938e4a01SJed Brown 1054938e4a01SJed Brown Collective 1055938e4a01SJed Brown 1056938e4a01SJed Brown Input Arguments: 1057938e4a01SJed Brown . snes - SNESFAS 1058938e4a01SJed Brown 1059938e4a01SJed Brown Output Arguments: 1060938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1061938e4a01SJed Brown 1062938e4a01SJed Brown Level: developer 1063938e4a01SJed Brown 1064938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1065938e4a01SJed Brown @*/ 1066938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1067938e4a01SJed Brown { 1068938e4a01SJed Brown PetscErrorCode ierr; 1069938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1070938e4a01SJed Brown 1071938e4a01SJed Brown PetscFunctionBegin; 1072938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1073938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1074938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1075938e4a01SJed Brown PetscFunctionReturn(0); 1076938e4a01SJed Brown } 1077938e4a01SJed Brown 1078e9923e8dSJed Brown #undef __FUNCT__ 1079e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1080e9923e8dSJed Brown /*@ 1081e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1082e9923e8dSJed Brown 1083e9923e8dSJed Brown Collective 1084e9923e8dSJed Brown 1085e9923e8dSJed Brown Input Arguments: 1086e9923e8dSJed Brown + fine - SNES from which to restrict 1087e9923e8dSJed Brown - Xfine - vector to restrict 1088e9923e8dSJed Brown 1089e9923e8dSJed Brown Output Arguments: 1090e9923e8dSJed Brown . Xcoarse - result of restriction 1091e9923e8dSJed Brown 1092e9923e8dSJed Brown Level: developer 1093e9923e8dSJed Brown 1094e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1095e9923e8dSJed Brown @*/ 1096e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1097e9923e8dSJed Brown { 1098e9923e8dSJed Brown PetscErrorCode ierr; 1099e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1100e9923e8dSJed Brown 1101e9923e8dSJed Brown PetscFunctionBegin; 1102e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1103e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1104e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1105e9923e8dSJed Brown if (fas->inject) { 1106e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1107e9923e8dSJed Brown } else { 1108e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1109e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1110e9923e8dSJed Brown } 1111e9923e8dSJed Brown PetscFunctionReturn(0); 1112e9923e8dSJed Brown } 1113e9923e8dSJed Brown 1114e9923e8dSJed Brown #undef __FUNCT__ 111539bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 111639bd7f45SPeter Brune /* 111739bd7f45SPeter Brune 111839bd7f45SPeter Brune Performs the FAS coarse correction as: 111939bd7f45SPeter Brune 112039bd7f45SPeter Brune fine problem: F(x) = 0 112139bd7f45SPeter Brune coarse problem: F^c(x) = b^c 112239bd7f45SPeter Brune 112339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 112439bd7f45SPeter Brune 112539bd7f45SPeter Brune with correction: 112639bd7f45SPeter Brune 112739bd7f45SPeter Brune 112839bd7f45SPeter Brune 112939bd7f45SPeter Brune */ 113039a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 113139bd7f45SPeter Brune PetscErrorCode ierr; 113239bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 113339bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 113439bd7f45SPeter Brune SNESConvergedReason reason; 113539bd7f45SPeter Brune PetscFunctionBegin; 1136fa9694d7SPeter Brune if (fas->next) { 1137c90fad12SPeter Brune X_c = fas->next->vec_sol; 1138293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1139c90fad12SPeter Brune F_c = fas->next->vec_func; 1140742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1141efe1f98aSPeter Brune 1142938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1143293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1144293a7e31SPeter Brune 1145293a7e31SPeter Brune /* restrict the defect */ 1146293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1147293a7e31SPeter Brune 1148c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1149ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1150c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1151742fe5e2SPeter Brune 1152293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1153c90fad12SPeter Brune 1154ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1155ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1156c90fad12SPeter Brune 1157c90fad12SPeter Brune /* recurse to the next level */ 1158f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1159742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1160742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1161742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1162742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1163742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1164742fe5e2SPeter Brune PetscFunctionReturn(0); 1165742fe5e2SPeter Brune } 1166742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1167ee78dd50SPeter Brune 1168fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1169fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 117039bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1171293a7e31SPeter Brune } 117239bd7f45SPeter Brune PetscFunctionReturn(0); 117339bd7f45SPeter Brune } 117439bd7f45SPeter Brune 117539bd7f45SPeter Brune #undef __FUNCT__ 117639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 117739bd7f45SPeter Brune /* 117839bd7f45SPeter Brune 117939bd7f45SPeter Brune The additive cycle looks like: 118039bd7f45SPeter Brune 118107144faaSPeter Brune xhat = x 118207144faaSPeter Brune xhat = dS(x, b) 118307144faaSPeter Brune x = coarsecorrection(xhat, b_d) 118407144faaSPeter Brune x = x + nu*(xhat - x); 118539bd7f45SPeter Brune (optional) x = uS(x, b) 118639bd7f45SPeter Brune 118739bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 118839bd7f45SPeter Brune 118939bd7f45SPeter Brune */ 119039bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 119107144faaSPeter Brune Vec F, B, Xhat; 1192ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 119339bd7f45SPeter Brune PetscErrorCode ierr; 119407144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 119507144faaSPeter Brune SNESConvergedReason reason; 1196ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1197ddebd997SPeter Brune PetscBool lssucceed; 119839bd7f45SPeter Brune PetscFunctionBegin; 119939bd7f45SPeter Brune 120039bd7f45SPeter Brune F = snes->vec_func; 120139bd7f45SPeter Brune B = snes->vec_rhs; 1202fde0ff24SPeter Brune Xhat = snes->work[3]; 1203fde0ff24SPeter Brune G = snes->work[1]; 1204fde0ff24SPeter Brune W = snes->work[2]; 120507144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 120607144faaSPeter Brune /* recurse first */ 120707144faaSPeter Brune if (fas->next) { 120807144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 120939bd7f45SPeter Brune 121007144faaSPeter Brune X_c = fas->next->vec_sol; 121107144faaSPeter Brune Xo_c = fas->next->work[0]; 121207144faaSPeter Brune F_c = fas->next->vec_func; 121307144faaSPeter Brune B_c = fas->next->vec_rhs; 121439bd7f45SPeter Brune 1215938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 121607144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 121707144faaSPeter Brune 121807144faaSPeter Brune /* restrict the defect */ 121907144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 122007144faaSPeter Brune 122107144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 122207144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 122307144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 122407144faaSPeter Brune 122507144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 122607144faaSPeter Brune 122707144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 122807144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 122907144faaSPeter Brune 123007144faaSPeter Brune /* recurse */ 123107144faaSPeter Brune fas->next->vec_rhs = B_c; 123207144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 123307144faaSPeter Brune 123407144faaSPeter Brune /* smooth on this level */ 123507144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 123607144faaSPeter Brune 123707144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 123807144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 123907144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 124007144faaSPeter Brune PetscFunctionReturn(0); 124107144faaSPeter Brune } 124207144faaSPeter Brune 124307144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 124407144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1245ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 124607144faaSPeter Brune 1247ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1248ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1249ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1250eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 1251ddebd997SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1252ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1253ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1254ddebd997SPeter Brune fnorm = gnorm; 125507144faaSPeter Brune } else { 125607144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 125707144faaSPeter Brune } 125839bd7f45SPeter Brune PetscFunctionReturn(0); 125939bd7f45SPeter Brune } 126039bd7f45SPeter Brune 126139bd7f45SPeter Brune #undef __FUNCT__ 126239bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 126339bd7f45SPeter Brune /* 126439bd7f45SPeter Brune 126539bd7f45SPeter Brune Defines the FAS cycle as: 126639bd7f45SPeter Brune 126739bd7f45SPeter Brune fine problem: F(x) = 0 126839bd7f45SPeter Brune coarse problem: F^c(x) = b^c 126939bd7f45SPeter Brune 127039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 127139bd7f45SPeter Brune 127239bd7f45SPeter Brune correction: 127339bd7f45SPeter Brune 127439bd7f45SPeter Brune x = x + I(x^c - Rx) 127539bd7f45SPeter Brune 127639bd7f45SPeter Brune */ 127739bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 127839bd7f45SPeter Brune 127939bd7f45SPeter Brune PetscErrorCode ierr; 128039bd7f45SPeter Brune Vec F,B; 128139bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 128239bd7f45SPeter Brune 128339bd7f45SPeter Brune PetscFunctionBegin; 128439bd7f45SPeter Brune F = snes->vec_func; 128539bd7f45SPeter Brune B = snes->vec_rhs; 128639bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 128739bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 128839bd7f45SPeter Brune 128939a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 129039bd7f45SPeter Brune 1291c90fad12SPeter Brune if (fas->level != 0) { 129239bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1293fe6f9142SPeter Brune } 1294fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1295fa9694d7SPeter Brune 1296fa9694d7SPeter Brune PetscFunctionReturn(0); 1297421d9b32SPeter Brune } 1298421d9b32SPeter Brune 1299421d9b32SPeter Brune #undef __FUNCT__ 1300421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1301421d9b32SPeter Brune 1302421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1303421d9b32SPeter Brune { 1304fa9694d7SPeter Brune PetscErrorCode ierr; 1305fe6f9142SPeter Brune PetscInt i, maxits; 1306ddb5aff1SPeter Brune Vec X, F; 1307fe6f9142SPeter Brune PetscReal fnorm; 1308*b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1309*b17ce1afSJed Brown DM dm; 1310*b17ce1afSJed Brown 1311421d9b32SPeter Brune PetscFunctionBegin; 1312fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1313fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1314fa9694d7SPeter Brune X = snes->vec_sol; 1315f5a6d4f9SBarry Smith F = snes->vec_func; 1316293a7e31SPeter Brune 1317293a7e31SPeter Brune /*norm setup */ 1318fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1319fe6f9142SPeter Brune snes->iter = 0; 1320fe6f9142SPeter Brune snes->norm = 0.; 1321fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1322fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1323fe6f9142SPeter Brune if (snes->domainerror) { 1324fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1325fe6f9142SPeter Brune PetscFunctionReturn(0); 1326fe6f9142SPeter Brune } 1327fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1328fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1329fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1330fe6f9142SPeter Brune snes->norm = fnorm; 1331fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1332fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1333fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1334fe6f9142SPeter Brune 1335fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1336fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1337fe6f9142SPeter Brune /* test convergence */ 1338fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1339fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1340*b17ce1afSJed Brown 1341*b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1342*b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1343*b17ce1afSJed Brown DM dmcoarse; 1344*b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1345*b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1346*b17ce1afSJed Brown dm = dmcoarse; 1347*b17ce1afSJed Brown } 1348*b17ce1afSJed Brown 1349fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1350fe6f9142SPeter Brune /* Call general purpose update function */ 1351646217ecSPeter Brune 1352fe6f9142SPeter Brune if (snes->ops->update) { 1353fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1354fe6f9142SPeter Brune } 135507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 135639bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 135707144faaSPeter Brune } else { 135807144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 135907144faaSPeter Brune } 1360742fe5e2SPeter Brune 1361742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1362742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1363742fe5e2SPeter Brune PetscFunctionReturn(0); 1364742fe5e2SPeter Brune } 1365c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1366c90fad12SPeter Brune /* Monitor convergence */ 1367c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1368c90fad12SPeter Brune snes->iter = i+1; 1369c90fad12SPeter Brune snes->norm = fnorm; 1370c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1371c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1372c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1373c90fad12SPeter Brune /* Test for convergence */ 1374c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1375c90fad12SPeter Brune if (snes->reason) break; 1376fe6f9142SPeter Brune } 1377fe6f9142SPeter Brune if (i == maxits) { 1378fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1379fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1380fe6f9142SPeter Brune } 1381421d9b32SPeter Brune PetscFunctionReturn(0); 1382421d9b32SPeter Brune } 1383