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 70*0e444f03SPeter Brune snes->max_funcs = 30000; 71*0e444f03SPeter Brune snes->max_its = 10000; 72*0e444f03SPeter 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 /*@ 101722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 102c79ef259SPeter Brune 103c79ef259SPeter Brune Input Parameter: 104c79ef259SPeter Brune . snes - the preconditioner 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 /*@ 125c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. 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; 177c79ef259SPeter Brune if (level > top_level) 178c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 179c79ef259SPeter Brune /* get to the correct level */ 180c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 181c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 182c79ef259SPeter Brune } 183c79ef259SPeter Brune if (fas->level != level) 184c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 185c79ef259SPeter Brune fas->n_cycles = cycles; 186c79ef259SPeter Brune PetscFunctionReturn(0); 187c79ef259SPeter Brune } 188c79ef259SPeter Brune 189c79ef259SPeter Brune #undef __FUNCT__ 190eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 191aeed3662SMatthew G Knepley /*@C 192c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 193c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 194c79ef259SPeter Brune and nonlinear preconditioners. 195c79ef259SPeter Brune 196c79ef259SPeter Brune Logically Collective on SNES 197c79ef259SPeter Brune 198c79ef259SPeter Brune Input Parameters: 199c79ef259SPeter Brune + snes - the multigrid context 200c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 201c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 202c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 203c79ef259SPeter Brune 204c79ef259SPeter Brune Level: advanced 205c79ef259SPeter Brune 206c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 207c79ef259SPeter Brune 208c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 209c79ef259SPeter Brune @*/ 210eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 211eff52c0eSPeter Brune PetscErrorCode ierr = 0; 212eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 213eff52c0eSPeter Brune PetscFunctionBegin; 214eff52c0eSPeter Brune 215eff52c0eSPeter Brune if (gsfunc) { 216eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 217eff52c0eSPeter Brune /* push the provided GS up the tree */ 218eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 219eff52c0eSPeter Brune } else if (snes->ops->computegs) { 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) { 223eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 224eff52c0eSPeter Brune } 225eff52c0eSPeter Brune PetscFunctionReturn(0); 226eff52c0eSPeter Brune } 227eff52c0eSPeter Brune 228eff52c0eSPeter Brune #undef __FUNCT__ 229eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 230aeed3662SMatthew G Knepley /*@C 231c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 232c79ef259SPeter Brune 233c79ef259SPeter Brune Logically Collective on SNES 234c79ef259SPeter Brune 235c79ef259SPeter Brune Input Parameters: 236c79ef259SPeter Brune + snes - the multigrid context 237c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 238c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 239c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 240c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 241c79ef259SPeter Brune 242c79ef259SPeter Brune Level: advanced 243c79ef259SPeter Brune 244c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 245c79ef259SPeter Brune 246c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 247c79ef259SPeter Brune @*/ 248eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 249eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 250eff52c0eSPeter Brune PetscErrorCode ierr; 251eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 252eff52c0eSPeter Brune SNES cur_snes = snes; 253eff52c0eSPeter Brune PetscFunctionBegin; 254eff52c0eSPeter Brune if (level > top_level) 255eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 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 } 261eff52c0eSPeter Brune if (fas->level != level) 262eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 263eff52c0eSPeter Brune if (gsfunc) { 2646273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 265eff52c0eSPeter Brune } 266eff52c0eSPeter Brune PetscFunctionReturn(0); 267eff52c0eSPeter Brune } 268eff52c0eSPeter Brune 269646217ecSPeter Brune #undef __FUNCT__ 270421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 271c79ef259SPeter Brune /*@ 272c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 273c79ef259SPeter Brune level of the FAS hierarchy. 274c79ef259SPeter Brune 275c79ef259SPeter Brune Input Parameters: 276c79ef259SPeter Brune + snes - the multigrid context 277c79ef259SPeter Brune level - the level to get 278c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 279c79ef259SPeter Brune 280c79ef259SPeter Brune Level: advanced 281c79ef259SPeter Brune 282c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 283c79ef259SPeter Brune 284c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 285c79ef259SPeter Brune @*/ 286421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 287421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 288421d9b32SPeter Brune PetscInt levels = fas->level; 289421d9b32SPeter Brune PetscInt i; 290421d9b32SPeter Brune PetscFunctionBegin; 291421d9b32SPeter Brune *lsnes = snes; 292421d9b32SPeter Brune if (fas->level < level) { 293421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 294421d9b32SPeter Brune } 295421d9b32SPeter Brune if (level > levels - 1) { 296421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 297421d9b32SPeter Brune } 298421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 299421d9b32SPeter Brune *lsnes = fas->next; 300421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 301421d9b32SPeter Brune } 302421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 303421d9b32SPeter Brune PetscFunctionReturn(0); 304421d9b32SPeter Brune } 305421d9b32SPeter Brune 306421d9b32SPeter Brune #undef __FUNCT__ 30707144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 30807144faaSPeter Brune /*@ 30907144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 31007144faaSPeter Brune e 31107144faaSPeter Brune 31207144faaSPeter Brune 31307144faaSPeter Brune @*/ 31407144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 31507144faaSPeter Brune { 31607144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31707144faaSPeter Brune 31807144faaSPeter Brune PetscFunctionBegin; 31907144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 32007144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 32107144faaSPeter Brune fas->fastype = fastype; 32207144faaSPeter Brune PetscFunctionReturn(0); 32307144faaSPeter Brune } 32407144faaSPeter Brune 32507144faaSPeter Brune 32607144faaSPeter Brune 32707144faaSPeter Brune #undef __FUNCT__ 328421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 329aeed3662SMatthew G Knepley /*@C 330c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 331c79ef259SPeter Brune Must be called before any other FAS routine. 332c79ef259SPeter Brune 333c79ef259SPeter Brune Input Parameters: 334c79ef259SPeter Brune + snes - the snes context 335c79ef259SPeter Brune . levels - the number of levels 336c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 337c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 338c79ef259SPeter Brune Fortran. 339c79ef259SPeter Brune 340c79ef259SPeter Brune Level: intermediate 341c79ef259SPeter Brune 342c79ef259SPeter Brune Notes: 343c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 344c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 345c79ef259SPeter Brune 346c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 347c79ef259SPeter Brune 348c79ef259SPeter Brune .seealso: SNESFASGetLevels() 349c79ef259SPeter Brune @*/ 350421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 351421d9b32SPeter Brune PetscErrorCode ierr; 352421d9b32SPeter Brune PetscInt i; 353421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3546273346dSPeter Brune SNES prevsnes; 355421d9b32SPeter Brune MPI_Comm comm; 356421d9b32SPeter Brune PetscFunctionBegin; 357ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 358ee1fd11aSPeter Brune if (levels == fas->levels) { 359ee1fd11aSPeter Brune if (!comms) 360ee1fd11aSPeter Brune PetscFunctionReturn(0); 361ee1fd11aSPeter Brune } 362421d9b32SPeter Brune /* user has changed the number of levels; reset */ 363421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 364421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 365421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 366ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3676273346dSPeter Brune fas->previous = PETSC_NULL; 3686273346dSPeter Brune prevsnes = snes; 369421d9b32SPeter Brune /* setup the finest level */ 370421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 371421d9b32SPeter Brune if (comms) comm = comms[i]; 372421d9b32SPeter Brune fas->level = i; 373421d9b32SPeter Brune fas->levels = levels; 374ee1fd11aSPeter Brune fas->next = PETSC_NULL; 375421d9b32SPeter Brune if (i > 0) { 376421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3774a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 378421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 379794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3806273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3816273346dSPeter Brune prevsnes = fas->next; 3826273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 383421d9b32SPeter Brune } 384421d9b32SPeter Brune } 385421d9b32SPeter Brune PetscFunctionReturn(0); 386421d9b32SPeter Brune } 387421d9b32SPeter Brune 388421d9b32SPeter Brune #undef __FUNCT__ 389c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 390c79ef259SPeter Brune /*@ 391c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 392c79ef259SPeter Brune use on all levels. 393c79ef259SPeter Brune 394fde0ff24SPeter Brune Logically Collective on SNES 395c79ef259SPeter Brune 396c79ef259SPeter Brune Input Parameters: 397c79ef259SPeter Brune + snes - the multigrid context 398c79ef259SPeter Brune - n - the number of smoothing steps 399c79ef259SPeter Brune 400c79ef259SPeter Brune Options Database Key: 401d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 402c79ef259SPeter Brune 403c79ef259SPeter Brune Level: advanced 404c79ef259SPeter Brune 405c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 406c79ef259SPeter Brune 407c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 408c79ef259SPeter Brune @*/ 409c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 410c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 411c79ef259SPeter Brune PetscErrorCode ierr = 0; 412c79ef259SPeter Brune PetscFunctionBegin; 413c79ef259SPeter Brune fas->max_up_it = n; 414c79ef259SPeter Brune if (fas->next) { 415c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 416c79ef259SPeter Brune } 417c79ef259SPeter Brune PetscFunctionReturn(0); 418c79ef259SPeter Brune } 419c79ef259SPeter Brune 420c79ef259SPeter Brune #undef __FUNCT__ 421c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 422c79ef259SPeter Brune /*@ 423c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 424c79ef259SPeter Brune use on all levels. 425c79ef259SPeter Brune 426fde0ff24SPeter Brune Logically Collective on SNES 427c79ef259SPeter Brune 428c79ef259SPeter Brune Input Parameters: 429c79ef259SPeter Brune + snes - the multigrid context 430c79ef259SPeter Brune - n - the number of smoothing steps 431c79ef259SPeter Brune 432c79ef259SPeter Brune Options Database Key: 433d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 434c79ef259SPeter Brune 435c79ef259SPeter Brune Level: advanced 436c79ef259SPeter Brune 437c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 438c79ef259SPeter Brune 439c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 440c79ef259SPeter Brune @*/ 441c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 442c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 443c79ef259SPeter Brune PetscErrorCode ierr = 0; 444c79ef259SPeter Brune PetscFunctionBegin; 445c79ef259SPeter Brune fas->max_down_it = n; 446c79ef259SPeter Brune if (fas->next) { 447c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 448c79ef259SPeter Brune } 449c79ef259SPeter Brune PetscFunctionReturn(0); 450c79ef259SPeter Brune } 451c79ef259SPeter Brune 452c79ef259SPeter Brune #undef __FUNCT__ 453421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 454c79ef259SPeter Brune /*@ 455c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 456c79ef259SPeter Brune interpolation from l-1 to the lth level 457c79ef259SPeter Brune 458c79ef259SPeter Brune Input Parameters: 459c79ef259SPeter Brune + snes - the multigrid context 460c79ef259SPeter Brune . mat - the interpolation operator 461c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 462c79ef259SPeter Brune 463c79ef259SPeter Brune Level: advanced 464c79ef259SPeter Brune 465c79ef259SPeter Brune Notes: 466c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 467c79ef259SPeter Brune for the same level. 468c79ef259SPeter Brune 469c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 470c79ef259SPeter Brune out from the matrix size which one it is. 471c79ef259SPeter Brune 472c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 473c79ef259SPeter Brune 474bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 475c79ef259SPeter Brune @*/ 476421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 477421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 478d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 479bccf9bb3SJed Brown PetscErrorCode ierr; 480421d9b32SPeter Brune 481421d9b32SPeter Brune PetscFunctionBegin; 482421d9b32SPeter Brune if (level > top_level) 483421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 484421d9b32SPeter Brune /* get to the correct level */ 485d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 486421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 487421d9b32SPeter Brune } 488421d9b32SPeter Brune if (fas->level != level) 489421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 490bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 491421d9b32SPeter Brune fas->interpolate = mat; 492421d9b32SPeter Brune PetscFunctionReturn(0); 493421d9b32SPeter Brune } 494421d9b32SPeter Brune 495421d9b32SPeter Brune #undef __FUNCT__ 496421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 497c79ef259SPeter Brune /*@ 498c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 499c79ef259SPeter Brune from level l to l-1. 500c79ef259SPeter Brune 501c79ef259SPeter Brune Input Parameters: 502c79ef259SPeter Brune + snes - the multigrid context 503c79ef259SPeter Brune . mat - the restriction matrix 504c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 505c79ef259SPeter Brune 506c79ef259SPeter Brune Level: advanced 507c79ef259SPeter Brune 508c79ef259SPeter Brune Notes: 509c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 510c79ef259SPeter Brune for the same level. 511c79ef259SPeter Brune 512c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 513c79ef259SPeter Brune out from the matrix size which one it is. 514c79ef259SPeter Brune 515fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 516c79ef259SPeter Brune is used. 517c79ef259SPeter Brune 518c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 519c79ef259SPeter Brune 520c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 521c79ef259SPeter Brune @*/ 522421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 523421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 524d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 525bccf9bb3SJed Brown PetscErrorCode ierr; 526421d9b32SPeter Brune 527421d9b32SPeter Brune PetscFunctionBegin; 528421d9b32SPeter Brune if (level > top_level) 529421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 530421d9b32SPeter Brune /* get to the correct level */ 531d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 532421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 533421d9b32SPeter Brune } 534421d9b32SPeter Brune if (fas->level != level) 535421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 536bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 537421d9b32SPeter Brune fas->restrct = mat; 538421d9b32SPeter Brune PetscFunctionReturn(0); 539421d9b32SPeter Brune } 540421d9b32SPeter Brune 541efe1f98aSPeter Brune 542421d9b32SPeter Brune #undef __FUNCT__ 543421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 544c79ef259SPeter Brune /*@ 545c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 546c79ef259SPeter Brune operator from level l to l-1. 547c79ef259SPeter Brune 548c79ef259SPeter Brune Input Parameters: 549c79ef259SPeter Brune + snes - the multigrid context 550c79ef259SPeter Brune . rscale - the restriction scaling 551c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 552c79ef259SPeter Brune 553c79ef259SPeter Brune Level: advanced 554c79ef259SPeter Brune 555c79ef259SPeter Brune Notes: 556c79ef259SPeter Brune This is only used in the case that the injection is not set. 557c79ef259SPeter Brune 558c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 559c79ef259SPeter Brune 560c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 561c79ef259SPeter Brune @*/ 562421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 563421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 564d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 565bccf9bb3SJed Brown PetscErrorCode ierr; 566421d9b32SPeter Brune 567421d9b32SPeter Brune PetscFunctionBegin; 568421d9b32SPeter Brune if (level > top_level) 569421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 570421d9b32SPeter Brune /* get to the correct level */ 571d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 572421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 573421d9b32SPeter Brune } 574421d9b32SPeter Brune if (fas->level != level) 575421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 576bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 577421d9b32SPeter Brune fas->rscale = rscale; 578421d9b32SPeter Brune PetscFunctionReturn(0); 579421d9b32SPeter Brune } 580421d9b32SPeter Brune 581efe1f98aSPeter Brune 582efe1f98aSPeter Brune #undef __FUNCT__ 583efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 584c79ef259SPeter Brune /*@ 585c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 586c79ef259SPeter Brune from level l to l-1. 587c79ef259SPeter Brune 588c79ef259SPeter Brune Input Parameters: 589c79ef259SPeter Brune + snes - the multigrid context 590c79ef259SPeter Brune . mat - the restriction matrix 591c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 592c79ef259SPeter Brune 593c79ef259SPeter Brune Level: advanced 594c79ef259SPeter Brune 595c79ef259SPeter Brune Notes: 596c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 597c79ef259SPeter Brune project the solution instead. 598c79ef259SPeter Brune 599c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 600c79ef259SPeter Brune 601c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 602c79ef259SPeter Brune @*/ 603efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 604efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 605efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 606bccf9bb3SJed Brown PetscErrorCode ierr; 607efe1f98aSPeter Brune 608efe1f98aSPeter Brune PetscFunctionBegin; 609efe1f98aSPeter Brune if (level > top_level) 610efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 611efe1f98aSPeter Brune /* get to the correct level */ 612efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 613efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 614efe1f98aSPeter Brune } 615efe1f98aSPeter Brune if (fas->level != level) 616efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 617bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 618efe1f98aSPeter Brune fas->inject = mat; 619efe1f98aSPeter Brune PetscFunctionReturn(0); 620efe1f98aSPeter Brune } 621efe1f98aSPeter Brune 622421d9b32SPeter Brune #undef __FUNCT__ 623421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 624421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 625421d9b32SPeter Brune { 62677df8cc4SPeter Brune PetscErrorCode ierr = 0; 627421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 628421d9b32SPeter Brune 629421d9b32SPeter Brune PetscFunctionBegin; 630bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 631bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 6323dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 633bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 634bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 635bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6363dccd265SPeter Brune 6373dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6383dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 639742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 640421d9b32SPeter Brune PetscFunctionReturn(0); 641421d9b32SPeter Brune } 642421d9b32SPeter Brune 643421d9b32SPeter Brune #undef __FUNCT__ 644421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 645421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 646421d9b32SPeter Brune { 647421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 648742fe5e2SPeter Brune PetscErrorCode ierr = 0; 649421d9b32SPeter Brune 650421d9b32SPeter Brune PetscFunctionBegin; 651421d9b32SPeter Brune /* recursively resets and then destroys */ 65279d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 653421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 654421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 655ee1fd11aSPeter Brune 656421d9b32SPeter Brune PetscFunctionReturn(0); 657421d9b32SPeter Brune } 658421d9b32SPeter Brune 659421d9b32SPeter Brune #undef __FUNCT__ 660421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 661421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 662421d9b32SPeter Brune { 66348bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 664421d9b32SPeter Brune PetscErrorCode ierr; 665efe1f98aSPeter Brune VecScatter injscatter; 666d1adcc6fSPeter Brune PetscInt dm_levels; 6673dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 668d1adcc6fSPeter Brune 669421d9b32SPeter Brune PetscFunctionBegin; 670eff52c0eSPeter Brune 671cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 672d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 673d1adcc6fSPeter Brune dm_levels++; 674cc05f883SPeter Brune if (dm_levels > fas->levels) { 6753dccd265SPeter Brune 6763dccd265SPeter Brune /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/ 6773dccd265SPeter Brune vec_sol = snes->vec_sol; 6783dccd265SPeter Brune vec_func = snes->vec_func; 6793dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6803dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6813dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6823dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6833dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6843dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6853dccd265SPeter Brune 6863dccd265SPeter Brune /* reset the number of levels */ 687d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 688cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6893dccd265SPeter Brune 6903dccd265SPeter Brune snes->vec_sol = vec_sol; 6913dccd265SPeter Brune snes->vec_func = vec_func; 6923dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6933dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 694d1adcc6fSPeter Brune } 695d1adcc6fSPeter Brune } 696d1adcc6fSPeter Brune 6973dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6983dccd265SPeter Brune 69907144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 700fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 70107144faaSPeter Brune } else { 702ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 70307144faaSPeter Brune } 704cc05f883SPeter Brune 70579d9a41aSPeter Brune if (snes->dm) { 70679d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 70779d9a41aSPeter Brune if (fas->next) { 70879d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 70979d9a41aSPeter Brune if (!fas->next->dm) { 71079d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 71179d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 71279d9a41aSPeter Brune } 71379d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 71479d9a41aSPeter Brune if (!fas->interpolate) { 71579d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 716bccf9bb3SJed Brown if (!fas->restrct) { 717bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 71879d9a41aSPeter Brune fas->restrct = fas->interpolate; 71979d9a41aSPeter Brune } 720bccf9bb3SJed Brown } 72179d9a41aSPeter Brune /* set the injection from the DM */ 72279d9a41aSPeter Brune if (!fas->inject) { 72379d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 72479d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 72579d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 72679d9a41aSPeter Brune } 72779d9a41aSPeter Brune } 72879d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 72979d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 73079d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 73179d9a41aSPeter Brune } 73279d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 73379d9a41aSPeter Brune 73479d9a41aSPeter Brune if (fas->next) { 73579d9a41aSPeter Brune if (fas->galerkin) { 73679d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 73779d9a41aSPeter Brune } else { 73879d9a41aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 73979d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 74079d9a41aSPeter Brune } 74179d9a41aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 74279d9a41aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 74379d9a41aSPeter Brune } 74479d9a41aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 74579d9a41aSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 74679d9a41aSPeter Brune } 74779d9a41aSPeter Brune } 74879d9a41aSPeter Brune } 74979d9a41aSPeter Brune 75079d9a41aSPeter Brune if (fas->next) { 75179d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 75279d9a41aSPeter Brune if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 75379d9a41aSPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 75479d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 75579d9a41aSPeter Brune } 75679d9a41aSPeter Brune 75748bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 75848bfdf8aSPeter Brune if (fas->upsmooth) { 75948bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 76048bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 76148bfdf8aSPeter Brune } 76248bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 76348bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 76448bfdf8aSPeter Brune } 76548bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 76648bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 76748bfdf8aSPeter Brune } 768fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 76948bfdf8aSPeter Brune } 77048bfdf8aSPeter Brune if (fas->downsmooth) { 77148bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 77248bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 77348bfdf8aSPeter Brune } 77448bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 77548bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 77648bfdf8aSPeter Brune } 77748bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 77848bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 77948bfdf8aSPeter Brune } 780fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 7811a266240SBarry Smith } 782cc05f883SPeter Brune 7836273346dSPeter Brune /* setup FAS work vectors */ 7846273346dSPeter Brune if (fas->galerkin) { 7856273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7866273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7876273346dSPeter Brune } 7886273346dSPeter Brune 78979d9a41aSPeter Brune 790fa9694d7SPeter Brune /* got to set them all up at once */ 791421d9b32SPeter Brune PetscFunctionReturn(0); 792421d9b32SPeter Brune } 793421d9b32SPeter Brune 794421d9b32SPeter Brune #undef __FUNCT__ 795421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 796421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 797421d9b32SPeter Brune { 798ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 799ee78dd50SPeter Brune PetscInt levels = 1; 800fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 801421d9b32SPeter Brune PetscErrorCode ierr; 802ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 80307144faaSPeter Brune SNESFASType fastype; 804fde0ff24SPeter Brune const char *optionsprefix; 805fde0ff24SPeter Brune const char *prefix; 806421d9b32SPeter Brune 807421d9b32SPeter Brune PetscFunctionBegin; 808c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 809ee78dd50SPeter Brune 810ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 811ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 812ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 813c732cbdbSBarry Smith if (!flg && snes->dm) { 814c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 815c732cbdbSBarry Smith levels++; 816d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 817c732cbdbSBarry Smith } 818ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 819ee78dd50SPeter Brune } 82007144faaSPeter Brune fastype = fas->fastype; 82107144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 82207144faaSPeter Brune if (flg) { 82307144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 82407144faaSPeter Brune } 825ee78dd50SPeter Brune 826fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 827fde0ff24SPeter Brune 828fde0ff24SPeter Brune /* smoother setup options */ 829fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 830fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 831fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 832fde0ff24SPeter Brune if (fas->level == 0) { 833fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 834fde0ff24SPeter Brune } 835ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 836fde0ff24SPeter Brune 837d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 838ee78dd50SPeter Brune 8396273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8406273346dSPeter Brune 841646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 842ee78dd50SPeter Brune 843ee78dd50SPeter Brune 844421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8458cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 846eff52c0eSPeter Brune 847fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 848fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 849fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 850fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 851fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 852fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 853fde0ff24SPeter Brune } 854fde0ff24SPeter Brune 855fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8568cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8578cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8588cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 859fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8608cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8618cc86e31SPeter Brune } 8628cc86e31SPeter Brune 863fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 86467339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 865ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 86667339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 867fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 868293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 869ee78dd50SPeter Brune } 870fde0ff24SPeter Brune 871fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 872fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 873fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 874fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 875fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 876fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 877fde0ff24SPeter Brune } 878fde0ff24SPeter Brune 879fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 880fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 881fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 882fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 883fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 884fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 885fde0ff24SPeter Brune } 886fde0ff24SPeter Brune 8871a266240SBarry Smith if (fas->upsmooth) { 888fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8891a266240SBarry Smith } 8901a266240SBarry Smith 8911a266240SBarry Smith if (fas->downsmooth) { 892fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8931a266240SBarry Smith } 894ee78dd50SPeter Brune 895742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 896fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 897742fe5e2SPeter Brune } 898742fe5e2SPeter Brune 899ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 900ce11c761SPeter Brune if (!fas->downsmooth) { 90179d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 90293dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 903ce11c761SPeter Brune if (fas->level == 0) { 90479d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 905ce11c761SPeter Brune } 906ce11c761SPeter Brune } 907ce11c761SPeter Brune 908ce11c761SPeter Brune if (!fas->upsmooth) { 90979d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 91093dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 911ce11c761SPeter Brune } 912ce11c761SPeter Brune 913ee78dd50SPeter Brune if (monflg) { 914646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 915794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 9162f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 917742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 918293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 919293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 920d28543b3SPeter Brune } else { 921d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 922d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 923d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 924d28543b3SPeter Brune } 925ee78dd50SPeter Brune } 926ee78dd50SPeter Brune 927ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 928d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 929421d9b32SPeter Brune PetscFunctionReturn(0); 930421d9b32SPeter Brune } 931421d9b32SPeter Brune 932421d9b32SPeter Brune #undef __FUNCT__ 933421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 934421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 935421d9b32SPeter Brune { 936421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 937421d9b32SPeter Brune PetscBool iascii; 938421d9b32SPeter Brune PetscErrorCode ierr; 939421d9b32SPeter Brune 940421d9b32SPeter Brune PetscFunctionBegin; 941421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 942421d9b32SPeter Brune if (iascii) { 943421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 944421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 945421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 946ee78dd50SPeter Brune if (fas->upsmooth) { 94739a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 948421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 949ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 950421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 951421d9b32SPeter Brune } else { 952ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 953421d9b32SPeter Brune } 954ee78dd50SPeter Brune if (fas->downsmooth) { 95539a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 956421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 957ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 958421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 959421d9b32SPeter Brune } else { 960ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 961421d9b32SPeter Brune } 962421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 963421d9b32SPeter Brune } else { 964421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 965421d9b32SPeter Brune } 966421d9b32SPeter Brune PetscFunctionReturn(0); 967421d9b32SPeter Brune } 968421d9b32SPeter Brune 969421d9b32SPeter Brune #undef __FUNCT__ 97039bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 97139bd7f45SPeter Brune /* 97239bd7f45SPeter Brune Defines the action of the downsmoother 97339bd7f45SPeter Brune */ 97439bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 97539bd7f45SPeter Brune PetscErrorCode ierr = 0; 976fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 977742fe5e2SPeter Brune SNESConvergedReason reason; 97839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 9794b32a720SPeter Brune Vec G, W, Y, FPC; 980fde0ff24SPeter Brune PetscBool lssuccess; 981fde0ff24SPeter Brune PetscInt k; 982421d9b32SPeter Brune PetscFunctionBegin; 983fde0ff24SPeter Brune G = snes->work[1]; 984fde0ff24SPeter Brune W = snes->work[2]; 985fde0ff24SPeter Brune Y = snes->work[3]; 986d1adcc6fSPeter Brune if (fas->downsmooth) { 987d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 988742fe5e2SPeter Brune /* check convergence reason for the smoother */ 989742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 990742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 991742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 992742fe5e2SPeter Brune PetscFunctionReturn(0); 993742fe5e2SPeter Brune } 9944b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9954b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 996fde0ff24SPeter Brune } else { 997fde0ff24SPeter Brune /* basic richardson smoother */ 998fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 999794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1000794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1001fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1002fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1003fde0ff24SPeter Brune if (!lssuccess) { 1004fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1005fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 10062f7ea302SPeter Brune PetscFunctionReturn(0); 10072f7ea302SPeter Brune } 1008fe6f9142SPeter Brune } 1009fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1010fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1011fde0ff24SPeter Brune fnorm = gnorm; 1012fde0ff24SPeter Brune } 1013fde0ff24SPeter Brune } 101439bd7f45SPeter Brune PetscFunctionReturn(0); 101539bd7f45SPeter Brune } 101639bd7f45SPeter Brune 101739bd7f45SPeter Brune 101839bd7f45SPeter Brune #undef __FUNCT__ 101939bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 102039bd7f45SPeter Brune /* 102107144faaSPeter Brune Defines the action of the upsmoother 102239bd7f45SPeter Brune */ 102339bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 102439bd7f45SPeter Brune PetscErrorCode ierr = 0; 1025fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 102639bd7f45SPeter Brune SNESConvergedReason reason; 102739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 10284b32a720SPeter Brune Vec G, W, Y, FPC; 1029fde0ff24SPeter Brune PetscBool lssuccess; 1030fde0ff24SPeter Brune PetscInt k; 103139bd7f45SPeter Brune PetscFunctionBegin; 1032fde0ff24SPeter Brune G = snes->work[1]; 1033fde0ff24SPeter Brune W = snes->work[2]; 1034fde0ff24SPeter Brune Y = snes->work[3]; 103539bd7f45SPeter Brune if (fas->upsmooth) { 1036fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 103739bd7f45SPeter Brune /* check convergence reason for the smoother */ 1038fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 103939bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 104039bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 104139bd7f45SPeter Brune PetscFunctionReturn(0); 104239bd7f45SPeter Brune } 10434b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 10444b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 1045fde0ff24SPeter Brune } else { 1046fde0ff24SPeter Brune /* basic richardson smoother */ 1047fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 104839bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 104939bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1050fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 1051fde0ff24SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1052fde0ff24SPeter Brune if (!lssuccess) { 1053fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1054fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 105539bd7f45SPeter Brune PetscFunctionReturn(0); 105639bd7f45SPeter Brune } 105739bd7f45SPeter Brune } 1058fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1059fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1060fde0ff24SPeter Brune fnorm = gnorm; 1061fde0ff24SPeter Brune } 1062fde0ff24SPeter Brune } 106339bd7f45SPeter Brune PetscFunctionReturn(0); 106439bd7f45SPeter Brune } 106539bd7f45SPeter Brune 106639bd7f45SPeter Brune #undef __FUNCT__ 106739bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 106839bd7f45SPeter Brune /* 106939bd7f45SPeter Brune 107039bd7f45SPeter Brune Performs the FAS coarse correction as: 107139bd7f45SPeter Brune 107239bd7f45SPeter Brune fine problem: F(x) = 0 107339bd7f45SPeter Brune coarse problem: F^c(x) = b^c 107439bd7f45SPeter Brune 107539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 107639bd7f45SPeter Brune 107739bd7f45SPeter Brune with correction: 107839bd7f45SPeter Brune 107939bd7f45SPeter Brune 108039bd7f45SPeter Brune 108139bd7f45SPeter Brune */ 108239a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 108339bd7f45SPeter Brune PetscErrorCode ierr; 108439bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 108539bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 108639bd7f45SPeter Brune SNESConvergedReason reason; 108739bd7f45SPeter Brune PetscFunctionBegin; 1088fa9694d7SPeter Brune if (fas->next) { 1089c90fad12SPeter Brune X_c = fas->next->vec_sol; 1090293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1091c90fad12SPeter Brune F_c = fas->next->vec_func; 1092742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1093efe1f98aSPeter Brune 1094efe1f98aSPeter Brune /* inject the solution */ 1095efe1f98aSPeter Brune if (fas->inject) { 1096a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 1097efe1f98aSPeter Brune } else { 1098a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 1099a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1100efe1f98aSPeter Brune } 1101293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1102293a7e31SPeter Brune 1103293a7e31SPeter Brune /* restrict the defect */ 1104293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1105293a7e31SPeter Brune 1106c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1107ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1108c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1109742fe5e2SPeter Brune 1110293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1111c90fad12SPeter Brune 1112ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1113ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1114c90fad12SPeter Brune 1115c90fad12SPeter Brune /* recurse to the next level */ 1116f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1117742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1118742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1119742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1120742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1121742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1122742fe5e2SPeter Brune PetscFunctionReturn(0); 1123742fe5e2SPeter Brune } 1124742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1125ee78dd50SPeter Brune 1126fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1127fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 112839bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1129293a7e31SPeter Brune } 113039bd7f45SPeter Brune PetscFunctionReturn(0); 113139bd7f45SPeter Brune } 113239bd7f45SPeter Brune 113339bd7f45SPeter Brune #undef __FUNCT__ 113439bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 113539bd7f45SPeter Brune /* 113639bd7f45SPeter Brune 113739bd7f45SPeter Brune The additive cycle looks like: 113839bd7f45SPeter Brune 113907144faaSPeter Brune xhat = x 114007144faaSPeter Brune xhat = dS(x, b) 114107144faaSPeter Brune x = coarsecorrection(xhat, b_d) 114207144faaSPeter Brune x = x + nu*(xhat - x); 114339bd7f45SPeter Brune (optional) x = uS(x, b) 114439bd7f45SPeter Brune 114539bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 114639bd7f45SPeter Brune 114739bd7f45SPeter Brune */ 114839bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 114907144faaSPeter Brune Vec F, B, Xhat; 1150ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 115139bd7f45SPeter Brune PetscErrorCode ierr; 115207144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 115307144faaSPeter Brune SNESConvergedReason reason; 1154ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1155ddebd997SPeter Brune PetscBool lssucceed; 115639bd7f45SPeter Brune PetscFunctionBegin; 115739bd7f45SPeter Brune 115839bd7f45SPeter Brune F = snes->vec_func; 115939bd7f45SPeter Brune B = snes->vec_rhs; 1160fde0ff24SPeter Brune Xhat = snes->work[3]; 1161fde0ff24SPeter Brune G = snes->work[1]; 1162fde0ff24SPeter Brune W = snes->work[2]; 116307144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 116407144faaSPeter Brune /* recurse first */ 116507144faaSPeter Brune if (fas->next) { 116607144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 116739bd7f45SPeter Brune 116807144faaSPeter Brune X_c = fas->next->vec_sol; 116907144faaSPeter Brune Xo_c = fas->next->work[0]; 117007144faaSPeter Brune F_c = fas->next->vec_func; 117107144faaSPeter Brune B_c = fas->next->vec_rhs; 117239bd7f45SPeter Brune 117307144faaSPeter Brune /* inject the solution */ 117407144faaSPeter Brune if (fas->inject) { 117507144faaSPeter Brune ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr); 117607144faaSPeter Brune } else { 117707144faaSPeter Brune ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr); 117807144faaSPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 117907144faaSPeter Brune } 118007144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 118107144faaSPeter Brune 118207144faaSPeter Brune /* restrict the defect */ 118307144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 118407144faaSPeter Brune 118507144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 118607144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 118707144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 118807144faaSPeter Brune 118907144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 119007144faaSPeter Brune 119107144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 119207144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 119307144faaSPeter Brune 119407144faaSPeter Brune /* recurse */ 119507144faaSPeter Brune fas->next->vec_rhs = B_c; 119607144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 119707144faaSPeter Brune 119807144faaSPeter Brune /* smooth on this level */ 119907144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 120007144faaSPeter Brune 120107144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 120207144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 120307144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 120407144faaSPeter Brune PetscFunctionReturn(0); 120507144faaSPeter Brune } 120607144faaSPeter Brune 120707144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 120807144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1209ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 121007144faaSPeter Brune 1211ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1212ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1213ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1214eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 1215ddebd997SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1216ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1217ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1218ddebd997SPeter Brune fnorm = gnorm; 121907144faaSPeter Brune } else { 122007144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 122107144faaSPeter Brune } 122239bd7f45SPeter Brune PetscFunctionReturn(0); 122339bd7f45SPeter Brune } 122439bd7f45SPeter Brune 122539bd7f45SPeter Brune #undef __FUNCT__ 122639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 122739bd7f45SPeter Brune /* 122839bd7f45SPeter Brune 122939bd7f45SPeter Brune Defines the FAS cycle as: 123039bd7f45SPeter Brune 123139bd7f45SPeter Brune fine problem: F(x) = 0 123239bd7f45SPeter Brune coarse problem: F^c(x) = b^c 123339bd7f45SPeter Brune 123439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 123539bd7f45SPeter Brune 123639bd7f45SPeter Brune correction: 123739bd7f45SPeter Brune 123839bd7f45SPeter Brune x = x + I(x^c - Rx) 123939bd7f45SPeter Brune 124039bd7f45SPeter Brune */ 124139bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 124239bd7f45SPeter Brune 124339bd7f45SPeter Brune PetscErrorCode ierr; 124439bd7f45SPeter Brune Vec F,B; 124539bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 124639bd7f45SPeter Brune 124739bd7f45SPeter Brune PetscFunctionBegin; 124839bd7f45SPeter Brune F = snes->vec_func; 124939bd7f45SPeter Brune B = snes->vec_rhs; 125039bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 125139bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 125239bd7f45SPeter Brune 125339a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 125439bd7f45SPeter Brune 1255c90fad12SPeter Brune if (fas->level != 0) { 125639bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1257fe6f9142SPeter Brune } 1258fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1259fa9694d7SPeter Brune 1260fa9694d7SPeter Brune PetscFunctionReturn(0); 1261421d9b32SPeter Brune } 1262421d9b32SPeter Brune 1263421d9b32SPeter Brune #undef __FUNCT__ 1264421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1265421d9b32SPeter Brune 1266421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1267421d9b32SPeter Brune { 1268fa9694d7SPeter Brune PetscErrorCode ierr; 1269fe6f9142SPeter Brune PetscInt i, maxits; 1270ddb5aff1SPeter Brune Vec X, F; 1271fe6f9142SPeter Brune PetscReal fnorm; 127207144faaSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 1273421d9b32SPeter Brune PetscFunctionBegin; 1274fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1275fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1276fa9694d7SPeter Brune X = snes->vec_sol; 1277f5a6d4f9SBarry Smith F = snes->vec_func; 1278293a7e31SPeter Brune 1279293a7e31SPeter Brune /*norm setup */ 1280fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1281fe6f9142SPeter Brune snes->iter = 0; 1282fe6f9142SPeter Brune snes->norm = 0.; 1283fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1284fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1285fe6f9142SPeter Brune if (snes->domainerror) { 1286fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1287fe6f9142SPeter Brune PetscFunctionReturn(0); 1288fe6f9142SPeter Brune } 1289fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1290fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1291fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1292fe6f9142SPeter Brune snes->norm = fnorm; 1293fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1294fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1295fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1296fe6f9142SPeter Brune 1297fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1298fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1299fe6f9142SPeter Brune /* test convergence */ 1300fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1301fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1302fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1303fe6f9142SPeter Brune /* Call general purpose update function */ 1304646217ecSPeter Brune 1305fe6f9142SPeter Brune if (snes->ops->update) { 1306fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1307fe6f9142SPeter Brune } 130807144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 130939bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 131007144faaSPeter Brune } else { 131107144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 131207144faaSPeter Brune } 1313742fe5e2SPeter Brune 1314742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1315742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1316742fe5e2SPeter Brune PetscFunctionReturn(0); 1317742fe5e2SPeter Brune } 1318c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1319c90fad12SPeter Brune /* Monitor convergence */ 1320c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1321c90fad12SPeter Brune snes->iter = i+1; 1322c90fad12SPeter Brune snes->norm = fnorm; 1323c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1324c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1325c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1326c90fad12SPeter Brune /* Test for convergence */ 1327c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1328c90fad12SPeter Brune if (snes->reason) break; 1329fe6f9142SPeter Brune } 1330fe6f9142SPeter Brune if (i == maxits) { 1331fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1332fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1333fe6f9142SPeter Brune } 1334421d9b32SPeter Brune PetscFunctionReturn(0); 1335421d9b32SPeter Brune } 1336