1421d9b32SPeter Brune /* Defines the basic SNES object */ 2421d9b32SPeter Brune #include <../src/snes/impls/fas/fasimpls.h> 3421d9b32SPeter Brune 4421d9b32SPeter Brune /*MC 5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver. 6421d9b32SPeter Brune 7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8421d9b32SPeter Brune 9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10421d9b32SPeter Brune M*/ 11421d9b32SPeter Brune 12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 18421d9b32SPeter Brune 19421d9b32SPeter Brune EXTERN_C_BEGIN 20421d9b32SPeter Brune 21421d9b32SPeter Brune #undef __FUNCT__ 22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 24421d9b32SPeter Brune { 25421d9b32SPeter Brune SNES_FAS * fas; 26421d9b32SPeter Brune PetscErrorCode ierr; 27421d9b32SPeter Brune 28421d9b32SPeter Brune PetscFunctionBegin; 29421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 30421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 31421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 33421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 34421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 35421d9b32SPeter Brune 36ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 37ed020824SBarry Smith snes->usespc = PETSC_FALSE; 38ed020824SBarry Smith 39421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 40421d9b32SPeter Brune snes->data = (void*) fas; 41421d9b32SPeter Brune fas->level = 0; 42293a7e31SPeter Brune fas->levels = 1; 43ee78dd50SPeter Brune fas->n_cycles = 1; 44ee78dd50SPeter Brune fas->max_up_it = 1; 45ee78dd50SPeter Brune fas->max_down_it = 1; 46ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 47ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 48421d9b32SPeter Brune fas->next = PETSC_NULL; 49421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 50421d9b32SPeter Brune fas->restrct = PETSC_NULL; 51efe1f98aSPeter Brune fas->inject = PETSC_NULL; 52cc05f883SPeter Brune fas->monitor = PETSC_NULL; 53cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 54eff52c0eSPeter Brune fas->useGS = PETSC_FALSE; 55421d9b32SPeter Brune PetscFunctionReturn(0); 56421d9b32SPeter Brune } 57421d9b32SPeter Brune EXTERN_C_END 58421d9b32SPeter Brune 59421d9b32SPeter Brune #undef __FUNCT__ 60421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 61c79ef259SPeter Brune /*@ 62722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 63c79ef259SPeter Brune 64c79ef259SPeter Brune Input Parameter: 65c79ef259SPeter Brune . snes - the preconditioner context 66c79ef259SPeter Brune 67c79ef259SPeter Brune Output parameter: 68c79ef259SPeter Brune . levels - the number of levels 69c79ef259SPeter Brune 70c79ef259SPeter Brune Level: advanced 71c79ef259SPeter Brune 72c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 73c79ef259SPeter Brune 74c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 75c79ef259SPeter Brune @*/ 76421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 77421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 78421d9b32SPeter Brune PetscFunctionBegin; 79ee1fd11aSPeter Brune *levels = fas->levels; 80421d9b32SPeter Brune PetscFunctionReturn(0); 81421d9b32SPeter Brune } 82421d9b32SPeter Brune 83421d9b32SPeter Brune #undef __FUNCT__ 84646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 85c79ef259SPeter Brune /*@ 86c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 87c79ef259SPeter Brune complicated cycling. 88c79ef259SPeter Brune 89c79ef259SPeter Brune Logically Collective on SNES 90c79ef259SPeter Brune 91c79ef259SPeter Brune Input Parameters: 92c79ef259SPeter Brune + snes - the multigrid context 93c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 94c79ef259SPeter Brune 95c79ef259SPeter Brune Options Database Key: 96c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 97c79ef259SPeter Brune 98c79ef259SPeter Brune Level: advanced 99c79ef259SPeter Brune 100c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 101c79ef259SPeter Brune 102c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 103c79ef259SPeter Brune @*/ 104646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 105646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 106646217ecSPeter Brune PetscErrorCode ierr; 107646217ecSPeter Brune PetscFunctionBegin; 108646217ecSPeter Brune fas->n_cycles = cycles; 109646217ecSPeter Brune if (fas->next) { 110646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 111646217ecSPeter Brune } 112646217ecSPeter Brune PetscFunctionReturn(0); 113646217ecSPeter Brune } 114646217ecSPeter Brune 115eff52c0eSPeter Brune #undef __FUNCT__ 116c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 117c79ef259SPeter Brune /*@ 118722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 119c79ef259SPeter Brune 120c79ef259SPeter Brune Logically Collective on SNES 121c79ef259SPeter Brune 122c79ef259SPeter Brune Input Parameters: 123c79ef259SPeter Brune + snes - the multigrid context 124c79ef259SPeter Brune . level - the level to set the number of cycles on 125c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 126c79ef259SPeter Brune 127c79ef259SPeter Brune Level: advanced 128c79ef259SPeter Brune 129c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 130c79ef259SPeter Brune 131c79ef259SPeter Brune .seealso: SNESFASSetCycles() 132c79ef259SPeter Brune @*/ 133c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 134c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 135c79ef259SPeter Brune PetscInt top_level = fas->level,i; 136c79ef259SPeter Brune 137c79ef259SPeter Brune PetscFunctionBegin; 138c79ef259SPeter Brune if (level > top_level) 139c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 140c79ef259SPeter Brune /* get to the correct level */ 141c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 142c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 143c79ef259SPeter Brune } 144c79ef259SPeter Brune if (fas->level != level) 145c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 146c79ef259SPeter Brune fas->n_cycles = cycles; 147c79ef259SPeter Brune PetscFunctionReturn(0); 148c79ef259SPeter Brune } 149c79ef259SPeter Brune 150c79ef259SPeter Brune #undef __FUNCT__ 151eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 152c79ef259SPeter Brune /*@ 153c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 154c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 155c79ef259SPeter Brune and nonlinear preconditioners. 156c79ef259SPeter Brune 157c79ef259SPeter Brune Logically Collective on SNES 158c79ef259SPeter Brune 159c79ef259SPeter Brune Input Parameters: 160c79ef259SPeter Brune + snes - the multigrid context 161c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 162c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 163c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 164c79ef259SPeter Brune 165c79ef259SPeter Brune Level: advanced 166c79ef259SPeter Brune 167c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 168c79ef259SPeter Brune 169c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 170c79ef259SPeter Brune @*/ 171eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 172eff52c0eSPeter Brune PetscErrorCode ierr = 0; 173eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 174eff52c0eSPeter Brune PetscFunctionBegin; 175eff52c0eSPeter Brune 176eff52c0eSPeter Brune /* use or don't use it according to user wishes*/ 177eff52c0eSPeter Brune fas->useGS = use_gs; 178eff52c0eSPeter Brune if (gsfunc) { 179eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 180eff52c0eSPeter Brune /* push the provided GS up the tree */ 181eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 182eff52c0eSPeter Brune } else if (snes->ops->computegs) { 183eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 184eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 185eff52c0eSPeter Brune } else if (use_gs) { 186eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 187eff52c0eSPeter Brune } 188eff52c0eSPeter Brune PetscFunctionReturn(0); 189eff52c0eSPeter Brune } 190eff52c0eSPeter Brune 191eff52c0eSPeter Brune #undef __FUNCT__ 192eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 193c79ef259SPeter Brune /*@ 194c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 195c79ef259SPeter Brune 196c79ef259SPeter Brune Logically Collective on SNES 197c79ef259SPeter Brune 198c79ef259SPeter Brune Input Parameters: 199c79ef259SPeter Brune + snes - the multigrid context 200c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 201c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 202c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 203c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 204c79ef259SPeter Brune 205c79ef259SPeter Brune Level: advanced 206c79ef259SPeter Brune 207c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 208c79ef259SPeter Brune 209c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 210c79ef259SPeter Brune @*/ 211eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 212eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 213eff52c0eSPeter Brune PetscErrorCode ierr; 214eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 215eff52c0eSPeter Brune SNES cur_snes = snes; 216eff52c0eSPeter Brune PetscFunctionBegin; 217eff52c0eSPeter Brune if (level > top_level) 218eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 219eff52c0eSPeter Brune /* get to the correct level */ 220eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 221eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 222eff52c0eSPeter Brune cur_snes = fas->next; 223eff52c0eSPeter Brune } 224eff52c0eSPeter Brune if (fas->level != level) 225eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 226eff52c0eSPeter Brune fas->useGS = use_gs; 227eff52c0eSPeter Brune if (gsfunc) { 228eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 229eff52c0eSPeter Brune } 230eff52c0eSPeter Brune PetscFunctionReturn(0); 231eff52c0eSPeter Brune } 232eff52c0eSPeter Brune 233646217ecSPeter Brune #undef __FUNCT__ 234421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 235c79ef259SPeter Brune /*@ 236c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 237c79ef259SPeter Brune level of the FAS hierarchy. 238c79ef259SPeter Brune 239c79ef259SPeter Brune Input Parameters: 240c79ef259SPeter Brune + snes - the multigrid context 241c79ef259SPeter Brune level - the level to get 242c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 243c79ef259SPeter Brune 244c79ef259SPeter Brune Level: advanced 245c79ef259SPeter Brune 246c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 247c79ef259SPeter Brune 248c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 249c79ef259SPeter Brune @*/ 250421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 251421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 252421d9b32SPeter Brune PetscInt levels = fas->level; 253421d9b32SPeter Brune PetscInt i; 254421d9b32SPeter Brune PetscFunctionBegin; 255421d9b32SPeter Brune *lsnes = snes; 256421d9b32SPeter Brune if (fas->level < level) { 257421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 258421d9b32SPeter Brune } 259421d9b32SPeter Brune if (level > levels - 1) { 260421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 261421d9b32SPeter Brune } 262421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 263421d9b32SPeter Brune *lsnes = fas->next; 264421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 265421d9b32SPeter Brune } 266421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 267421d9b32SPeter Brune PetscFunctionReturn(0); 268421d9b32SPeter Brune } 269421d9b32SPeter Brune 270421d9b32SPeter Brune #undef __FUNCT__ 271421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 272c79ef259SPeter Brune /*@ 273c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 274c79ef259SPeter Brune Must be called before any other FAS routine. 275c79ef259SPeter Brune 276c79ef259SPeter Brune Input Parameters: 277c79ef259SPeter Brune + snes - the snes context 278c79ef259SPeter Brune . levels - the number of levels 279c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 280c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 281c79ef259SPeter Brune Fortran. 282c79ef259SPeter Brune 283c79ef259SPeter Brune Level: intermediate 284c79ef259SPeter Brune 285c79ef259SPeter Brune Notes: 286c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 287c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 288c79ef259SPeter Brune 289c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 290c79ef259SPeter Brune 291c79ef259SPeter Brune .seealso: SNESFASGetLevels() 292c79ef259SPeter Brune @*/ 293421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 294421d9b32SPeter Brune PetscErrorCode ierr; 295421d9b32SPeter Brune PetscInt i; 296421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 297421d9b32SPeter Brune MPI_Comm comm; 298421d9b32SPeter Brune PetscFunctionBegin; 299ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 300ee1fd11aSPeter Brune if (levels == fas->levels) { 301ee1fd11aSPeter Brune if (!comms) 302ee1fd11aSPeter Brune PetscFunctionReturn(0); 303ee1fd11aSPeter Brune } 304421d9b32SPeter Brune /* user has changed the number of levels; reset */ 305421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 306421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 307421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 308ee1fd11aSPeter Brune fas->next = PETSC_NULL; 309421d9b32SPeter Brune /* setup the finest level */ 310421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 311421d9b32SPeter Brune if (comms) comm = comms[i]; 312421d9b32SPeter Brune fas->level = i; 313421d9b32SPeter Brune fas->levels = levels; 314ee1fd11aSPeter Brune fas->next = PETSC_NULL; 315421d9b32SPeter Brune if (i > 0) { 316421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3174a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 318421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 319794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 320421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 321421d9b32SPeter Brune } 322421d9b32SPeter Brune } 323eff52c0eSPeter Brune /* by default set the GS smoothers up the chain if one exists on the finest level */ 324eff52c0eSPeter Brune if (snes->ops->computegs) 325eff52c0eSPeter Brune ierr = SNESFASSetGS(snes, snes->ops->computegs, snes->gsP, fas->useGS);CHKERRQ(ierr); 326421d9b32SPeter Brune PetscFunctionReturn(0); 327421d9b32SPeter Brune } 328421d9b32SPeter Brune 329421d9b32SPeter Brune #undef __FUNCT__ 330c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 331c79ef259SPeter Brune /*@ 332c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 333c79ef259SPeter Brune use on all levels. 334c79ef259SPeter Brune 335c79ef259SPeter Brune Logically Collective on PC 336c79ef259SPeter Brune 337c79ef259SPeter Brune Input Parameters: 338c79ef259SPeter Brune + snes - the multigrid context 339c79ef259SPeter Brune - n - the number of smoothing steps 340c79ef259SPeter Brune 341c79ef259SPeter Brune Options Database Key: 342c79ef259SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 343c79ef259SPeter Brune 344c79ef259SPeter Brune Level: advanced 345c79ef259SPeter Brune 346c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 347c79ef259SPeter Brune 348c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 349c79ef259SPeter Brune @*/ 350c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 351c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 352c79ef259SPeter Brune PetscErrorCode ierr = 0; 353c79ef259SPeter Brune PetscFunctionBegin; 354c79ef259SPeter Brune fas->max_up_it = n; 355c79ef259SPeter Brune if (fas->next) { 356c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 357c79ef259SPeter Brune } 358c79ef259SPeter Brune PetscFunctionReturn(0); 359c79ef259SPeter Brune } 360c79ef259SPeter Brune 361c79ef259SPeter Brune #undef __FUNCT__ 362c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 363c79ef259SPeter Brune /*@ 364c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 365c79ef259SPeter Brune use on all levels. 366c79ef259SPeter Brune 367c79ef259SPeter Brune Logically Collective on PC 368c79ef259SPeter Brune 369c79ef259SPeter Brune Input Parameters: 370c79ef259SPeter Brune + snes - the multigrid context 371c79ef259SPeter Brune - n - the number of smoothing steps 372c79ef259SPeter Brune 373c79ef259SPeter Brune Options Database Key: 374c79ef259SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 375c79ef259SPeter Brune 376c79ef259SPeter Brune Level: advanced 377c79ef259SPeter Brune 378c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 379c79ef259SPeter Brune 380c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 381c79ef259SPeter Brune @*/ 382c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 383c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 384c79ef259SPeter Brune PetscErrorCode ierr = 0; 385c79ef259SPeter Brune PetscFunctionBegin; 386c79ef259SPeter Brune fas->max_down_it = n; 387c79ef259SPeter Brune if (fas->next) { 388c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 389c79ef259SPeter Brune } 390c79ef259SPeter Brune PetscFunctionReturn(0); 391c79ef259SPeter Brune } 392c79ef259SPeter Brune 393c79ef259SPeter Brune #undef __FUNCT__ 394421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 395c79ef259SPeter Brune /*@ 396c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 397c79ef259SPeter Brune interpolation from l-1 to the lth level 398c79ef259SPeter Brune 399c79ef259SPeter Brune Input Parameters: 400c79ef259SPeter Brune + snes - the multigrid context 401c79ef259SPeter Brune . mat - the interpolation operator 402c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 403c79ef259SPeter Brune 404c79ef259SPeter Brune Level: advanced 405c79ef259SPeter Brune 406c79ef259SPeter Brune Notes: 407c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 408c79ef259SPeter Brune for the same level. 409c79ef259SPeter Brune 410c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 411c79ef259SPeter Brune out from the matrix size which one it is. 412c79ef259SPeter Brune 413c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 414c79ef259SPeter Brune 415c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 416c79ef259SPeter Brune @*/ 417421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 418421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 419d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 420421d9b32SPeter Brune 421421d9b32SPeter Brune PetscFunctionBegin; 422421d9b32SPeter Brune if (level > top_level) 423421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 424421d9b32SPeter Brune /* get to the correct level */ 425d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 426421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 427421d9b32SPeter Brune } 428421d9b32SPeter Brune if (fas->level != level) 429421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 430421d9b32SPeter Brune fas->interpolate = mat; 431421d9b32SPeter Brune PetscFunctionReturn(0); 432421d9b32SPeter Brune } 433421d9b32SPeter Brune 434421d9b32SPeter Brune #undef __FUNCT__ 435421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 436c79ef259SPeter Brune /*@ 437c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 438c79ef259SPeter Brune from level l to l-1. 439c79ef259SPeter Brune 440c79ef259SPeter Brune Input Parameters: 441c79ef259SPeter Brune + snes - the multigrid context 442c79ef259SPeter Brune . mat - the restriction matrix 443c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 444c79ef259SPeter Brune 445c79ef259SPeter Brune Level: advanced 446c79ef259SPeter Brune 447c79ef259SPeter Brune Notes: 448c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 449c79ef259SPeter Brune for the same level. 450c79ef259SPeter Brune 451c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 452c79ef259SPeter Brune out from the matrix size which one it is. 453c79ef259SPeter Brune 454c79ef259SPeter Brune If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 455c79ef259SPeter Brune is used. 456c79ef259SPeter Brune 457c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 458c79ef259SPeter Brune 459c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 460c79ef259SPeter Brune @*/ 461421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 462421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 463d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 464421d9b32SPeter Brune 465421d9b32SPeter Brune PetscFunctionBegin; 466421d9b32SPeter Brune if (level > top_level) 467421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 468421d9b32SPeter Brune /* get to the correct level */ 469d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 470421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 471421d9b32SPeter Brune } 472421d9b32SPeter Brune if (fas->level != level) 473421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 474421d9b32SPeter Brune fas->restrct = mat; 475421d9b32SPeter Brune PetscFunctionReturn(0); 476421d9b32SPeter Brune } 477421d9b32SPeter Brune 478efe1f98aSPeter Brune 479421d9b32SPeter Brune #undef __FUNCT__ 480421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 481c79ef259SPeter Brune /*@ 482c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 483c79ef259SPeter Brune operator from level l to l-1. 484c79ef259SPeter Brune 485c79ef259SPeter Brune Input Parameters: 486c79ef259SPeter Brune + snes - the multigrid context 487c79ef259SPeter Brune . rscale - the restriction scaling 488c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 489c79ef259SPeter Brune 490c79ef259SPeter Brune Level: advanced 491c79ef259SPeter Brune 492c79ef259SPeter Brune Notes: 493c79ef259SPeter Brune This is only used in the case that the injection is not set. 494c79ef259SPeter Brune 495c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 496c79ef259SPeter Brune 497c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 498c79ef259SPeter Brune @*/ 499421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 500421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 501d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 502421d9b32SPeter Brune 503421d9b32SPeter Brune PetscFunctionBegin; 504421d9b32SPeter Brune if (level > top_level) 505421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 506421d9b32SPeter Brune /* get to the correct level */ 507d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 508421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 509421d9b32SPeter Brune } 510421d9b32SPeter Brune if (fas->level != level) 511421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 512421d9b32SPeter Brune fas->rscale = rscale; 513421d9b32SPeter Brune PetscFunctionReturn(0); 514421d9b32SPeter Brune } 515421d9b32SPeter Brune 516efe1f98aSPeter Brune 517efe1f98aSPeter Brune #undef __FUNCT__ 518efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 519c79ef259SPeter Brune /*@ 520c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 521c79ef259SPeter Brune from level l to l-1. 522c79ef259SPeter Brune 523c79ef259SPeter Brune Input Parameters: 524c79ef259SPeter Brune + snes - the multigrid context 525c79ef259SPeter Brune . mat - the restriction matrix 526c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 527c79ef259SPeter Brune 528c79ef259SPeter Brune Level: advanced 529c79ef259SPeter Brune 530c79ef259SPeter Brune Notes: 531c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 532c79ef259SPeter Brune project the solution instead. 533c79ef259SPeter Brune 534c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 535c79ef259SPeter Brune 536c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 537c79ef259SPeter Brune @*/ 538efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 539efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 540efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 541efe1f98aSPeter Brune 542efe1f98aSPeter Brune PetscFunctionBegin; 543efe1f98aSPeter Brune if (level > top_level) 544efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 545efe1f98aSPeter Brune /* get to the correct level */ 546efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 547efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 548efe1f98aSPeter Brune } 549efe1f98aSPeter Brune if (fas->level != level) 550efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 551efe1f98aSPeter Brune fas->inject = mat; 552efe1f98aSPeter Brune PetscFunctionReturn(0); 553efe1f98aSPeter Brune } 554efe1f98aSPeter Brune 555421d9b32SPeter Brune #undef __FUNCT__ 556421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 557421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 558421d9b32SPeter Brune { 55977df8cc4SPeter Brune PetscErrorCode ierr = 0; 560421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 561421d9b32SPeter Brune 562421d9b32SPeter Brune PetscFunctionBegin; 563421d9b32SPeter Brune /* destroy local data created in SNESSetup_FAS */ 56451e86f29SPeter Brune #if 0 565293a7e31SPeter Brune /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 56651e86f29SPeter Brune #endif 567ee1fd11aSPeter Brune if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr); 56851e86f29SPeter Brune #if 0 56951e86f29SPeter Brune #endif 570421d9b32SPeter Brune PetscFunctionReturn(0); 571421d9b32SPeter Brune } 572421d9b32SPeter Brune 573421d9b32SPeter Brune #undef __FUNCT__ 574421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 575421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 576421d9b32SPeter Brune { 577421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 578421d9b32SPeter Brune PetscErrorCode ierr; 579421d9b32SPeter Brune 580421d9b32SPeter Brune PetscFunctionBegin; 581421d9b32SPeter Brune /* recursively resets and then destroys */ 582421d9b32SPeter Brune ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 58351e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 58451e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 585efe1f98aSPeter Brune if (fas->inject) { 586efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 587efe1f98aSPeter Brune } 58851e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 58951e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 59051e86f29SPeter Brune fas->restrct = PETSC_NULL; 59151e86f29SPeter Brune } else { 59251e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 59351e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 59451e86f29SPeter Brune } 59551e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 59651e86f29SPeter Brune if (snes->work) ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 597421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 598421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 599ee1fd11aSPeter Brune 600421d9b32SPeter Brune PetscFunctionReturn(0); 601421d9b32SPeter Brune } 602421d9b32SPeter Brune 603421d9b32SPeter Brune #undef __FUNCT__ 604421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 605421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 606421d9b32SPeter Brune { 6071a266240SBarry Smith SNES_FAS *fas = (SNES_FAS *) snes->data,*tmp; 608421d9b32SPeter Brune PetscErrorCode ierr; 609efe1f98aSPeter Brune VecScatter injscatter; 610d1adcc6fSPeter Brune PetscInt dm_levels; 611d1adcc6fSPeter Brune 612421d9b32SPeter Brune 613421d9b32SPeter Brune PetscFunctionBegin; 614eff52c0eSPeter Brune 615cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 616d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 617d1adcc6fSPeter Brune dm_levels++; 618cc05f883SPeter Brune if (dm_levels > fas->levels) { 619d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 620cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 621d1adcc6fSPeter Brune } 622d1adcc6fSPeter Brune } 623d1adcc6fSPeter Brune 624cc05f883SPeter Brune if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */ 625cc05f883SPeter Brune 626d1adcc6fSPeter Brune /* should call the SNESSetFromOptions() only when appropriate */ 6271a266240SBarry Smith tmp = fas; 6281a266240SBarry Smith while (tmp) { 629*6bed07a3SBarry Smith if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);} 630*6bed07a3SBarry Smith if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);} 6311a266240SBarry Smith tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0; 6321a266240SBarry Smith } 6331a266240SBarry Smith 634421d9b32SPeter Brune /* gets the solver ready for solution */ 635646217ecSPeter Brune if (fas->next) { 636646217ecSPeter Brune if (snes->ops->computegs) { 637646217ecSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 638646217ecSPeter Brune } 639646217ecSPeter Brune } 640421d9b32SPeter Brune if (snes->dm) { 641421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 642421d9b32SPeter Brune if (fas->next) { 643421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 644ee78dd50SPeter Brune if (!fas->next->dm) { 645ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 646ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 647ee78dd50SPeter Brune } 648421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 649ee78dd50SPeter Brune if (!fas->interpolate) { 650421d9b32SPeter Brune ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 651421d9b32SPeter Brune fas->restrct = fas->interpolate; 652421d9b32SPeter Brune } 653efe1f98aSPeter Brune /* set the injection from the DM */ 654efe1f98aSPeter Brune if (!fas->inject) { 655efe1f98aSPeter Brune ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 656efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 657efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 658efe1f98aSPeter Brune } 659421d9b32SPeter Brune } 660ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 661*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 662*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 663421d9b32SPeter Brune } 664cc05f883SPeter Brune 665fa9694d7SPeter Brune if (fas->next) { 666fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 667*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 668fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 669fa9694d7SPeter Brune } 670fa9694d7SPeter Brune /* got to set them all up at once */ 671421d9b32SPeter Brune PetscFunctionReturn(0); 672421d9b32SPeter Brune } 673421d9b32SPeter Brune 674421d9b32SPeter Brune #undef __FUNCT__ 675421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 676421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 677421d9b32SPeter Brune { 678ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 679ee78dd50SPeter Brune PetscInt levels = 1; 680d1adcc6fSPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 681421d9b32SPeter Brune PetscErrorCode ierr; 682ee78dd50SPeter Brune const char * def_smooth = SNESNRICHARDSON; 683ee78dd50SPeter Brune char pre_type[256]; 684ee78dd50SPeter Brune char post_type[256]; 685ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 686421d9b32SPeter Brune 687421d9b32SPeter Brune PetscFunctionBegin; 688c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 689ee78dd50SPeter Brune 690ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 691ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 692ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 693c732cbdbSBarry Smith if (!flg && snes->dm) { 694c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 695c732cbdbSBarry Smith levels++; 696d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 697c732cbdbSBarry Smith } 698ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 699ee78dd50SPeter Brune } 700ee78dd50SPeter Brune 701ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 702ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 703ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 704d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 705d1adcc6fSPeter Brune if (smoothflg) { 706ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 707ee78dd50SPeter Brune } else { 708d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 709d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 710ee78dd50SPeter Brune } 711ee78dd50SPeter Brune 712ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 713d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 714d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 715d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 716ee78dd50SPeter Brune 717646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 718ee78dd50SPeter Brune 719eff52c0eSPeter Brune ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 720eff52c0eSPeter Brune 721ee78dd50SPeter Brune /* other options for the coarsest level */ 722ee78dd50SPeter Brune if (fas->level == 0) { 723d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 724eff52c0eSPeter Brune ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 725ee78dd50SPeter Brune } 726ee78dd50SPeter Brune 727421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7288cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 729eff52c0eSPeter Brune 730eff52c0eSPeter Brune if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) { 7318cc86e31SPeter Brune const char *prefix; 7328cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 7338cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 7348cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 7358cc86e31SPeter Brune if (fas->level || (fas->levels == 1)) { 736eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 7378cc86e31SPeter Brune } else { 7388cc86e31SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 7398cc86e31SPeter Brune } 7408cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 7418cc86e31SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 7428cc86e31SPeter Brune if (snes->ops->computefunction) { 7438cc86e31SPeter Brune ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 7448cc86e31SPeter Brune } 7458cc86e31SPeter Brune } 7468cc86e31SPeter Brune 747eff52c0eSPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) { 74867339d5cSBarry Smith const char *prefix; 74967339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 750ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 75167339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 752eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 753293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 754ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 75569bed9afSBarry Smith if (snes->ops->computefunction) { 75669bed9afSBarry Smith ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 75769bed9afSBarry Smith } 758ee78dd50SPeter Brune } 7591a266240SBarry Smith if (fas->upsmooth) { 7601a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 7611a266240SBarry Smith } 7621a266240SBarry Smith 7631a266240SBarry Smith if (fas->downsmooth) { 7641a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 7651a266240SBarry Smith } 766ee78dd50SPeter Brune 767ee78dd50SPeter Brune if (monflg) { 768646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 769794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 770293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 771293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 772ee78dd50SPeter Brune } 773ee78dd50SPeter Brune 774ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 7751a266240SBarry Smith if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 776421d9b32SPeter Brune PetscFunctionReturn(0); 777421d9b32SPeter Brune } 778421d9b32SPeter Brune 779421d9b32SPeter Brune #undef __FUNCT__ 780421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 781421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 782421d9b32SPeter Brune { 783421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 784421d9b32SPeter Brune PetscBool iascii; 785421d9b32SPeter Brune PetscErrorCode ierr; 786421d9b32SPeter Brune PetscInt levels = fas->levels; 787421d9b32SPeter Brune PetscInt i; 788421d9b32SPeter Brune 789421d9b32SPeter Brune PetscFunctionBegin; 790421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 791421d9b32SPeter Brune if (iascii) { 792421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 793421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 794421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 795421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 796ee78dd50SPeter Brune if (fas->upsmooth) { 797ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 798421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 799ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 800421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 801421d9b32SPeter Brune } else { 802ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 803421d9b32SPeter Brune } 804ee78dd50SPeter Brune if (fas->downsmooth) { 805ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 806421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 807ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 808421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 809421d9b32SPeter Brune } else { 810ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 811421d9b32SPeter Brune } 812eff52c0eSPeter Brune if (fas->useGS) { 813eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 814eff52c0eSPeter Brune } 815421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 816421d9b32SPeter Brune } 817421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 818421d9b32SPeter Brune } else { 819421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 820421d9b32SPeter Brune } 821421d9b32SPeter Brune PetscFunctionReturn(0); 822421d9b32SPeter Brune } 823421d9b32SPeter Brune 824421d9b32SPeter Brune /* 825fa9694d7SPeter Brune 826fa9694d7SPeter Brune Defines the FAS cycle as: 827fa9694d7SPeter Brune 828fa9694d7SPeter Brune fine problem: F(x) = 0 829fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 830fa9694d7SPeter Brune 831fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 832fa9694d7SPeter Brune 833fa9694d7SPeter Brune correction: 834fa9694d7SPeter Brune 835fa9694d7SPeter Brune x = x + I(x^c - Rx) 836fa9694d7SPeter Brune 837421d9b32SPeter Brune */ 838fa9694d7SPeter Brune 839421d9b32SPeter Brune #undef __FUNCT__ 840421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 841646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 842fa9694d7SPeter Brune 843fa9694d7SPeter Brune PetscErrorCode ierr; 844646217ecSPeter Brune Vec X_c, Xo_c, F_c, B_c,F,B; 845fa9694d7SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 846646217ecSPeter Brune PetscInt i, k; 8472d15c84fSPeter Brune PetscReal fnorm; 848421d9b32SPeter Brune 849421d9b32SPeter Brune PetscFunctionBegin; 850f5a6d4f9SBarry Smith F = snes->vec_func; 851646217ecSPeter Brune B = snes->vec_rhs; 852fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 853d1adcc6fSPeter Brune if (fas->downsmooth) { 854d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 855d1adcc6fSPeter Brune } else if (snes->ops->computegs) { 856794bee33SPeter Brune if (fas->monitor) { 857794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 858794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 859d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 860eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 861d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 862794bee33SPeter Brune } 863646217ecSPeter Brune for (k = 0; k < fas->max_up_it; k++) { 864646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 865cc05f883SPeter Brune if (fas->monitor) { 866794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 867794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 868d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 869eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 870d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 871794bee33SPeter Brune } 872646217ecSPeter Brune } 873c90fad12SPeter Brune } else if (snes->pc) { 874c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 875fe6f9142SPeter Brune } 876fa9694d7SPeter Brune if (fas->next) { 877293a7e31SPeter Brune for (i = 0; i < fas->n_cycles; i++) { 878ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 879794bee33SPeter Brune 880c90fad12SPeter Brune X_c = fas->next->vec_sol; 881293a7e31SPeter Brune Xo_c = fas->next->work[0]; 882c90fad12SPeter Brune F_c = fas->next->vec_func; 883293a7e31SPeter Brune B_c = fas->next->work[1]; 884efe1f98aSPeter Brune 885efe1f98aSPeter Brune /* inject the solution */ 886efe1f98aSPeter Brune if (fas->inject) { 887a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 888efe1f98aSPeter Brune } else { 889a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 890a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 891efe1f98aSPeter Brune } 892293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 893293a7e31SPeter Brune 894293a7e31SPeter Brune /* restrict the defect */ 895293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 896293a7e31SPeter Brune 897c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 898ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 899c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 900293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 901c90fad12SPeter Brune 902ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 903ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 904c90fad12SPeter Brune 905c90fad12SPeter Brune /* recurse to the next level */ 906f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 907646217ecSPeter Brune ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); 908f5a6d4f9SBarry Smith fas->next->vec_rhs = PETSC_NULL; 909ee78dd50SPeter Brune 910fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 911fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 912ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 913293a7e31SPeter Brune } 914fa9694d7SPeter Brune } 915d1adcc6fSPeter Brune /* up-smooth -- just update using the down-smoother */ 916c90fad12SPeter Brune if (fas->level != 0) { 917d1adcc6fSPeter Brune if (fas->upsmooth) { 918d1adcc6fSPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 919d1adcc6fSPeter Brune } else if (snes->ops->computegs) { 920794bee33SPeter Brune if (fas->monitor) { 921794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 922794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 923d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 924eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 925d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 926794bee33SPeter Brune } 927646217ecSPeter Brune for (k = 0; k < fas->max_down_it; k++) { 928646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 929794bee33SPeter Brune if (fas->monitor) { 930794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 931794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 932d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 933eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 934d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 935794bee33SPeter Brune } 936646217ecSPeter Brune } 937c90fad12SPeter Brune } else if (snes->pc) { 938c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 939c90fad12SPeter Brune } 940fe6f9142SPeter Brune } 941fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 942fa9694d7SPeter Brune 943fa9694d7SPeter Brune PetscFunctionReturn(0); 944421d9b32SPeter Brune } 945421d9b32SPeter Brune 946421d9b32SPeter Brune #undef __FUNCT__ 947421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 948421d9b32SPeter Brune 949421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 950421d9b32SPeter Brune { 951fa9694d7SPeter Brune PetscErrorCode ierr; 952fe6f9142SPeter Brune PetscInt i, maxits; 953f5a6d4f9SBarry Smith Vec X, B,F; 954fe6f9142SPeter Brune PetscReal fnorm; 955421d9b32SPeter Brune PetscFunctionBegin; 956fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 957fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 958fa9694d7SPeter Brune X = snes->vec_sol; 959fe6f9142SPeter Brune B = snes->vec_rhs; 960f5a6d4f9SBarry Smith F = snes->vec_func; 961293a7e31SPeter Brune 962293a7e31SPeter Brune /*norm setup */ 963fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 964fe6f9142SPeter Brune snes->iter = 0; 965fe6f9142SPeter Brune snes->norm = 0.; 966fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 967fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 968fe6f9142SPeter Brune if (snes->domainerror) { 969fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 970fe6f9142SPeter Brune PetscFunctionReturn(0); 971fe6f9142SPeter Brune } 972fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 973fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 974fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 975fe6f9142SPeter Brune snes->norm = fnorm; 976fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 977fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 978fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 979fe6f9142SPeter Brune 980fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 981fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 982fe6f9142SPeter Brune /* test convergence */ 983fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 984fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 985fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 986fe6f9142SPeter Brune /* Call general purpose update function */ 987646217ecSPeter Brune 988fe6f9142SPeter Brune if (snes->ops->update) { 989fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 990fe6f9142SPeter Brune } 991646217ecSPeter Brune ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 992c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 993c90fad12SPeter Brune /* Monitor convergence */ 994c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 995c90fad12SPeter Brune snes->iter = i+1; 996c90fad12SPeter Brune snes->norm = fnorm; 997c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 998c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 999c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1000c90fad12SPeter Brune /* Test for convergence */ 1001c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1002c90fad12SPeter Brune if (snes->reason) break; 1003fe6f9142SPeter Brune } 1004fe6f9142SPeter Brune if (i == maxits) { 1005fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1006fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1007fe6f9142SPeter Brune } 1008421d9b32SPeter Brune PetscFunctionReturn(0); 1009421d9b32SPeter Brune } 1010