1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 507144faaSPeter Brune 6421d9b32SPeter Brune /*MC 7ef536925SPeter Brune 8ef536925SPeter Brune SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 9421d9b32SPeter Brune 10421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 11421d9b32SPeter Brune 12421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 13421d9b32SPeter Brune M*/ 14421d9b32SPeter Brune 15421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 16421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 18421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 19421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 20421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 216273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 22421d9b32SPeter Brune 23421d9b32SPeter Brune EXTERN_C_BEGIN 24ddebd997SPeter Brune #undef __FUNCT__ 25ddebd997SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_FAS" 26ddebd997SPeter Brune PetscErrorCode SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type) 27ddebd997SPeter Brune { 28ddebd997SPeter Brune PetscErrorCode ierr; 29ddebd997SPeter Brune PetscFunctionBegin; 30ddebd997SPeter Brune 31ddebd997SPeter Brune switch (type) { 32ddebd997SPeter Brune case SNES_LS_BASIC: 33ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 34ddebd997SPeter Brune break; 35ddebd997SPeter Brune case SNES_LS_BASIC_NONORMS: 36ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 37ddebd997SPeter Brune break; 38ddebd997SPeter Brune case SNES_LS_QUADRATIC: 39ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 40ddebd997SPeter Brune break; 41ddebd997SPeter Brune default: 42ddebd997SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 43ddebd997SPeter Brune break; 44ddebd997SPeter Brune } 45ddebd997SPeter Brune snes->ls_type = type; 46ddebd997SPeter Brune PetscFunctionReturn(0); 47ddebd997SPeter Brune } 48ddebd997SPeter Brune EXTERN_C_END 49ddebd997SPeter Brune 50ddebd997SPeter Brune EXTERN_C_BEGIN 51421d9b32SPeter Brune 52421d9b32SPeter Brune #undef __FUNCT__ 53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 54421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 55421d9b32SPeter Brune { 56421d9b32SPeter Brune SNES_FAS * fas; 57421d9b32SPeter Brune PetscErrorCode ierr; 58421d9b32SPeter Brune 59421d9b32SPeter Brune PetscFunctionBegin; 60421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 61421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 62421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 63421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 64421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 65421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 66421d9b32SPeter Brune 67ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 68ed020824SBarry Smith snes->usespc = PETSC_FALSE; 69ed020824SBarry Smith 70421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 71421d9b32SPeter Brune snes->data = (void*) fas; 72421d9b32SPeter Brune fas->level = 0; 73293a7e31SPeter Brune fas->levels = 1; 74ee78dd50SPeter Brune fas->n_cycles = 1; 75ee78dd50SPeter Brune fas->max_up_it = 1; 76ee78dd50SPeter Brune fas->max_down_it = 1; 77ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 78ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 79421d9b32SPeter Brune fas->next = PETSC_NULL; 806273346dSPeter Brune fas->previous = PETSC_NULL; 81421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 82421d9b32SPeter Brune fas->restrct = PETSC_NULL; 83efe1f98aSPeter Brune fas->inject = PETSC_NULL; 84cc05f883SPeter Brune fas->monitor = PETSC_NULL; 85cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 86ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 87ddebd997SPeter Brune 88ddebd997SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 89ddebd997SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 90ddebd997SPeter Brune 91421d9b32SPeter Brune PetscFunctionReturn(0); 92421d9b32SPeter Brune } 93421d9b32SPeter Brune EXTERN_C_END 94421d9b32SPeter Brune 95421d9b32SPeter Brune #undef __FUNCT__ 96421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 97c79ef259SPeter Brune /*@ 98722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 99c79ef259SPeter Brune 100c79ef259SPeter Brune Input Parameter: 101c79ef259SPeter Brune . snes - the preconditioner context 102c79ef259SPeter Brune 103c79ef259SPeter Brune Output parameter: 104c79ef259SPeter Brune . levels - the number of levels 105c79ef259SPeter Brune 106c79ef259SPeter Brune Level: advanced 107c79ef259SPeter Brune 108c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 109c79ef259SPeter Brune 110c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 111c79ef259SPeter Brune @*/ 112421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 113421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 114421d9b32SPeter Brune PetscFunctionBegin; 115ee1fd11aSPeter Brune *levels = fas->levels; 116421d9b32SPeter Brune PetscFunctionReturn(0); 117421d9b32SPeter Brune } 118421d9b32SPeter Brune 119421d9b32SPeter Brune #undef __FUNCT__ 120646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 121c79ef259SPeter Brune /*@ 122c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 123c79ef259SPeter Brune complicated cycling. 124c79ef259SPeter Brune 125c79ef259SPeter Brune Logically Collective on SNES 126c79ef259SPeter Brune 127c79ef259SPeter Brune Input Parameters: 128c79ef259SPeter Brune + snes - the multigrid context 129c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 130c79ef259SPeter Brune 131c79ef259SPeter Brune Options Database Key: 132c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 133c79ef259SPeter Brune 134c79ef259SPeter Brune Level: advanced 135c79ef259SPeter Brune 136c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 137c79ef259SPeter Brune 138c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 139c79ef259SPeter Brune @*/ 140646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 141646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 142646217ecSPeter Brune PetscErrorCode ierr; 143646217ecSPeter Brune PetscFunctionBegin; 144646217ecSPeter Brune fas->n_cycles = cycles; 145646217ecSPeter Brune if (fas->next) { 146646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 147646217ecSPeter Brune } 148646217ecSPeter Brune PetscFunctionReturn(0); 149646217ecSPeter Brune } 150646217ecSPeter Brune 151eff52c0eSPeter Brune #undef __FUNCT__ 152c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 153c79ef259SPeter Brune /*@ 154722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 155c79ef259SPeter Brune 156c79ef259SPeter Brune Logically Collective on SNES 157c79ef259SPeter Brune 158c79ef259SPeter Brune Input Parameters: 159c79ef259SPeter Brune + snes - the multigrid context 160c79ef259SPeter Brune . level - the level to set the number of cycles on 161c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 162c79ef259SPeter Brune 163c79ef259SPeter Brune Level: advanced 164c79ef259SPeter Brune 165c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 166c79ef259SPeter Brune 167c79ef259SPeter Brune .seealso: SNESFASSetCycles() 168c79ef259SPeter Brune @*/ 169c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 170c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 171c79ef259SPeter Brune PetscInt top_level = fas->level,i; 172c79ef259SPeter Brune 173c79ef259SPeter Brune PetscFunctionBegin; 174c79ef259SPeter Brune if (level > top_level) 175c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 176c79ef259SPeter Brune /* get to the correct level */ 177c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 178c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 179c79ef259SPeter Brune } 180c79ef259SPeter Brune if (fas->level != level) 181c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 182c79ef259SPeter Brune fas->n_cycles = cycles; 183c79ef259SPeter Brune PetscFunctionReturn(0); 184c79ef259SPeter Brune } 185c79ef259SPeter Brune 186c79ef259SPeter Brune #undef __FUNCT__ 187eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 188aeed3662SMatthew G Knepley /*@C 189c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 190c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 191c79ef259SPeter Brune and nonlinear preconditioners. 192c79ef259SPeter Brune 193c79ef259SPeter Brune Logically Collective on SNES 194c79ef259SPeter Brune 195c79ef259SPeter Brune Input Parameters: 196c79ef259SPeter Brune + snes - the multigrid context 197c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 198c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 199c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 200c79ef259SPeter Brune 201c79ef259SPeter Brune Level: advanced 202c79ef259SPeter Brune 203c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 204c79ef259SPeter Brune 205c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 206c79ef259SPeter Brune @*/ 207eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 208eff52c0eSPeter Brune PetscErrorCode ierr = 0; 209eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 210eff52c0eSPeter Brune PetscFunctionBegin; 211eff52c0eSPeter Brune 212eff52c0eSPeter Brune /* use or don't use it according to user wishes*/ 213d28543b3SPeter Brune snes->usegs = use_gs; 214eff52c0eSPeter Brune if (gsfunc) { 215eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 216eff52c0eSPeter Brune /* push the provided GS up the tree */ 217eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 218eff52c0eSPeter Brune } else if (snes->ops->computegs) { 219eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 220eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 221eff52c0eSPeter Brune } else if (use_gs) { 222eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 223eff52c0eSPeter Brune } 224eff52c0eSPeter Brune PetscFunctionReturn(0); 225eff52c0eSPeter Brune } 226eff52c0eSPeter Brune 227eff52c0eSPeter Brune #undef __FUNCT__ 228eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 229aeed3662SMatthew G Knepley /*@C 230c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 231c79ef259SPeter Brune 232c79ef259SPeter Brune Logically Collective on SNES 233c79ef259SPeter Brune 234c79ef259SPeter Brune Input Parameters: 235c79ef259SPeter Brune + snes - the multigrid context 236c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 237c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 238c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 239c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 240c79ef259SPeter Brune 241c79ef259SPeter Brune Level: advanced 242c79ef259SPeter Brune 243c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 244c79ef259SPeter Brune 245c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 246c79ef259SPeter Brune @*/ 247eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 248eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 249eff52c0eSPeter Brune PetscErrorCode ierr; 250eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 251eff52c0eSPeter Brune SNES cur_snes = snes; 252eff52c0eSPeter Brune PetscFunctionBegin; 253eff52c0eSPeter Brune if (level > top_level) 254eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 255eff52c0eSPeter Brune /* get to the correct level */ 256eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 257eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 258eff52c0eSPeter Brune cur_snes = fas->next; 259eff52c0eSPeter Brune } 260eff52c0eSPeter Brune if (fas->level != level) 261eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 262d28543b3SPeter Brune snes->usegs = use_gs; 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 394c79ef259SPeter Brune Logically Collective on PC 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 426c79ef259SPeter Brune Logically Collective on PC 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 474c79ef259SPeter Brune .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; 479421d9b32SPeter Brune 480421d9b32SPeter Brune PetscFunctionBegin; 481421d9b32SPeter Brune if (level > top_level) 482421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 483421d9b32SPeter Brune /* get to the correct level */ 484d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 485421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 486421d9b32SPeter Brune } 487421d9b32SPeter Brune if (fas->level != level) 488421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 489421d9b32SPeter Brune fas->interpolate = mat; 490421d9b32SPeter Brune PetscFunctionReturn(0); 491421d9b32SPeter Brune } 492421d9b32SPeter Brune 493421d9b32SPeter Brune #undef __FUNCT__ 494421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 495c79ef259SPeter Brune /*@ 496c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 497c79ef259SPeter Brune from level l to l-1. 498c79ef259SPeter Brune 499c79ef259SPeter Brune Input Parameters: 500c79ef259SPeter Brune + snes - the multigrid context 501c79ef259SPeter Brune . mat - the restriction matrix 502c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 503c79ef259SPeter Brune 504c79ef259SPeter Brune Level: advanced 505c79ef259SPeter Brune 506c79ef259SPeter Brune Notes: 507c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 508c79ef259SPeter Brune for the same level. 509c79ef259SPeter Brune 510c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 511c79ef259SPeter Brune out from the matrix size which one it is. 512c79ef259SPeter Brune 513c79ef259SPeter Brune If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 514c79ef259SPeter Brune is used. 515c79ef259SPeter Brune 516c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 517c79ef259SPeter Brune 518c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 519c79ef259SPeter Brune @*/ 520421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 521421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 522d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 523421d9b32SPeter Brune 524421d9b32SPeter Brune PetscFunctionBegin; 525421d9b32SPeter Brune if (level > top_level) 526421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 527421d9b32SPeter Brune /* get to the correct level */ 528d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 529421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 530421d9b32SPeter Brune } 531421d9b32SPeter Brune if (fas->level != level) 532421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 533421d9b32SPeter Brune fas->restrct = mat; 534421d9b32SPeter Brune PetscFunctionReturn(0); 535421d9b32SPeter Brune } 536421d9b32SPeter Brune 537efe1f98aSPeter Brune 538421d9b32SPeter Brune #undef __FUNCT__ 539421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 540c79ef259SPeter Brune /*@ 541c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 542c79ef259SPeter Brune operator from level l to l-1. 543c79ef259SPeter Brune 544c79ef259SPeter Brune Input Parameters: 545c79ef259SPeter Brune + snes - the multigrid context 546c79ef259SPeter Brune . rscale - the restriction scaling 547c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 548c79ef259SPeter Brune 549c79ef259SPeter Brune Level: advanced 550c79ef259SPeter Brune 551c79ef259SPeter Brune Notes: 552c79ef259SPeter Brune This is only used in the case that the injection is not set. 553c79ef259SPeter Brune 554c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 555c79ef259SPeter Brune 556c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 557c79ef259SPeter Brune @*/ 558421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 559421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 560d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 561421d9b32SPeter Brune 562421d9b32SPeter Brune PetscFunctionBegin; 563421d9b32SPeter Brune if (level > top_level) 564421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 565421d9b32SPeter Brune /* get to the correct level */ 566d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 567421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 568421d9b32SPeter Brune } 569421d9b32SPeter Brune if (fas->level != level) 570421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 571421d9b32SPeter Brune fas->rscale = rscale; 572421d9b32SPeter Brune PetscFunctionReturn(0); 573421d9b32SPeter Brune } 574421d9b32SPeter Brune 575efe1f98aSPeter Brune 576efe1f98aSPeter Brune #undef __FUNCT__ 577efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 578c79ef259SPeter Brune /*@ 579c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 580c79ef259SPeter Brune from level l to l-1. 581c79ef259SPeter Brune 582c79ef259SPeter Brune Input Parameters: 583c79ef259SPeter Brune + snes - the multigrid context 584c79ef259SPeter Brune . mat - the restriction matrix 585c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 586c79ef259SPeter Brune 587c79ef259SPeter Brune Level: advanced 588c79ef259SPeter Brune 589c79ef259SPeter Brune Notes: 590c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 591c79ef259SPeter Brune project the solution instead. 592c79ef259SPeter Brune 593c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 594c79ef259SPeter Brune 595c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 596c79ef259SPeter Brune @*/ 597efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 598efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 599efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 600efe1f98aSPeter Brune 601efe1f98aSPeter Brune PetscFunctionBegin; 602efe1f98aSPeter Brune if (level > top_level) 603efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 604efe1f98aSPeter Brune /* get to the correct level */ 605efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 606efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 607efe1f98aSPeter Brune } 608efe1f98aSPeter Brune if (fas->level != level) 609efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 610efe1f98aSPeter Brune fas->inject = mat; 611efe1f98aSPeter Brune PetscFunctionReturn(0); 612efe1f98aSPeter Brune } 613efe1f98aSPeter Brune 614421d9b32SPeter Brune #undef __FUNCT__ 615421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 616421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 617421d9b32SPeter Brune { 61877df8cc4SPeter Brune PetscErrorCode ierr = 0; 619421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 620421d9b32SPeter Brune 621421d9b32SPeter Brune PetscFunctionBegin; 622742fe5e2SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 623742fe5e2SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 624742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 625421d9b32SPeter Brune PetscFunctionReturn(0); 626421d9b32SPeter Brune } 627421d9b32SPeter Brune 628421d9b32SPeter Brune #undef __FUNCT__ 629421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 630421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 631421d9b32SPeter Brune { 632421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 633742fe5e2SPeter Brune PetscErrorCode ierr = 0; 634421d9b32SPeter Brune 635421d9b32SPeter Brune PetscFunctionBegin; 636421d9b32SPeter Brune /* recursively resets and then destroys */ 63751e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 63851e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 639efe1f98aSPeter Brune if (fas->inject) { 640efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 641efe1f98aSPeter Brune } 64251e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 64351e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 64451e86f29SPeter Brune fas->restrct = PETSC_NULL; 64551e86f29SPeter Brune } else { 64651e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 64751e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 64851e86f29SPeter Brune } 64951e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 650421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 651421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 652ee1fd11aSPeter Brune 653421d9b32SPeter Brune PetscFunctionReturn(0); 654421d9b32SPeter Brune } 655421d9b32SPeter Brune 656421d9b32SPeter Brune #undef __FUNCT__ 657421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 658421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 659421d9b32SPeter Brune { 66048bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 661421d9b32SPeter Brune PetscErrorCode ierr; 662efe1f98aSPeter Brune VecScatter injscatter; 663d1adcc6fSPeter Brune PetscInt dm_levels; 664d1adcc6fSPeter Brune 665421d9b32SPeter Brune PetscFunctionBegin; 666eff52c0eSPeter Brune 667cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 668d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 669d1adcc6fSPeter Brune dm_levels++; 670cc05f883SPeter Brune if (dm_levels > fas->levels) { 671d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 672cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 673d1adcc6fSPeter Brune } 674d1adcc6fSPeter Brune } 675d1adcc6fSPeter Brune 67607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 6772f7ea302SPeter Brune ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 67807144faaSPeter Brune } else { 679ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 68007144faaSPeter Brune } 681cc05f883SPeter Brune 68248bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 68348bfdf8aSPeter Brune if (fas->upsmooth) { 68448bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 68548bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 68648bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 68748bfdf8aSPeter Brune } 68848bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 68948bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 69048bfdf8aSPeter Brune } 69148bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 69248bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 69348bfdf8aSPeter Brune } 69448bfdf8aSPeter Brune } 69548bfdf8aSPeter Brune if (fas->downsmooth) { 69648bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 69748bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 69848bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 69948bfdf8aSPeter Brune } 70048bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 70148bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 70248bfdf8aSPeter Brune } 70348bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 70448bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 70548bfdf8aSPeter Brune } 7061a266240SBarry Smith } 70748bfdf8aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 708646217ecSPeter Brune if (fas->next) { 7096273346dSPeter Brune if (fas->galerkin) { 7106273346dSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 7116273346dSPeter Brune } else { 71248bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 71348bfdf8aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 71448bfdf8aSPeter Brune } 71548bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 71648bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 71748bfdf8aSPeter Brune } 71848bfdf8aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 719646217ecSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 720646217ecSPeter Brune } 721646217ecSPeter Brune } 7226273346dSPeter Brune } 723421d9b32SPeter Brune if (snes->dm) { 724421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 725421d9b32SPeter Brune if (fas->next) { 726421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 727ee78dd50SPeter Brune if (!fas->next->dm) { 728ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 729ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 730ee78dd50SPeter Brune } 731421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 732ee78dd50SPeter Brune if (!fas->interpolate) { 733e727c939SJed Brown ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 734421d9b32SPeter Brune fas->restrct = fas->interpolate; 735421d9b32SPeter Brune } 736efe1f98aSPeter Brune /* set the injection from the DM */ 737efe1f98aSPeter Brune if (!fas->inject) { 738e727c939SJed Brown ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 739efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 740efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 741efe1f98aSPeter Brune } 742421d9b32SPeter Brune } 743ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 744*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 745*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 746421d9b32SPeter Brune } 747cc05f883SPeter Brune 7486273346dSPeter Brune /* setup FAS work vectors */ 7496273346dSPeter Brune if (fas->galerkin) { 7506273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7516273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7526273346dSPeter Brune } 7536273346dSPeter Brune 754fa9694d7SPeter Brune if (fas->next) { 755fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 756*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 757742fe5e2SPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 758fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 759fa9694d7SPeter Brune } 760fa9694d7SPeter Brune /* got to set them all up at once */ 761421d9b32SPeter Brune PetscFunctionReturn(0); 762421d9b32SPeter Brune } 763421d9b32SPeter Brune 764421d9b32SPeter Brune #undef __FUNCT__ 765421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 766421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 767421d9b32SPeter Brune { 768ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 769ee78dd50SPeter Brune PetscInt levels = 1; 770d1adcc6fSPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 771421d9b32SPeter Brune PetscErrorCode ierr; 772ee78dd50SPeter Brune const char *def_smooth = SNESNRICHARDSON; 773ee78dd50SPeter Brune char pre_type[256]; 774ee78dd50SPeter Brune char post_type[256]; 775ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 77607144faaSPeter Brune SNESFASType fastype; 777421d9b32SPeter Brune 778421d9b32SPeter Brune PetscFunctionBegin; 779c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 780ee78dd50SPeter Brune 781ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 782ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 783ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 784c732cbdbSBarry Smith if (!flg && snes->dm) { 785c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 786c732cbdbSBarry Smith levels++; 787d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 788c732cbdbSBarry Smith } 789ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 790ee78dd50SPeter Brune } 791ee78dd50SPeter Brune 792ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 793ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 794ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 79507144faaSPeter Brune fastype = fas->fastype; 79607144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 79707144faaSPeter Brune if (flg) { 79807144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 79907144faaSPeter Brune } 800d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 801d1adcc6fSPeter Brune if (smoothflg) { 802ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 803ee78dd50SPeter Brune } else { 804d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 805d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 806ee78dd50SPeter Brune } 807ee78dd50SPeter Brune 808ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 809d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 810d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 811d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 812ee78dd50SPeter Brune 8136273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8146273346dSPeter Brune 815646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 816ee78dd50SPeter Brune 817ee78dd50SPeter Brune /* other options for the coarsest level */ 818ee78dd50SPeter Brune if (fas->level == 0) { 819d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 820ee78dd50SPeter Brune } 821ee78dd50SPeter Brune 822421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8238cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 824eff52c0eSPeter Brune 825d28543b3SPeter Brune if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) { 8268cc86e31SPeter Brune const char *prefix; 8278cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8288cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8298cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 8308cc86e31SPeter Brune if (fas->level || (fas->levels == 1)) { 831eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 8328cc86e31SPeter Brune } else { 8338cc86e31SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 8348cc86e31SPeter Brune } 8358cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8368cc86e31SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 8378cc86e31SPeter Brune } 8388cc86e31SPeter Brune 839d28543b3SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) { 84067339d5cSBarry Smith const char *prefix; 84167339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 842ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 84367339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 844eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 845293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 846ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 847ee78dd50SPeter Brune } 8481a266240SBarry Smith if (fas->upsmooth) { 8491a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 8501a266240SBarry Smith } 8511a266240SBarry Smith 8521a266240SBarry Smith if (fas->downsmooth) { 8531a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 8541a266240SBarry Smith } 855ee78dd50SPeter Brune 856742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 857742fe5e2SPeter Brune ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr); 858742fe5e2SPeter Brune } 859742fe5e2SPeter Brune 860ee78dd50SPeter Brune if (monflg) { 861646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 862794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 8632f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 864742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 865293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 866293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 867d28543b3SPeter Brune } else { 868d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 869d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 870d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 871d28543b3SPeter Brune } 872ee78dd50SPeter Brune } 873ee78dd50SPeter Brune 874ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 875d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 876421d9b32SPeter Brune PetscFunctionReturn(0); 877421d9b32SPeter Brune } 878421d9b32SPeter Brune 879421d9b32SPeter Brune #undef __FUNCT__ 880421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 881421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 882421d9b32SPeter Brune { 883421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 884421d9b32SPeter Brune PetscBool iascii; 885421d9b32SPeter Brune PetscErrorCode ierr; 886421d9b32SPeter Brune 887421d9b32SPeter Brune PetscFunctionBegin; 88839a4b097SPeter Brune if (fas->level != fas->levels - 1) PetscFunctionReturn(0); 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 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 894ee78dd50SPeter Brune if (fas->upsmooth) { 89539a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 896421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 897ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 898421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 899421d9b32SPeter Brune } else { 900ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 901421d9b32SPeter Brune } 902ee78dd50SPeter Brune if (fas->downsmooth) { 90339a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 904421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 905ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 906421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 907421d9b32SPeter Brune } else { 908ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 909421d9b32SPeter Brune } 910d28543b3SPeter Brune if (snes->usegs) { 91139a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D -- smoothdown=%D, smoothup=%D\n", 91239a4b097SPeter Brune fas->level, fas->max_down_it, fas->max_up_it);CHKERRQ(ierr); 913421d9b32SPeter Brune } 914421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 915421d9b32SPeter Brune } else { 916421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 917421d9b32SPeter Brune } 918421d9b32SPeter Brune PetscFunctionReturn(0); 919421d9b32SPeter Brune } 920421d9b32SPeter Brune 921421d9b32SPeter Brune #undef __FUNCT__ 92239bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 92339bd7f45SPeter Brune /* 92439bd7f45SPeter Brune Defines the action of the downsmoother 92539bd7f45SPeter Brune */ 92639bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 92739bd7f45SPeter Brune PetscErrorCode ierr = 0; 9282d15c84fSPeter Brune PetscReal fnorm; 929742fe5e2SPeter Brune SNESConvergedReason reason; 93039bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 931421d9b32SPeter Brune PetscFunctionBegin; 932d1adcc6fSPeter Brune if (fas->downsmooth) { 933d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 934742fe5e2SPeter Brune /* check convergence reason for the smoother */ 935742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 936742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 937742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 938742fe5e2SPeter Brune PetscFunctionReturn(0); 939742fe5e2SPeter Brune } 940d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 941794bee33SPeter Brune if (fas->monitor) { 942794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 943794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 944d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 945eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 946d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 947794bee33SPeter Brune } 9489c44eddcSPeter Brune ierr = SNESSetGSSweeps(snes, fas->max_down_it);CHKERRQ(ierr); 949646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 950cc05f883SPeter Brune if (fas->monitor) { 951794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 952794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 953d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 9549c44eddcSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "1 SNES GS Function norm %14.12e\n", fnorm);CHKERRQ(ierr); 955d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 956794bee33SPeter Brune } 957c90fad12SPeter Brune } else if (snes->pc) { 958c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 9592f7ea302SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 9602f7ea302SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 9612f7ea302SPeter Brune snes->reason = SNES_DIVERGED_INNER; 9622f7ea302SPeter Brune PetscFunctionReturn(0); 9632f7ea302SPeter Brune } 964fe6f9142SPeter Brune } 96539bd7f45SPeter Brune PetscFunctionReturn(0); 96639bd7f45SPeter Brune } 96739bd7f45SPeter Brune 96839bd7f45SPeter Brune 96939bd7f45SPeter Brune #undef __FUNCT__ 97039bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 97139bd7f45SPeter Brune /* 97207144faaSPeter Brune Defines the action of the upsmoother 97339bd7f45SPeter Brune */ 97439bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 97539bd7f45SPeter Brune PetscErrorCode ierr = 0; 97639bd7f45SPeter Brune PetscReal fnorm; 97739bd7f45SPeter Brune SNESConvergedReason reason; 97839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 97939bd7f45SPeter Brune PetscFunctionBegin; 98039bd7f45SPeter Brune if (fas->upsmooth) { 98139bd7f45SPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 98239bd7f45SPeter Brune /* check convergence reason for the smoother */ 98339bd7f45SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 98439bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 98539bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 98639bd7f45SPeter Brune PetscFunctionReturn(0); 98739bd7f45SPeter Brune } 98839bd7f45SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 98939bd7f45SPeter Brune if (fas->monitor) { 99039bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 99139bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 99239bd7f45SPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 99339bd7f45SPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 99439bd7f45SPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 99539bd7f45SPeter Brune } 9969c44eddcSPeter Brune ierr = SNESSetGSSweeps(snes, fas->max_up_it);CHKERRQ(ierr); 99739bd7f45SPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 99839bd7f45SPeter Brune if (fas->monitor) { 99939bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 100039bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 100139bd7f45SPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 10029c44eddcSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "1 SNES GS Function norm %14.12e\n", fnorm);CHKERRQ(ierr); 100339bd7f45SPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 100439bd7f45SPeter Brune } 100539bd7f45SPeter Brune } else if (snes->pc) { 100639bd7f45SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 100739bd7f45SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 100839bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 100939bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 101039bd7f45SPeter Brune PetscFunctionReturn(0); 101139bd7f45SPeter Brune } 101239bd7f45SPeter Brune } 101339bd7f45SPeter Brune PetscFunctionReturn(0); 101439bd7f45SPeter Brune } 101539bd7f45SPeter Brune 101639bd7f45SPeter Brune #undef __FUNCT__ 101739bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 101839bd7f45SPeter Brune /* 101939bd7f45SPeter Brune 102039bd7f45SPeter Brune Performs the FAS coarse correction as: 102139bd7f45SPeter Brune 102239bd7f45SPeter Brune fine problem: F(x) = 0 102339bd7f45SPeter Brune coarse problem: F^c(x) = b^c 102439bd7f45SPeter Brune 102539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 102639bd7f45SPeter Brune 102739bd7f45SPeter Brune with correction: 102839bd7f45SPeter Brune 102939bd7f45SPeter Brune 103039bd7f45SPeter Brune 103139bd7f45SPeter Brune */ 103239a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 103339bd7f45SPeter Brune PetscErrorCode ierr; 103439bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 103539bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 103639bd7f45SPeter Brune SNESConvergedReason reason; 103739bd7f45SPeter Brune PetscFunctionBegin; 1038fa9694d7SPeter Brune if (fas->next) { 1039ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1040794bee33SPeter Brune 1041c90fad12SPeter Brune X_c = fas->next->vec_sol; 1042293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1043c90fad12SPeter Brune F_c = fas->next->vec_func; 1044742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1045efe1f98aSPeter Brune 1046efe1f98aSPeter Brune /* inject the solution */ 1047efe1f98aSPeter Brune if (fas->inject) { 1048a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 1049efe1f98aSPeter Brune } else { 1050a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 1051a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1052efe1f98aSPeter Brune } 1053293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1054293a7e31SPeter Brune 1055293a7e31SPeter Brune /* restrict the defect */ 1056293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1057293a7e31SPeter Brune 1058c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1059ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1060c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1061742fe5e2SPeter Brune 1062293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1063c90fad12SPeter Brune 1064ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1065ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1066c90fad12SPeter Brune 1067c90fad12SPeter Brune /* recurse to the next level */ 1068f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1069742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1070742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1071742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1072742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1073742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1074742fe5e2SPeter Brune PetscFunctionReturn(0); 1075742fe5e2SPeter Brune } 1076742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1077ee78dd50SPeter Brune 1078fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1079fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 108039bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1081293a7e31SPeter Brune } 108239bd7f45SPeter Brune PetscFunctionReturn(0); 108339bd7f45SPeter Brune } 108439bd7f45SPeter Brune 108539bd7f45SPeter Brune #undef __FUNCT__ 108639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 108739bd7f45SPeter Brune /* 108839bd7f45SPeter Brune 108939bd7f45SPeter Brune The additive cycle looks like: 109039bd7f45SPeter Brune 109107144faaSPeter Brune xhat = x 109207144faaSPeter Brune xhat = dS(x, b) 109307144faaSPeter Brune x = coarsecorrection(xhat, b_d) 109407144faaSPeter Brune x = x + nu*(xhat - x); 109539bd7f45SPeter Brune (optional) x = uS(x, b) 109639bd7f45SPeter Brune 109739bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 109839bd7f45SPeter Brune 109939bd7f45SPeter Brune */ 110039bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 110107144faaSPeter Brune Vec F, B, Xhat; 1102ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 110339bd7f45SPeter Brune PetscErrorCode ierr; 110407144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 110507144faaSPeter Brune SNESConvergedReason reason; 1106ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1107ddebd997SPeter Brune PetscBool lssucceed; 110839bd7f45SPeter Brune PetscFunctionBegin; 110939bd7f45SPeter Brune 111039bd7f45SPeter Brune F = snes->vec_func; 111139bd7f45SPeter Brune B = snes->vec_rhs; 111207144faaSPeter Brune Xhat = snes->work[1]; 1113ddebd997SPeter Brune G = snes->work[2]; 1114ddebd997SPeter Brune W = snes->work[3]; 111507144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 111607144faaSPeter Brune /* recurse first */ 111707144faaSPeter Brune if (fas->next) { 111807144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 111939bd7f45SPeter Brune 112007144faaSPeter Brune X_c = fas->next->vec_sol; 112107144faaSPeter Brune Xo_c = fas->next->work[0]; 112207144faaSPeter Brune F_c = fas->next->vec_func; 112307144faaSPeter Brune B_c = fas->next->vec_rhs; 112439bd7f45SPeter Brune 112507144faaSPeter Brune /* inject the solution */ 112607144faaSPeter Brune if (fas->inject) { 112707144faaSPeter Brune ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr); 112807144faaSPeter Brune } else { 112907144faaSPeter Brune ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr); 113007144faaSPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 113107144faaSPeter Brune } 113207144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 113307144faaSPeter Brune 113407144faaSPeter Brune /* restrict the defect */ 113507144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 113607144faaSPeter Brune 113707144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 113807144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 113907144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 114007144faaSPeter Brune 114107144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 114207144faaSPeter Brune 114307144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 114407144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 114507144faaSPeter Brune 114607144faaSPeter Brune /* recurse */ 114707144faaSPeter Brune fas->next->vec_rhs = B_c; 114807144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 114907144faaSPeter Brune 115007144faaSPeter Brune /* smooth on this level */ 115107144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 115207144faaSPeter Brune 115307144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 115407144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 115507144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 115607144faaSPeter Brune PetscFunctionReturn(0); 115707144faaSPeter Brune } 115807144faaSPeter Brune 115907144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 116007144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1161ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 116207144faaSPeter Brune 1163ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1164ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1165ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1166eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 1167ddebd997SPeter Brune ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1168ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1169ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1170ddebd997SPeter Brune fnorm = gnorm; 117107144faaSPeter Brune } else { 117207144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 117307144faaSPeter Brune } 117439bd7f45SPeter Brune PetscFunctionReturn(0); 117539bd7f45SPeter Brune } 117639bd7f45SPeter Brune 117739bd7f45SPeter Brune #undef __FUNCT__ 117839bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 117939bd7f45SPeter Brune /* 118039bd7f45SPeter Brune 118139bd7f45SPeter Brune Defines the FAS cycle as: 118239bd7f45SPeter Brune 118339bd7f45SPeter Brune fine problem: F(x) = 0 118439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 118539bd7f45SPeter Brune 118639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 118739bd7f45SPeter Brune 118839bd7f45SPeter Brune correction: 118939bd7f45SPeter Brune 119039bd7f45SPeter Brune x = x + I(x^c - Rx) 119139bd7f45SPeter Brune 119239bd7f45SPeter Brune */ 119339bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 119439bd7f45SPeter Brune 119539bd7f45SPeter Brune PetscErrorCode ierr; 119639bd7f45SPeter Brune Vec F,B; 119739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 119839bd7f45SPeter Brune 119939bd7f45SPeter Brune PetscFunctionBegin; 120039bd7f45SPeter Brune F = snes->vec_func; 120139bd7f45SPeter Brune B = snes->vec_rhs; 120239bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 120339bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 120439bd7f45SPeter Brune 120539a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 120639bd7f45SPeter Brune 1207c90fad12SPeter Brune if (fas->level != 0) { 120839bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1209fe6f9142SPeter Brune } 1210fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1211fa9694d7SPeter Brune 1212fa9694d7SPeter Brune PetscFunctionReturn(0); 1213421d9b32SPeter Brune } 1214421d9b32SPeter Brune 1215421d9b32SPeter Brune #undef __FUNCT__ 1216421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1217421d9b32SPeter Brune 1218421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1219421d9b32SPeter Brune { 1220fa9694d7SPeter Brune PetscErrorCode ierr; 1221fe6f9142SPeter Brune PetscInt i, maxits; 1222ddb5aff1SPeter Brune Vec X, F; 1223fe6f9142SPeter Brune PetscReal fnorm; 122407144faaSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 1225421d9b32SPeter Brune PetscFunctionBegin; 1226fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1227fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1228fa9694d7SPeter Brune X = snes->vec_sol; 1229f5a6d4f9SBarry Smith F = snes->vec_func; 1230293a7e31SPeter Brune 1231293a7e31SPeter Brune /*norm setup */ 1232fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1233fe6f9142SPeter Brune snes->iter = 0; 1234fe6f9142SPeter Brune snes->norm = 0.; 1235fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1236fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1237fe6f9142SPeter Brune if (snes->domainerror) { 1238fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1239fe6f9142SPeter Brune PetscFunctionReturn(0); 1240fe6f9142SPeter Brune } 1241fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1242fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1243fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1244fe6f9142SPeter Brune snes->norm = fnorm; 1245fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1246fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1247fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1248fe6f9142SPeter Brune 1249fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1250fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1251fe6f9142SPeter Brune /* test convergence */ 1252fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1253fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1254fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1255fe6f9142SPeter Brune /* Call general purpose update function */ 1256646217ecSPeter Brune 1257fe6f9142SPeter Brune if (snes->ops->update) { 1258fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1259fe6f9142SPeter Brune } 126007144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 126139bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 126207144faaSPeter Brune } else { 126307144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 126407144faaSPeter Brune } 1265742fe5e2SPeter Brune 1266742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1267742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1268742fe5e2SPeter Brune PetscFunctionReturn(0); 1269742fe5e2SPeter Brune } 1270c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1271c90fad12SPeter Brune /* Monitor convergence */ 1272c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1273c90fad12SPeter Brune snes->iter = i+1; 1274c90fad12SPeter Brune snes->norm = fnorm; 1275c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1276c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1277c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1278c90fad12SPeter Brune /* Test for convergence */ 1279c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1280c90fad12SPeter Brune if (snes->reason) break; 1281fe6f9142SPeter Brune } 1282fe6f9142SPeter Brune if (i == maxits) { 1283fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1284fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1285fe6f9142SPeter Brune } 1286421d9b32SPeter Brune PetscFunctionReturn(0); 1287421d9b32SPeter Brune } 1288