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 7421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver. 8421d9b32SPeter Brune 9421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 10421d9b32SPeter Brune 11421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 12421d9b32SPeter Brune M*/ 13421d9b32SPeter Brune 14421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 15421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 16421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 18421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 19421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 206273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 21421d9b32SPeter Brune 22421d9b32SPeter Brune EXTERN_C_BEGIN 23ddebd997SPeter Brune #undef __FUNCT__ 24ddebd997SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_FAS" 25ddebd997SPeter Brune PetscErrorCode SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type) 26ddebd997SPeter Brune { 27ddebd997SPeter Brune PetscErrorCode ierr; 28ddebd997SPeter Brune PetscFunctionBegin; 29ddebd997SPeter Brune 30ddebd997SPeter Brune switch (type) { 31ddebd997SPeter Brune case SNES_LS_BASIC: 32ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 33ddebd997SPeter Brune break; 34ddebd997SPeter Brune case SNES_LS_BASIC_NONORMS: 35ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 36ddebd997SPeter Brune break; 37ddebd997SPeter Brune case SNES_LS_QUADRATIC: 38ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 39ddebd997SPeter Brune break; 40ddebd997SPeter Brune default: 41ddebd997SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 42ddebd997SPeter Brune break; 43ddebd997SPeter Brune } 44ddebd997SPeter Brune snes->ls_type = type; 45ddebd997SPeter Brune PetscFunctionReturn(0); 46ddebd997SPeter Brune } 47ddebd997SPeter Brune EXTERN_C_END 48ddebd997SPeter Brune 49ddebd997SPeter Brune EXTERN_C_BEGIN 50421d9b32SPeter Brune 51421d9b32SPeter Brune #undef __FUNCT__ 52421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 53421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 54421d9b32SPeter Brune { 55421d9b32SPeter Brune SNES_FAS * fas; 56421d9b32SPeter Brune PetscErrorCode ierr; 57421d9b32SPeter Brune 58421d9b32SPeter Brune PetscFunctionBegin; 59421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 60421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 61421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 62421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 63421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 64421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 65421d9b32SPeter Brune 66ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 67ed020824SBarry Smith snes->usespc = PETSC_FALSE; 68ed020824SBarry Smith 69421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 70421d9b32SPeter Brune snes->data = (void*) fas; 71421d9b32SPeter Brune fas->level = 0; 72293a7e31SPeter Brune fas->levels = 1; 73ee78dd50SPeter Brune fas->n_cycles = 1; 74ee78dd50SPeter Brune fas->max_up_it = 1; 75ee78dd50SPeter Brune fas->max_down_it = 1; 76ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 77ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 78421d9b32SPeter Brune fas->next = PETSC_NULL; 796273346dSPeter Brune fas->previous = PETSC_NULL; 80421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 81421d9b32SPeter Brune fas->restrct = PETSC_NULL; 82efe1f98aSPeter Brune fas->inject = PETSC_NULL; 83cc05f883SPeter Brune fas->monitor = PETSC_NULL; 84cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 85ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 86ddebd997SPeter Brune 87ddebd997SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 88ddebd997SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 89ddebd997SPeter Brune 90421d9b32SPeter Brune PetscFunctionReturn(0); 91421d9b32SPeter Brune } 92421d9b32SPeter Brune EXTERN_C_END 93421d9b32SPeter Brune 94421d9b32SPeter Brune #undef __FUNCT__ 95421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 96c79ef259SPeter Brune /*@ 97722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 98c79ef259SPeter Brune 99c79ef259SPeter Brune Input Parameter: 100c79ef259SPeter Brune . snes - the preconditioner context 101c79ef259SPeter Brune 102c79ef259SPeter Brune Output parameter: 103c79ef259SPeter Brune . levels - the number of levels 104c79ef259SPeter Brune 105c79ef259SPeter Brune Level: advanced 106c79ef259SPeter Brune 107c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 108c79ef259SPeter Brune 109c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 110c79ef259SPeter Brune @*/ 111421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 112421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 113421d9b32SPeter Brune PetscFunctionBegin; 114ee1fd11aSPeter Brune *levels = fas->levels; 115421d9b32SPeter Brune PetscFunctionReturn(0); 116421d9b32SPeter Brune } 117421d9b32SPeter Brune 118421d9b32SPeter Brune #undef __FUNCT__ 119646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 120c79ef259SPeter Brune /*@ 121c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 122c79ef259SPeter Brune complicated cycling. 123c79ef259SPeter Brune 124c79ef259SPeter Brune Logically Collective on SNES 125c79ef259SPeter Brune 126c79ef259SPeter Brune Input Parameters: 127c79ef259SPeter Brune + snes - the multigrid context 128c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 129c79ef259SPeter Brune 130c79ef259SPeter Brune Options Database Key: 131c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 132c79ef259SPeter Brune 133c79ef259SPeter Brune Level: advanced 134c79ef259SPeter Brune 135c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 136c79ef259SPeter Brune 137c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 138c79ef259SPeter Brune @*/ 139646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 140646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 141646217ecSPeter Brune PetscErrorCode ierr; 142646217ecSPeter Brune PetscFunctionBegin; 143646217ecSPeter Brune fas->n_cycles = cycles; 144646217ecSPeter Brune if (fas->next) { 145646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 146646217ecSPeter Brune } 147646217ecSPeter Brune PetscFunctionReturn(0); 148646217ecSPeter Brune } 149646217ecSPeter Brune 150eff52c0eSPeter Brune #undef __FUNCT__ 151c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 152c79ef259SPeter Brune /*@ 153722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 154c79ef259SPeter Brune 155c79ef259SPeter Brune Logically Collective on SNES 156c79ef259SPeter Brune 157c79ef259SPeter Brune Input Parameters: 158c79ef259SPeter Brune + snes - the multigrid context 159c79ef259SPeter Brune . level - the level to set the number of cycles on 160c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 161c79ef259SPeter Brune 162c79ef259SPeter Brune Level: advanced 163c79ef259SPeter Brune 164c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 165c79ef259SPeter Brune 166c79ef259SPeter Brune .seealso: SNESFASSetCycles() 167c79ef259SPeter Brune @*/ 168c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 169c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 170c79ef259SPeter Brune PetscInt top_level = fas->level,i; 171c79ef259SPeter Brune 172c79ef259SPeter Brune PetscFunctionBegin; 173c79ef259SPeter Brune if (level > top_level) 174c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 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 } 179c79ef259SPeter Brune if (fas->level != level) 180c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 181c79ef259SPeter Brune fas->n_cycles = cycles; 182c79ef259SPeter Brune PetscFunctionReturn(0); 183c79ef259SPeter Brune } 184c79ef259SPeter Brune 185c79ef259SPeter Brune #undef __FUNCT__ 186eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 187aeed3662SMatthew G Knepley /*@C 188c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 189c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 190c79ef259SPeter Brune and nonlinear preconditioners. 191c79ef259SPeter Brune 192c79ef259SPeter Brune Logically Collective on SNES 193c79ef259SPeter Brune 194c79ef259SPeter Brune Input Parameters: 195c79ef259SPeter Brune + snes - the multigrid context 196c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 197c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 198c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 199c79ef259SPeter Brune 200c79ef259SPeter Brune Level: advanced 201c79ef259SPeter Brune 202c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 203c79ef259SPeter Brune 204c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 205c79ef259SPeter Brune @*/ 206eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 207eff52c0eSPeter Brune PetscErrorCode ierr = 0; 208eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 209eff52c0eSPeter Brune PetscFunctionBegin; 210eff52c0eSPeter Brune 211eff52c0eSPeter Brune /* use or don't use it according to user wishes*/ 212d28543b3SPeter Brune snes->usegs = use_gs; 213eff52c0eSPeter Brune if (gsfunc) { 214eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 215eff52c0eSPeter Brune /* push the provided GS up the tree */ 216eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 217eff52c0eSPeter Brune } else if (snes->ops->computegs) { 218eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 219eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 220eff52c0eSPeter Brune } else if (use_gs) { 221eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 222eff52c0eSPeter Brune } 223eff52c0eSPeter Brune PetscFunctionReturn(0); 224eff52c0eSPeter Brune } 225eff52c0eSPeter Brune 226eff52c0eSPeter Brune #undef __FUNCT__ 227eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 228aeed3662SMatthew G Knepley /*@C 229c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 230c79ef259SPeter Brune 231c79ef259SPeter Brune Logically Collective on SNES 232c79ef259SPeter Brune 233c79ef259SPeter Brune Input Parameters: 234c79ef259SPeter Brune + snes - the multigrid context 235c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 236c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 237c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 238c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 239c79ef259SPeter Brune 240c79ef259SPeter Brune Level: advanced 241c79ef259SPeter Brune 242c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 243c79ef259SPeter Brune 244c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 245c79ef259SPeter Brune @*/ 246eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 247eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 248eff52c0eSPeter Brune PetscErrorCode ierr; 249eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 250eff52c0eSPeter Brune SNES cur_snes = snes; 251eff52c0eSPeter Brune PetscFunctionBegin; 252eff52c0eSPeter Brune if (level > top_level) 253eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 254eff52c0eSPeter Brune /* get to the correct level */ 255eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 256eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 257eff52c0eSPeter Brune cur_snes = fas->next; 258eff52c0eSPeter Brune } 259eff52c0eSPeter Brune if (fas->level != level) 260eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 261d28543b3SPeter Brune snes->usegs = use_gs; 262eff52c0eSPeter Brune if (gsfunc) { 2636273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 264eff52c0eSPeter Brune } 265eff52c0eSPeter Brune PetscFunctionReturn(0); 266eff52c0eSPeter Brune } 267eff52c0eSPeter Brune 268646217ecSPeter Brune #undef __FUNCT__ 269421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 270c79ef259SPeter Brune /*@ 271c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 272c79ef259SPeter Brune level of the FAS hierarchy. 273c79ef259SPeter Brune 274c79ef259SPeter Brune Input Parameters: 275c79ef259SPeter Brune + snes - the multigrid context 276c79ef259SPeter Brune level - the level to get 277c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 278c79ef259SPeter Brune 279c79ef259SPeter Brune Level: advanced 280c79ef259SPeter Brune 281c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 282c79ef259SPeter Brune 283c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 284c79ef259SPeter Brune @*/ 285421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 286421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 287421d9b32SPeter Brune PetscInt levels = fas->level; 288421d9b32SPeter Brune PetscInt i; 289421d9b32SPeter Brune PetscFunctionBegin; 290421d9b32SPeter Brune *lsnes = snes; 291421d9b32SPeter Brune if (fas->level < level) { 292421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 293421d9b32SPeter Brune } 294421d9b32SPeter Brune if (level > levels - 1) { 295421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 296421d9b32SPeter Brune } 297421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 298421d9b32SPeter Brune *lsnes = fas->next; 299421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 300421d9b32SPeter Brune } 301421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 302421d9b32SPeter Brune PetscFunctionReturn(0); 303421d9b32SPeter Brune } 304421d9b32SPeter Brune 305421d9b32SPeter Brune #undef __FUNCT__ 30607144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 30707144faaSPeter Brune /*@ 30807144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 30907144faaSPeter Brune e 31007144faaSPeter Brune 31107144faaSPeter Brune 31207144faaSPeter Brune @*/ 31307144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 31407144faaSPeter Brune { 31507144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31607144faaSPeter Brune 31707144faaSPeter Brune PetscFunctionBegin; 31807144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31907144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 32007144faaSPeter Brune fas->fastype = fastype; 32107144faaSPeter Brune PetscFunctionReturn(0); 32207144faaSPeter Brune } 32307144faaSPeter Brune 32407144faaSPeter Brune 32507144faaSPeter Brune 32607144faaSPeter Brune #undef __FUNCT__ 327421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 328aeed3662SMatthew G Knepley /*@C 329c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 330c79ef259SPeter Brune Must be called before any other FAS routine. 331c79ef259SPeter Brune 332c79ef259SPeter Brune Input Parameters: 333c79ef259SPeter Brune + snes - the snes context 334c79ef259SPeter Brune . levels - the number of levels 335c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 336c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 337c79ef259SPeter Brune Fortran. 338c79ef259SPeter Brune 339c79ef259SPeter Brune Level: intermediate 340c79ef259SPeter Brune 341c79ef259SPeter Brune Notes: 342c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 343c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 344c79ef259SPeter Brune 345c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 346c79ef259SPeter Brune 347c79ef259SPeter Brune .seealso: SNESFASGetLevels() 348c79ef259SPeter Brune @*/ 349421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 350421d9b32SPeter Brune PetscErrorCode ierr; 351421d9b32SPeter Brune PetscInt i; 352421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3536273346dSPeter Brune SNES prevsnes; 354421d9b32SPeter Brune MPI_Comm comm; 355421d9b32SPeter Brune PetscFunctionBegin; 356ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 357ee1fd11aSPeter Brune if (levels == fas->levels) { 358ee1fd11aSPeter Brune if (!comms) 359ee1fd11aSPeter Brune PetscFunctionReturn(0); 360ee1fd11aSPeter Brune } 361421d9b32SPeter Brune /* user has changed the number of levels; reset */ 362421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 363421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 364421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 365ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3666273346dSPeter Brune fas->previous = PETSC_NULL; 3676273346dSPeter Brune prevsnes = snes; 368421d9b32SPeter Brune /* setup the finest level */ 369421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 370421d9b32SPeter Brune if (comms) comm = comms[i]; 371421d9b32SPeter Brune fas->level = i; 372421d9b32SPeter Brune fas->levels = levels; 373ee1fd11aSPeter Brune fas->next = PETSC_NULL; 374421d9b32SPeter Brune if (i > 0) { 375421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3764a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 377421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 378794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3796273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3806273346dSPeter Brune prevsnes = fas->next; 3816273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 382421d9b32SPeter Brune } 383421d9b32SPeter Brune } 384421d9b32SPeter Brune PetscFunctionReturn(0); 385421d9b32SPeter Brune } 386421d9b32SPeter Brune 387421d9b32SPeter Brune #undef __FUNCT__ 388c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 389c79ef259SPeter Brune /*@ 390c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 391c79ef259SPeter Brune use on all levels. 392c79ef259SPeter Brune 393c79ef259SPeter Brune Logically Collective on PC 394c79ef259SPeter Brune 395c79ef259SPeter Brune Input Parameters: 396c79ef259SPeter Brune + snes - the multigrid context 397c79ef259SPeter Brune - n - the number of smoothing steps 398c79ef259SPeter Brune 399c79ef259SPeter Brune Options Database Key: 400d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 401c79ef259SPeter Brune 402c79ef259SPeter Brune Level: advanced 403c79ef259SPeter Brune 404c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 405c79ef259SPeter Brune 406c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 407c79ef259SPeter Brune @*/ 408c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 409c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 410c79ef259SPeter Brune PetscErrorCode ierr = 0; 411c79ef259SPeter Brune PetscFunctionBegin; 412c79ef259SPeter Brune fas->max_up_it = n; 413c79ef259SPeter Brune if (fas->next) { 414c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 415c79ef259SPeter Brune } 416c79ef259SPeter Brune PetscFunctionReturn(0); 417c79ef259SPeter Brune } 418c79ef259SPeter Brune 419c79ef259SPeter Brune #undef __FUNCT__ 420c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 421c79ef259SPeter Brune /*@ 422c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 423c79ef259SPeter Brune use on all levels. 424c79ef259SPeter Brune 425c79ef259SPeter Brune Logically Collective on PC 426c79ef259SPeter Brune 427c79ef259SPeter Brune Input Parameters: 428c79ef259SPeter Brune + snes - the multigrid context 429c79ef259SPeter Brune - n - the number of smoothing steps 430c79ef259SPeter Brune 431c79ef259SPeter Brune Options Database Key: 432d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 433c79ef259SPeter Brune 434c79ef259SPeter Brune Level: advanced 435c79ef259SPeter Brune 436c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 437c79ef259SPeter Brune 438c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 439c79ef259SPeter Brune @*/ 440c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 441c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 442c79ef259SPeter Brune PetscErrorCode ierr = 0; 443c79ef259SPeter Brune PetscFunctionBegin; 444c79ef259SPeter Brune fas->max_down_it = n; 445c79ef259SPeter Brune if (fas->next) { 446c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 447c79ef259SPeter Brune } 448c79ef259SPeter Brune PetscFunctionReturn(0); 449c79ef259SPeter Brune } 450c79ef259SPeter Brune 451c79ef259SPeter Brune #undef __FUNCT__ 452421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 453c79ef259SPeter Brune /*@ 454c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 455c79ef259SPeter Brune interpolation from l-1 to the lth level 456c79ef259SPeter Brune 457c79ef259SPeter Brune Input Parameters: 458c79ef259SPeter Brune + snes - the multigrid context 459c79ef259SPeter Brune . mat - the interpolation operator 460c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 461c79ef259SPeter Brune 462c79ef259SPeter Brune Level: advanced 463c79ef259SPeter Brune 464c79ef259SPeter Brune Notes: 465c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 466c79ef259SPeter Brune for the same level. 467c79ef259SPeter Brune 468c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 469c79ef259SPeter Brune out from the matrix size which one it is. 470c79ef259SPeter Brune 471c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 472c79ef259SPeter Brune 473c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 474c79ef259SPeter Brune @*/ 475421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 476421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 477d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 478421d9b32SPeter Brune 479421d9b32SPeter Brune PetscFunctionBegin; 480421d9b32SPeter Brune if (level > top_level) 481421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 482421d9b32SPeter Brune /* get to the correct level */ 483d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 484421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 485421d9b32SPeter Brune } 486421d9b32SPeter Brune if (fas->level != level) 487421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 488421d9b32SPeter Brune fas->interpolate = mat; 489421d9b32SPeter Brune PetscFunctionReturn(0); 490421d9b32SPeter Brune } 491421d9b32SPeter Brune 492421d9b32SPeter Brune #undef __FUNCT__ 493421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 494c79ef259SPeter Brune /*@ 495c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 496c79ef259SPeter Brune from level l to l-1. 497c79ef259SPeter Brune 498c79ef259SPeter Brune Input Parameters: 499c79ef259SPeter Brune + snes - the multigrid context 500c79ef259SPeter Brune . mat - the restriction matrix 501c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 502c79ef259SPeter Brune 503c79ef259SPeter Brune Level: advanced 504c79ef259SPeter Brune 505c79ef259SPeter Brune Notes: 506c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 507c79ef259SPeter Brune for the same level. 508c79ef259SPeter Brune 509c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 510c79ef259SPeter Brune out from the matrix size which one it is. 511c79ef259SPeter Brune 512c79ef259SPeter Brune If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 513c79ef259SPeter Brune is used. 514c79ef259SPeter Brune 515c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 516c79ef259SPeter Brune 517c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 518c79ef259SPeter Brune @*/ 519421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 520421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 521d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 522421d9b32SPeter Brune 523421d9b32SPeter Brune PetscFunctionBegin; 524421d9b32SPeter Brune if (level > top_level) 525421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 526421d9b32SPeter Brune /* get to the correct level */ 527d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 528421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 529421d9b32SPeter Brune } 530421d9b32SPeter Brune if (fas->level != level) 531421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 532421d9b32SPeter Brune fas->restrct = mat; 533421d9b32SPeter Brune PetscFunctionReturn(0); 534421d9b32SPeter Brune } 535421d9b32SPeter Brune 536efe1f98aSPeter Brune 537421d9b32SPeter Brune #undef __FUNCT__ 538421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 539c79ef259SPeter Brune /*@ 540c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 541c79ef259SPeter Brune operator from level l to l-1. 542c79ef259SPeter Brune 543c79ef259SPeter Brune Input Parameters: 544c79ef259SPeter Brune + snes - the multigrid context 545c79ef259SPeter Brune . rscale - the restriction scaling 546c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 547c79ef259SPeter Brune 548c79ef259SPeter Brune Level: advanced 549c79ef259SPeter Brune 550c79ef259SPeter Brune Notes: 551c79ef259SPeter Brune This is only used in the case that the injection is not set. 552c79ef259SPeter Brune 553c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 554c79ef259SPeter Brune 555c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 556c79ef259SPeter Brune @*/ 557421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 558421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 559d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 560421d9b32SPeter Brune 561421d9b32SPeter Brune PetscFunctionBegin; 562421d9b32SPeter Brune if (level > top_level) 563421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 564421d9b32SPeter Brune /* get to the correct level */ 565d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 566421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 567421d9b32SPeter Brune } 568421d9b32SPeter Brune if (fas->level != level) 569421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 570421d9b32SPeter Brune fas->rscale = rscale; 571421d9b32SPeter Brune PetscFunctionReturn(0); 572421d9b32SPeter Brune } 573421d9b32SPeter Brune 574efe1f98aSPeter Brune 575efe1f98aSPeter Brune #undef __FUNCT__ 576efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 577c79ef259SPeter Brune /*@ 578c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 579c79ef259SPeter Brune from level l to l-1. 580c79ef259SPeter Brune 581c79ef259SPeter Brune Input Parameters: 582c79ef259SPeter Brune + snes - the multigrid context 583c79ef259SPeter Brune . mat - the restriction matrix 584c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 585c79ef259SPeter Brune 586c79ef259SPeter Brune Level: advanced 587c79ef259SPeter Brune 588c79ef259SPeter Brune Notes: 589c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 590c79ef259SPeter Brune project the solution instead. 591c79ef259SPeter Brune 592c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 593c79ef259SPeter Brune 594c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 595c79ef259SPeter Brune @*/ 596efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 597efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 598efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 599efe1f98aSPeter Brune 600efe1f98aSPeter Brune PetscFunctionBegin; 601efe1f98aSPeter Brune if (level > top_level) 602efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 603efe1f98aSPeter Brune /* get to the correct level */ 604efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 605efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 606efe1f98aSPeter Brune } 607efe1f98aSPeter Brune if (fas->level != level) 608efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 609efe1f98aSPeter Brune fas->inject = mat; 610efe1f98aSPeter Brune PetscFunctionReturn(0); 611efe1f98aSPeter Brune } 612efe1f98aSPeter Brune 613421d9b32SPeter Brune #undef __FUNCT__ 614421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 615421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 616421d9b32SPeter Brune { 61777df8cc4SPeter Brune PetscErrorCode ierr = 0; 618421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 619421d9b32SPeter Brune 620421d9b32SPeter Brune PetscFunctionBegin; 621742fe5e2SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 622742fe5e2SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 623742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 624421d9b32SPeter Brune PetscFunctionReturn(0); 625421d9b32SPeter Brune } 626421d9b32SPeter Brune 627421d9b32SPeter Brune #undef __FUNCT__ 628421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 629421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 630421d9b32SPeter Brune { 631421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 632742fe5e2SPeter Brune PetscErrorCode ierr = 0; 633421d9b32SPeter Brune 634421d9b32SPeter Brune PetscFunctionBegin; 635421d9b32SPeter Brune /* recursively resets and then destroys */ 63651e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 63751e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 638efe1f98aSPeter Brune if (fas->inject) { 639efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 640efe1f98aSPeter Brune } 64151e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 64251e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 64351e86f29SPeter Brune fas->restrct = PETSC_NULL; 64451e86f29SPeter Brune } else { 64551e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 64651e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 64751e86f29SPeter Brune } 64851e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 649421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 650421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 651ee1fd11aSPeter Brune 652421d9b32SPeter Brune PetscFunctionReturn(0); 653421d9b32SPeter Brune } 654421d9b32SPeter Brune 655421d9b32SPeter Brune #undef __FUNCT__ 656421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 657421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 658421d9b32SPeter Brune { 65948bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 660421d9b32SPeter Brune PetscErrorCode ierr; 661efe1f98aSPeter Brune VecScatter injscatter; 662d1adcc6fSPeter Brune PetscInt dm_levels; 663d1adcc6fSPeter Brune 664421d9b32SPeter Brune PetscFunctionBegin; 665eff52c0eSPeter Brune 666cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 667d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 668d1adcc6fSPeter Brune dm_levels++; 669cc05f883SPeter Brune if (dm_levels > fas->levels) { 670d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 671cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 672d1adcc6fSPeter Brune } 673d1adcc6fSPeter Brune } 674d1adcc6fSPeter Brune 67507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 6762f7ea302SPeter Brune ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 67707144faaSPeter Brune } else { 678ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 67907144faaSPeter Brune } 680cc05f883SPeter Brune 68148bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 68248bfdf8aSPeter Brune if (fas->upsmooth) { 68348bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 68448bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 68548bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 68648bfdf8aSPeter Brune } 68748bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 68848bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 68948bfdf8aSPeter Brune } 69048bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 69148bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 69248bfdf8aSPeter Brune } 69348bfdf8aSPeter Brune } 69448bfdf8aSPeter Brune if (fas->downsmooth) { 69548bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 69648bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 69748bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 69848bfdf8aSPeter Brune } 69948bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 70048bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 70148bfdf8aSPeter Brune } 70248bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 70348bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 70448bfdf8aSPeter Brune } 7051a266240SBarry Smith } 70648bfdf8aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 707646217ecSPeter Brune if (fas->next) { 7086273346dSPeter Brune if (fas->galerkin) { 7096273346dSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 7106273346dSPeter Brune } else { 71148bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 71248bfdf8aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 71348bfdf8aSPeter Brune } 71448bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 71548bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 71648bfdf8aSPeter Brune } 71748bfdf8aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 718646217ecSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 719646217ecSPeter Brune } 720646217ecSPeter Brune } 7216273346dSPeter Brune } 722421d9b32SPeter Brune if (snes->dm) { 723421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 724421d9b32SPeter Brune if (fas->next) { 725421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 726ee78dd50SPeter Brune if (!fas->next->dm) { 727ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 728ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 729ee78dd50SPeter Brune } 730421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 731ee78dd50SPeter Brune if (!fas->interpolate) { 732e727c939SJed Brown ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 733421d9b32SPeter Brune fas->restrct = fas->interpolate; 734421d9b32SPeter Brune } 735efe1f98aSPeter Brune /* set the injection from the DM */ 736efe1f98aSPeter Brune if (!fas->inject) { 737e727c939SJed Brown ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 738efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 739efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 740efe1f98aSPeter Brune } 741421d9b32SPeter Brune } 742ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 743*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 744*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 745421d9b32SPeter Brune } 746cc05f883SPeter Brune 7476273346dSPeter Brune /* setup FAS work vectors */ 7486273346dSPeter Brune if (fas->galerkin) { 7496273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7506273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7516273346dSPeter Brune } 7526273346dSPeter Brune 753fa9694d7SPeter Brune if (fas->next) { 754fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 755*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 756742fe5e2SPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 757fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 758fa9694d7SPeter Brune } 759fa9694d7SPeter Brune /* got to set them all up at once */ 760421d9b32SPeter Brune PetscFunctionReturn(0); 761421d9b32SPeter Brune } 762421d9b32SPeter Brune 763421d9b32SPeter Brune #undef __FUNCT__ 764421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 765421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 766421d9b32SPeter Brune { 767ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 768ee78dd50SPeter Brune PetscInt levels = 1; 769d1adcc6fSPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 770421d9b32SPeter Brune PetscErrorCode ierr; 771ee78dd50SPeter Brune const char *def_smooth = SNESNRICHARDSON; 772ee78dd50SPeter Brune char pre_type[256]; 773ee78dd50SPeter Brune char post_type[256]; 774ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 77507144faaSPeter Brune SNESFASType fastype; 776421d9b32SPeter Brune 777421d9b32SPeter Brune PetscFunctionBegin; 778c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 779ee78dd50SPeter Brune 780ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 781ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 782ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 783c732cbdbSBarry Smith if (!flg && snes->dm) { 784c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 785c732cbdbSBarry Smith levels++; 786d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 787c732cbdbSBarry Smith } 788ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 789ee78dd50SPeter Brune } 790ee78dd50SPeter Brune 791ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 792ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 793ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 79407144faaSPeter Brune fastype = fas->fastype; 79507144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 79607144faaSPeter Brune if (flg) { 79707144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 79807144faaSPeter Brune } 799d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 800d1adcc6fSPeter Brune if (smoothflg) { 801ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 802ee78dd50SPeter Brune } else { 803d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 804d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 805ee78dd50SPeter Brune } 806ee78dd50SPeter Brune 807ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 808d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 809d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 810d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 811ee78dd50SPeter Brune 8126273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8136273346dSPeter Brune 814646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 815ee78dd50SPeter Brune 816ee78dd50SPeter Brune /* other options for the coarsest level */ 817ee78dd50SPeter Brune if (fas->level == 0) { 818d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 819ee78dd50SPeter Brune } 820ee78dd50SPeter Brune 821421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8228cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 823eff52c0eSPeter Brune 824d28543b3SPeter Brune if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) { 8258cc86e31SPeter Brune const char *prefix; 8268cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8278cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8288cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 8298cc86e31SPeter Brune if (fas->level || (fas->levels == 1)) { 830eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 8318cc86e31SPeter Brune } else { 8328cc86e31SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 8338cc86e31SPeter Brune } 8348cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8358cc86e31SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 8368cc86e31SPeter Brune } 8378cc86e31SPeter Brune 838d28543b3SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) { 83967339d5cSBarry Smith const char *prefix; 84067339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 841ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 84267339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 843eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 844293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 845ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 846ee78dd50SPeter Brune } 8471a266240SBarry Smith if (fas->upsmooth) { 8481a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 8491a266240SBarry Smith } 8501a266240SBarry Smith 8511a266240SBarry Smith if (fas->downsmooth) { 8521a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 8531a266240SBarry Smith } 854ee78dd50SPeter Brune 855742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 856742fe5e2SPeter Brune ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr); 857742fe5e2SPeter Brune } 858742fe5e2SPeter Brune 859ee78dd50SPeter Brune if (monflg) { 860646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 861794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8622f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 863742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 864293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 865293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 866d28543b3SPeter Brune } else { 867d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 868d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 869d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 870d28543b3SPeter Brune } 871ee78dd50SPeter Brune } 872ee78dd50SPeter Brune 873ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 874d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 875421d9b32SPeter Brune PetscFunctionReturn(0); 876421d9b32SPeter Brune } 877421d9b32SPeter Brune 878421d9b32SPeter Brune #undef __FUNCT__ 879421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 880421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 881421d9b32SPeter Brune { 882421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 883421d9b32SPeter Brune PetscBool iascii; 884421d9b32SPeter Brune PetscErrorCode ierr; 885421d9b32SPeter Brune PetscInt levels = fas->levels; 886421d9b32SPeter Brune PetscInt i; 887421d9b32SPeter Brune 888421d9b32SPeter Brune PetscFunctionBegin; 889421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 890421d9b32SPeter Brune if (iascii) { 891421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 892421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 893421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 894421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 895ee78dd50SPeter Brune if (fas->upsmooth) { 896ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 897421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 898ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 899421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 900421d9b32SPeter Brune } else { 901ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 902421d9b32SPeter Brune } 903ee78dd50SPeter Brune if (fas->downsmooth) { 904ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 905421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 906ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 907421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 908421d9b32SPeter Brune } else { 909ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 910421d9b32SPeter Brune } 911d28543b3SPeter Brune if (snes->usegs) { 912eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 913eff52c0eSPeter Brune } 914421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 915421d9b32SPeter Brune } 916421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 917421d9b32SPeter Brune } else { 918421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 919421d9b32SPeter Brune } 920421d9b32SPeter Brune PetscFunctionReturn(0); 921421d9b32SPeter Brune } 922421d9b32SPeter Brune 923421d9b32SPeter Brune #undef __FUNCT__ 92439bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92539bd7f45SPeter Brune /* 92639bd7f45SPeter Brune Defines the action of the downsmoother 92739bd7f45SPeter Brune */ 92839bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 92939bd7f45SPeter Brune PetscErrorCode ierr = 0; 9302d15c84fSPeter Brune PetscReal fnorm; 931742fe5e2SPeter Brune SNESConvergedReason reason; 93239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 93339bd7f45SPeter Brune PetscInt k; 934421d9b32SPeter Brune PetscFunctionBegin; 935d1adcc6fSPeter Brune if (fas->downsmooth) { 936d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 937742fe5e2SPeter Brune /* check convergence reason for the smoother */ 938742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 939742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 940742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 941742fe5e2SPeter Brune PetscFunctionReturn(0); 942742fe5e2SPeter Brune } 943d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 944794bee33SPeter Brune if (fas->monitor) { 945794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 946794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 947d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 948eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 949d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 950794bee33SPeter Brune } 95107144faaSPeter Brune for (k = 0; k < fas->max_down_it; k++) { 952646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 953cc05f883SPeter Brune if (fas->monitor) { 954794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 955794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 956d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 957eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 958d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 959794bee33SPeter Brune } 960646217ecSPeter Brune } 961c90fad12SPeter Brune } else if (snes->pc) { 962c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 9632f7ea302SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 9642f7ea302SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 9652f7ea302SPeter Brune snes->reason = SNES_DIVERGED_INNER; 9662f7ea302SPeter Brune PetscFunctionReturn(0); 9672f7ea302SPeter Brune } 968fe6f9142SPeter Brune } 96939bd7f45SPeter Brune PetscFunctionReturn(0); 97039bd7f45SPeter Brune } 97139bd7f45SPeter Brune 97239bd7f45SPeter Brune 97339bd7f45SPeter Brune #undef __FUNCT__ 97439bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 97539bd7f45SPeter Brune /* 97607144faaSPeter Brune Defines the action of the upsmoother 97739bd7f45SPeter Brune */ 97839bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 97939bd7f45SPeter Brune PetscErrorCode ierr = 0; 98039bd7f45SPeter Brune PetscReal fnorm; 98139bd7f45SPeter Brune SNESConvergedReason reason; 98239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 98339bd7f45SPeter Brune PetscInt k; 98439bd7f45SPeter Brune PetscFunctionBegin; 98539bd7f45SPeter Brune if (fas->upsmooth) { 98639bd7f45SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 98739bd7f45SPeter Brune /* check convergence reason for the smoother */ 98839bd7f45SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 98939bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 99039bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 99139bd7f45SPeter Brune PetscFunctionReturn(0); 99239bd7f45SPeter Brune } 99339bd7f45SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 99439bd7f45SPeter Brune if (fas->monitor) { 99539bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 99639bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 99739bd7f45SPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 99839bd7f45SPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 99939bd7f45SPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 100039bd7f45SPeter Brune } 100107144faaSPeter Brune for (k = 0; k < fas->max_up_it; k++) { 100239bd7f45SPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 100339bd7f45SPeter Brune if (fas->monitor) { 100439bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 100539bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 100639bd7f45SPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 100739bd7f45SPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 100839bd7f45SPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 100939bd7f45SPeter Brune } 101039bd7f45SPeter Brune } 101139bd7f45SPeter Brune } else if (snes->pc) { 101239bd7f45SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 101339bd7f45SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 101439bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 101539bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 101639bd7f45SPeter Brune PetscFunctionReturn(0); 101739bd7f45SPeter Brune } 101839bd7f45SPeter Brune } 101939bd7f45SPeter Brune PetscFunctionReturn(0); 102039bd7f45SPeter Brune } 102139bd7f45SPeter Brune 102239bd7f45SPeter Brune #undef __FUNCT__ 102339bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 102439bd7f45SPeter Brune /* 102539bd7f45SPeter Brune 102639bd7f45SPeter Brune Performs the FAS coarse correction as: 102739bd7f45SPeter Brune 102839bd7f45SPeter Brune fine problem: F(x) = 0 102939bd7f45SPeter Brune coarse problem: F^c(x) = b^c 103039bd7f45SPeter Brune 103139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 103239bd7f45SPeter Brune 103339bd7f45SPeter Brune with correction: 103439bd7f45SPeter Brune 103539bd7f45SPeter Brune 103639bd7f45SPeter Brune 103739bd7f45SPeter Brune */ 103839bd7f45SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec B, Vec X, Vec F, Vec X_new) { 103939bd7f45SPeter Brune PetscErrorCode ierr; 104039bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 104139bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 104239bd7f45SPeter Brune SNESConvergedReason reason; 104339bd7f45SPeter Brune PetscFunctionBegin; 1044fa9694d7SPeter Brune if (fas->next) { 1045ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1046794bee33SPeter Brune 1047c90fad12SPeter Brune X_c = fas->next->vec_sol; 1048293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1049c90fad12SPeter Brune F_c = fas->next->vec_func; 1050742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1051efe1f98aSPeter Brune 1052efe1f98aSPeter Brune /* inject the solution */ 1053efe1f98aSPeter Brune if (fas->inject) { 1054a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 1055efe1f98aSPeter Brune } else { 1056a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 1057a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1058efe1f98aSPeter Brune } 1059293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1060293a7e31SPeter Brune 1061293a7e31SPeter Brune /* restrict the defect */ 1062293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1063293a7e31SPeter Brune 1064c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1065ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1066c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1067742fe5e2SPeter Brune 1068293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1069c90fad12SPeter Brune 1070ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1071ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1072c90fad12SPeter Brune 1073c90fad12SPeter Brune /* recurse to the next level */ 1074f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1075742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1076742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1077742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1078742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1079742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1080742fe5e2SPeter Brune PetscFunctionReturn(0); 1081742fe5e2SPeter Brune } 1082742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1083ee78dd50SPeter Brune 1084fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1085fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 108639bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1087293a7e31SPeter Brune } 108839bd7f45SPeter Brune PetscFunctionReturn(0); 108939bd7f45SPeter Brune } 109039bd7f45SPeter Brune 109139bd7f45SPeter Brune #undef __FUNCT__ 109239bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 109339bd7f45SPeter Brune /* 109439bd7f45SPeter Brune 109539bd7f45SPeter Brune The additive cycle looks like: 109639bd7f45SPeter Brune 109707144faaSPeter Brune xhat = x 109807144faaSPeter Brune xhat = dS(x, b) 109907144faaSPeter Brune x = coarsecorrection(xhat, b_d) 110007144faaSPeter Brune x = x + nu*(xhat - x); 110139bd7f45SPeter Brune (optional) x = uS(x, b) 110239bd7f45SPeter Brune 110339bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 110439bd7f45SPeter Brune 110539bd7f45SPeter Brune */ 110639bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 110707144faaSPeter Brune Vec F, B, Xhat; 1108ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 110939bd7f45SPeter Brune PetscErrorCode ierr; 111007144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 111107144faaSPeter Brune SNESConvergedReason reason; 1112ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1113ddebd997SPeter Brune PetscBool lssucceed; 111439bd7f45SPeter Brune PetscFunctionBegin; 111539bd7f45SPeter Brune 111639bd7f45SPeter Brune F = snes->vec_func; 111739bd7f45SPeter Brune B = snes->vec_rhs; 111807144faaSPeter Brune Xhat = snes->work[1]; 1119ddebd997SPeter Brune G = snes->work[2]; 1120ddebd997SPeter Brune W = snes->work[3]; 112107144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 112207144faaSPeter Brune /* recurse first */ 112307144faaSPeter Brune if (fas->next) { 112407144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 112539bd7f45SPeter Brune 112607144faaSPeter Brune X_c = fas->next->vec_sol; 112707144faaSPeter Brune Xo_c = fas->next->work[0]; 112807144faaSPeter Brune F_c = fas->next->vec_func; 112907144faaSPeter Brune B_c = fas->next->vec_rhs; 113039bd7f45SPeter Brune 113107144faaSPeter Brune /* inject the solution */ 113207144faaSPeter Brune if (fas->inject) { 113307144faaSPeter Brune ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr); 113407144faaSPeter Brune } else { 113507144faaSPeter Brune ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr); 113607144faaSPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 113707144faaSPeter Brune } 113807144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 113907144faaSPeter Brune 114007144faaSPeter Brune /* restrict the defect */ 114107144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 114207144faaSPeter Brune 114307144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 114407144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 114507144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 114607144faaSPeter Brune 114707144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 114807144faaSPeter Brune 114907144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 115007144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 115107144faaSPeter Brune 115207144faaSPeter Brune /* recurse */ 115307144faaSPeter Brune fas->next->vec_rhs = B_c; 115407144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 115507144faaSPeter Brune 115607144faaSPeter Brune /* smooth on this level */ 115707144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 115807144faaSPeter Brune 115907144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 116007144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 116107144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 116207144faaSPeter Brune PetscFunctionReturn(0); 116307144faaSPeter Brune } 116407144faaSPeter Brune 116507144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 116607144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1167ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 116807144faaSPeter Brune 1169ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1170ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1171ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1172ddebd997SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1173ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1174ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1175ddebd997SPeter Brune fnorm = gnorm; 117607144faaSPeter Brune } else { 117707144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 117807144faaSPeter Brune } 117939bd7f45SPeter Brune PetscFunctionReturn(0); 118039bd7f45SPeter Brune } 118139bd7f45SPeter Brune 118239bd7f45SPeter Brune #undef __FUNCT__ 118339bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 118439bd7f45SPeter Brune /* 118539bd7f45SPeter Brune 118639bd7f45SPeter Brune Defines the FAS cycle as: 118739bd7f45SPeter Brune 118839bd7f45SPeter Brune fine problem: F(x) = 0 118939bd7f45SPeter Brune coarse problem: F^c(x) = b^c 119039bd7f45SPeter Brune 119139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 119239bd7f45SPeter Brune 119339bd7f45SPeter Brune correction: 119439bd7f45SPeter Brune 119539bd7f45SPeter Brune x = x + I(x^c - Rx) 119639bd7f45SPeter Brune 119739bd7f45SPeter Brune */ 119839bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 119939bd7f45SPeter Brune 120039bd7f45SPeter Brune PetscErrorCode ierr; 120139bd7f45SPeter Brune Vec F,B; 120239bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 120339bd7f45SPeter Brune 120439bd7f45SPeter Brune PetscFunctionBegin; 120539bd7f45SPeter Brune F = snes->vec_func; 120639bd7f45SPeter Brune B = snes->vec_rhs; 120739bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 120839bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 120939bd7f45SPeter Brune 121039bd7f45SPeter Brune ierr = FASCoarseCorrection(snes, B, X, F, X);CHKERRQ(ierr); 121139bd7f45SPeter Brune 1212c90fad12SPeter Brune if (fas->level != 0) { 121339bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1214fe6f9142SPeter Brune } 1215fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1216fa9694d7SPeter Brune 1217fa9694d7SPeter Brune PetscFunctionReturn(0); 1218421d9b32SPeter Brune } 1219421d9b32SPeter Brune 1220421d9b32SPeter Brune #undef __FUNCT__ 1221421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1222421d9b32SPeter Brune 1223421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1224421d9b32SPeter Brune { 1225fa9694d7SPeter Brune PetscErrorCode ierr; 1226fe6f9142SPeter Brune PetscInt i, maxits; 1227f5a6d4f9SBarry Smith Vec X, B, F; 1228fe6f9142SPeter Brune PetscReal fnorm; 122907144faaSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 1230421d9b32SPeter Brune PetscFunctionBegin; 1231fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1232fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1233fa9694d7SPeter Brune X = snes->vec_sol; 1234fe6f9142SPeter Brune B = snes->vec_rhs; 1235f5a6d4f9SBarry Smith F = snes->vec_func; 1236293a7e31SPeter Brune 1237293a7e31SPeter Brune /*norm setup */ 1238fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1239fe6f9142SPeter Brune snes->iter = 0; 1240fe6f9142SPeter Brune snes->norm = 0.; 1241fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1242fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1243fe6f9142SPeter Brune if (snes->domainerror) { 1244fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1245fe6f9142SPeter Brune PetscFunctionReturn(0); 1246fe6f9142SPeter Brune } 1247fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1248fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1249fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1250fe6f9142SPeter Brune snes->norm = fnorm; 1251fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1252fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1253fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1254fe6f9142SPeter Brune 1255fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1256fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1257fe6f9142SPeter Brune /* test convergence */ 1258fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1259fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1260fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1261fe6f9142SPeter Brune /* Call general purpose update function */ 1262646217ecSPeter Brune 1263fe6f9142SPeter Brune if (snes->ops->update) { 1264fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1265fe6f9142SPeter Brune } 126607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 126739bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 126807144faaSPeter Brune } else { 126907144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 127007144faaSPeter Brune } 1271742fe5e2SPeter Brune 1272742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1273742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1274742fe5e2SPeter Brune PetscFunctionReturn(0); 1275742fe5e2SPeter Brune } 1276c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1277c90fad12SPeter Brune /* Monitor convergence */ 1278c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1279c90fad12SPeter Brune snes->iter = i+1; 1280c90fad12SPeter Brune snes->norm = fnorm; 1281c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1282c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1283c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1284c90fad12SPeter Brune /* Test for convergence */ 1285c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1286c90fad12SPeter Brune if (snes->reason) break; 1287fe6f9142SPeter Brune } 1288fe6f9142SPeter Brune if (i == maxits) { 1289fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1290fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1291fe6f9142SPeter Brune } 1292421d9b32SPeter Brune PetscFunctionReturn(0); 1293421d9b32SPeter Brune } 1294