1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h" I*/ 26cd56b8bSPeter Brune 3b6266c6eSPeter Brune #undef __FUNCT__ 4be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSSetTolerances" 5b6266c6eSPeter Brune /*@ 6be95d8f1SBarry Smith SNESNGSSetTolerances - Sets various parameters used in convergence tests. 7b6266c6eSPeter Brune 8b6266c6eSPeter Brune Logically Collective on SNES 9b6266c6eSPeter Brune 10b6266c6eSPeter Brune Input Parameters: 11b6266c6eSPeter Brune + snes - the SNES context 12b6266c6eSPeter Brune . abstol - absolute convergence tolerance 13b6266c6eSPeter Brune . rtol - relative convergence tolerance 14b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 15b6266c6eSPeter Brune - maxit - maximum number of iterations 16b6266c6eSPeter Brune 17b6266c6eSPeter Brune 18b6266c6eSPeter Brune Options Database Keys: 19be95d8f1SBarry Smith + -snes_ngs_atol <abstol> - Sets abstol 20be95d8f1SBarry Smith . -snes_ngs_rtol <rtol> - Sets rtol 21be95d8f1SBarry Smith . -snes_ngs_stol <stol> - Sets stol 22b6266c6eSPeter Brune - -snes_max_it <maxit> - Sets maxit 23b6266c6eSPeter Brune 24b6266c6eSPeter Brune Level: intermediate 25b6266c6eSPeter Brune 26b6266c6eSPeter Brune .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances 27b6266c6eSPeter Brune 28b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance() 29b6266c6eSPeter Brune @*/ 30be95d8f1SBarry Smith PetscErrorCode SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit) 31b6266c6eSPeter Brune { 32be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 33b6266c6eSPeter Brune 34b6266c6eSPeter Brune PetscFunctionBegin; 35b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 36b6266c6eSPeter Brune 37b6266c6eSPeter Brune if (abstol != PETSC_DEFAULT) { 3857622a8eSBarry Smith if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol); 39b6266c6eSPeter Brune gs->abstol = abstol; 40b6266c6eSPeter Brune } 41b6266c6eSPeter Brune if (rtol != PETSC_DEFAULT) { 4257622a8eSBarry Smith if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol); 43b6266c6eSPeter Brune gs->rtol = rtol; 44b6266c6eSPeter Brune } 45b6266c6eSPeter Brune if (stol != PETSC_DEFAULT) { 4657622a8eSBarry Smith if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol); 47b6266c6eSPeter Brune gs->stol = stol; 48b6266c6eSPeter Brune } 49b6266c6eSPeter Brune if (maxit != PETSC_DEFAULT) { 50ce94432eSBarry Smith if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit); 51b6266c6eSPeter Brune gs->max_its = maxit; 52b6266c6eSPeter Brune } 53b6266c6eSPeter Brune PetscFunctionReturn(0); 54b6266c6eSPeter Brune } 55b6266c6eSPeter Brune 56b6266c6eSPeter Brune #undef __FUNCT__ 57be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSGetTolerances" 58b6266c6eSPeter Brune /*@ 59be95d8f1SBarry Smith SNESNGSGetTolerances - Gets various parameters used in convergence tests. 60b6266c6eSPeter Brune 61b6266c6eSPeter Brune Not Collective 62b6266c6eSPeter Brune 63b6266c6eSPeter Brune Input Parameters: 64b6266c6eSPeter Brune + snes - the SNES context 65b6266c6eSPeter Brune . atol - absolute convergence tolerance 66b6266c6eSPeter Brune . rtol - relative convergence tolerance 67b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm 68b6266c6eSPeter Brune of the change in the solution between steps 69b6266c6eSPeter Brune - maxit - maximum number of iterations 70b6266c6eSPeter Brune 71b6266c6eSPeter Brune Notes: 720298fd71SBarry Smith The user can specify NULL for any parameter that is not needed. 73b6266c6eSPeter Brune 74b6266c6eSPeter Brune Level: intermediate 75b6266c6eSPeter Brune 76b6266c6eSPeter Brune .keywords: SNES, nonlinear, get, convergence, tolerances 77b6266c6eSPeter Brune 78b6266c6eSPeter Brune .seealso: SNESSetTolerances() 79b6266c6eSPeter Brune @*/ 80be95d8f1SBarry Smith PetscErrorCode SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit) 81b6266c6eSPeter Brune { 82be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 83b6266c6eSPeter Brune 84b6266c6eSPeter Brune PetscFunctionBegin; 85b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 86b6266c6eSPeter Brune if (atol) *atol = gs->abstol; 87b6266c6eSPeter Brune if (rtol) *rtol = gs->rtol; 88b6266c6eSPeter Brune if (stol) *stol = gs->stol; 89b6266c6eSPeter Brune if (maxit) *maxit = gs->max_its; 90b6266c6eSPeter Brune PetscFunctionReturn(0); 91b6266c6eSPeter Brune } 926cd56b8bSPeter Brune 936cd56b8bSPeter Brune #undef __FUNCT__ 94be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSSetSweeps" 956cd56b8bSPeter Brune /*@ 96be95d8f1SBarry Smith SNESNGSSetSweeps - Sets the number of sweeps of GS to use. 976cd56b8bSPeter Brune 986cd56b8bSPeter Brune Input Parameters: 996cd56b8bSPeter Brune + snes - the SNES context 1006cd56b8bSPeter Brune - sweeps - the number of sweeps of GS to perform. 1016cd56b8bSPeter Brune 1026cd56b8bSPeter Brune Level: intermediate 1036cd56b8bSPeter Brune 1046cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 1056cd56b8bSPeter Brune 106be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSGetSweeps() 1076cd56b8bSPeter Brune @*/ 1086cd56b8bSPeter Brune 109be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) 1106cd56b8bSPeter Brune { 111be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 1126cd56b8bSPeter Brune 1136cd56b8bSPeter Brune PetscFunctionBegin; 1146cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1156cd56b8bSPeter Brune gs->sweeps = sweeps; 1166cd56b8bSPeter Brune PetscFunctionReturn(0); 1176cd56b8bSPeter Brune } 1186cd56b8bSPeter Brune 1196cd56b8bSPeter Brune #undef __FUNCT__ 120be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSGetSweeps" 1216cd56b8bSPeter Brune /*@ 122be95d8f1SBarry Smith SNESNGSGetSweeps - Gets the number of sweeps GS will use. 1236cd56b8bSPeter Brune 1246cd56b8bSPeter Brune Input Parameters: 1256cd56b8bSPeter Brune . snes - the SNES context 1266cd56b8bSPeter Brune 1276cd56b8bSPeter Brune Output Parameters: 1286cd56b8bSPeter Brune . sweeps - the number of sweeps of GS to perform. 1296cd56b8bSPeter Brune 1306cd56b8bSPeter Brune Level: intermediate 1316cd56b8bSPeter Brune 1326cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 1336cd56b8bSPeter Brune 134be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSSetSweeps() 1356cd56b8bSPeter Brune @*/ 136be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps) 1376cd56b8bSPeter Brune { 138be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 1396cd56b8bSPeter Brune 1406cd56b8bSPeter Brune PetscFunctionBegin; 1416cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1426cd56b8bSPeter Brune *sweeps = gs->sweeps; 1436cd56b8bSPeter Brune PetscFunctionReturn(0); 1446cd56b8bSPeter Brune } 1456cd56b8bSPeter Brune 146c79d6ae8SPeter Brune #undef __FUNCT__ 147be95d8f1SBarry Smith #define __FUNCT__ "SNESReset_NGS" 148be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes) 1493542a6bcSPeter Brune { 1503542a6bcSPeter Brune PetscFunctionBegin; 1513542a6bcSPeter Brune PetscFunctionReturn(0); 1523542a6bcSPeter Brune } 1533542a6bcSPeter Brune 1543542a6bcSPeter Brune #undef __FUNCT__ 155be95d8f1SBarry Smith #define __FUNCT__ "SNESDestroy_NGS" 156be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes) 1573542a6bcSPeter Brune { 1583542a6bcSPeter Brune PetscErrorCode ierr; 1593542a6bcSPeter Brune 1603542a6bcSPeter Brune PetscFunctionBegin; 161be95d8f1SBarry Smith ierr = SNESReset_NGS(snes);CHKERRQ(ierr); 16222d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 1633542a6bcSPeter Brune PetscFunctionReturn(0); 1643542a6bcSPeter Brune } 1653542a6bcSPeter Brune 1663542a6bcSPeter Brune #undef __FUNCT__ 167be95d8f1SBarry Smith #define __FUNCT__ "SNESSetUp_NGS" 168be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes) 1693542a6bcSPeter Brune { 170ef107cf5SPeter Brune PetscErrorCode ierr; 171ef107cf5SPeter Brune PetscErrorCode (*f)(SNES,Vec,Vec,void*); 172ef107cf5SPeter Brune 1733542a6bcSPeter Brune PetscFunctionBegin; 174be95d8f1SBarry Smith ierr = SNESGetNGS(snes,&f,NULL);CHKERRQ(ierr); 175ef107cf5SPeter Brune if (!f) { 176be95d8f1SBarry Smith ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr); 177ef107cf5SPeter Brune } 1783542a6bcSPeter Brune PetscFunctionReturn(0); 1793542a6bcSPeter Brune } 1803542a6bcSPeter Brune 1813542a6bcSPeter Brune #undef __FUNCT__ 182be95d8f1SBarry Smith #define __FUNCT__ "SNESSetFromOptions_NGS" 1834416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes) 1843542a6bcSPeter Brune { 185be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 1863542a6bcSPeter Brune PetscErrorCode ierr; 187b6266c6eSPeter Brune PetscInt sweeps,max_its=PETSC_DEFAULT; 188b6266c6eSPeter Brune PetscReal rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT; 189b6266c6eSPeter Brune PetscBool flg,flg1,flg2,flg3; 1903542a6bcSPeter Brune 1913542a6bcSPeter Brune PetscFunctionBegin; 192e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNES GS options");CHKERRQ(ierr); 1936cd56b8bSPeter Brune /* GS Options */ 194be95d8f1SBarry Smith ierr = PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr); 195b6266c6eSPeter Brune if (flg) { 196be95d8f1SBarry Smith ierr = SNESNGSSetSweeps(snes,sweeps);CHKERRQ(ierr); 197b6266c6eSPeter Brune } 19851b6f3e6SMatthew G. Knepley ierr = PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr); 19951b6f3e6SMatthew G. Knepley ierr = PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr); 20051b6f3e6SMatthew G. Knepley ierr = PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr); 20151b6f3e6SMatthew G. Knepley ierr = PetscOptionsInt("-snes_ngs_max_it","Maximum number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr); 202b6266c6eSPeter Brune if (flg || flg1 || flg2 || flg3) { 203be95d8f1SBarry Smith ierr = SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr); 204b6266c6eSPeter Brune } 205ef107cf5SPeter Brune flg = PETSC_FALSE; 206*98272f38SBarry Smith ierr = PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL);CHKERRQ(ierr); 207ef107cf5SPeter Brune if (flg) { 208be95d8f1SBarry Smith ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr); 209*98272f38SBarry Smith ierr = PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n");CHKERRQ(ierr); 210ef107cf5SPeter Brune } 211be95d8f1SBarry Smith ierr = PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr); 212*98272f38SBarry Smith ierr = PetscOptionsBool("-snes_ngs_secant_mat_coloring","Use the graph coloring of the Jacobian for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg);CHKERRQ(ierr); 213ef107cf5SPeter Brune 2144b32a720SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2153542a6bcSPeter Brune PetscFunctionReturn(0); 2163542a6bcSPeter Brune } 2173542a6bcSPeter Brune 2183542a6bcSPeter Brune #undef __FUNCT__ 219be95d8f1SBarry Smith #define __FUNCT__ "SNESView_NGS" 220be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) 2213542a6bcSPeter Brune { 222*98272f38SBarry Smith PetscErrorCode ierr; 223*98272f38SBarry Smith PetscErrorCode (*f)(SNES,Vec,Vec,void*); 224*98272f38SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 225*98272f38SBarry Smith PetscBool iascii; 226*98272f38SBarry Smith 2273542a6bcSPeter Brune PetscFunctionBegin; 228*98272f38SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 229*98272f38SBarry Smith if (iascii) { 230*98272f38SBarry Smith ierr = DMSNESGetNGS(snes->dm,&f,NULL);CHKERRQ(ierr); 231*98272f38SBarry Smith if (f == SNESComputeNGSDefaultSecant) { 232*98272f38SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," NGS: Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h);CHKERRQ(ierr); 233*98272f38SBarry Smith } 234*98272f38SBarry Smith } 2353542a6bcSPeter Brune PetscFunctionReturn(0); 2363542a6bcSPeter Brune } 2373542a6bcSPeter Brune 2383542a6bcSPeter Brune #undef __FUNCT__ 239be95d8f1SBarry Smith #define __FUNCT__ "SNESSolve_NGS" 240be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes) 2413542a6bcSPeter Brune { 2423542a6bcSPeter Brune Vec F; 2433542a6bcSPeter Brune Vec X; 2443542a6bcSPeter Brune Vec B; 2453542a6bcSPeter Brune PetscInt i; 24606e07b1aSPeter Brune PetscReal fnorm; 2473542a6bcSPeter Brune PetscErrorCode ierr; 248365a6726SPeter Brune SNESNormSchedule normschedule; 2493542a6bcSPeter Brune 2503542a6bcSPeter Brune PetscFunctionBegin; 251c579b300SPatrick Farrell 2526c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 253c579b300SPatrick Farrell 254fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 2553542a6bcSPeter Brune X = snes->vec_sol; 2563542a6bcSPeter Brune F = snes->vec_func; 2573542a6bcSPeter Brune B = snes->vec_rhs; 2583542a6bcSPeter Brune 259e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 260534ebe21SPeter Brune snes->iter = 0; 261534ebe21SPeter Brune snes->norm = 0.; 262e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 263526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 264534ebe21SPeter Brune 265365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 266365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 2673542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 268e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2693542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2701aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 271e4ed7901SPeter Brune 2723542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 273422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 274e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 275e4ed7901SPeter Brune snes->iter = 0; 2763542a6bcSPeter Brune snes->norm = fnorm; 277e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 278a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 279e4ed7901SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 2803542a6bcSPeter Brune 2813542a6bcSPeter Brune /* test convergence */ 2823542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2833542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 284534ebe21SPeter Brune } else { 285e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 286a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 287534ebe21SPeter Brune } 288534ebe21SPeter Brune 2893542a6bcSPeter Brune /* Call general purpose update function */ 2903542a6bcSPeter Brune if (snes->ops->update) { 2913542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 2923542a6bcSPeter Brune } 293534ebe21SPeter Brune 2943542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 295be95d8f1SBarry Smith ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr); 2964b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 297365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 2983542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2993542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 300422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 3013542a6bcSPeter Brune /* Monitor convergence */ 302e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 3033542a6bcSPeter Brune snes->iter = i+1; 3043542a6bcSPeter Brune snes->norm = fnorm; 305e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 306a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 3073542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 308ac4a6415SMatthew G. Knepley } 3093542a6bcSPeter Brune /* Test for convergence */ 310365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3113542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 3123542a6bcSPeter Brune /* Call general purpose update function */ 3133542a6bcSPeter Brune if (snes->ops->update) { 3143542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 3153542a6bcSPeter Brune } 3163542a6bcSPeter Brune } 317365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 318534ebe21SPeter Brune if (i == snes->max_its) { 319534ebe21SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 320534ebe21SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 321534ebe21SPeter Brune } 3221aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 3233542a6bcSPeter Brune PetscFunctionReturn(0); 3243542a6bcSPeter Brune } 3253542a6bcSPeter Brune 3263542a6bcSPeter Brune /*MC 327*98272f38SBarry Smith SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation 328*98272f38SBarry Smith using coloring. 3293542a6bcSPeter Brune 3303542a6bcSPeter Brune Level: advanced 3313542a6bcSPeter Brune 332*98272f38SBarry Smith Options Database: 333*98272f38SBarry Smith + -snes_ngs_sweeps - Number of sweeps of GS to apply 334*98272f38SBarry Smith . -snes_ngs_atol - Absolute residual tolerance for GS iteration 335*98272f38SBarry Smith . -snes_ngs_rtol - Relative residual tolerance for GS iteration 336*98272f38SBarry Smith . -snes_ngs_stol - Absolute update tolerance for GS iteration 337*98272f38SBarry Smith . -snes_ngs_max_it - Maximum number of sweeps of GS to apply 338*98272f38SBarry Smith . -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is 339*98272f38SBarry Smith used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring 340*98272f38SBarry Smith is available or a Jacobian sparse matrix is provided (from which to get the coloring). 341*98272f38SBarry Smith . -snes_ngs_secant_h - Differencing parameter for secant approximation 342*98272f38SBarry Smith - -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available. 343*98272f38SBarry Smith 344*98272f38SBarry Smith 3453542a6bcSPeter Brune Notes: 3463542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 3473542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 3483542a6bcSPeter Brune 3494f02bc6aSBarry Smith References: 35096a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 3514f02bc6aSBarry Smith SIAM Review, 57(4), 2015 3524f02bc6aSBarry Smith 353*98272f38SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances() 3543542a6bcSPeter Brune M*/ 3553542a6bcSPeter Brune 3563542a6bcSPeter Brune #undef __FUNCT__ 357be95d8f1SBarry Smith #define __FUNCT__ "SNESCreate_NGS" 358be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) 3593542a6bcSPeter Brune { 360be95d8f1SBarry Smith SNES_NGS *gs; 3613542a6bcSPeter Brune PetscErrorCode ierr; 3623542a6bcSPeter Brune 3633542a6bcSPeter Brune PetscFunctionBegin; 364be95d8f1SBarry Smith snes->ops->destroy = SNESDestroy_NGS; 365be95d8f1SBarry Smith snes->ops->setup = SNESSetUp_NGS; 366be95d8f1SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NGS; 367be95d8f1SBarry Smith snes->ops->view = SNESView_NGS; 368be95d8f1SBarry Smith snes->ops->solve = SNESSolve_NGS; 369be95d8f1SBarry Smith snes->ops->reset = SNESReset_NGS; 3703542a6bcSPeter Brune 3713542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 3723542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 3733542a6bcSPeter Brune 37488976e71SPeter Brune if (!snes->tolerancesset) { 3750e444f03SPeter Brune snes->max_its = 10000; 3760e444f03SPeter Brune snes->max_funcs = 10000; 37788976e71SPeter Brune } 3780e444f03SPeter Brune 379b00a9115SJed Brown ierr = PetscNewLog(snes,&gs);CHKERRQ(ierr); 3806cd56b8bSPeter Brune 3816cd56b8bSPeter Brune gs->sweeps = 1; 382b6266c6eSPeter Brune gs->rtol = 1e-5; 383b6266c6eSPeter Brune gs->abstol = 1e-15; 384b6266c6eSPeter Brune gs->stol = 1e-12; 385b6266c6eSPeter Brune gs->max_its = 50; 386737a7e12SPeter Brune gs->h = 1e-8; 3876cd56b8bSPeter Brune 3883542a6bcSPeter Brune snes->data = (void*) gs; 3893542a6bcSPeter Brune PetscFunctionReturn(0); 3903542a6bcSPeter Brune } 391