1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h" I*/ 26cd56b8bSPeter Brune 3b6266c6eSPeter Brune /*@ 4be95d8f1SBarry Smith SNESNGSSetTolerances - Sets various parameters used in convergence tests. 5b6266c6eSPeter Brune 6b6266c6eSPeter Brune Logically Collective on SNES 7b6266c6eSPeter Brune 8b6266c6eSPeter Brune Input Parameters: 9b6266c6eSPeter Brune + snes - the SNES context 10b6266c6eSPeter Brune . abstol - absolute convergence tolerance 11b6266c6eSPeter Brune . rtol - relative convergence tolerance 12b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 13b6266c6eSPeter Brune - maxit - maximum number of iterations 14b6266c6eSPeter Brune 15b6266c6eSPeter Brune Options Database Keys: 16be95d8f1SBarry Smith + -snes_ngs_atol <abstol> - Sets abstol 17be95d8f1SBarry Smith . -snes_ngs_rtol <rtol> - Sets rtol 18be95d8f1SBarry Smith . -snes_ngs_stol <stol> - Sets stol 19b6266c6eSPeter Brune - -snes_max_it <maxit> - Sets maxit 20b6266c6eSPeter Brune 21b6266c6eSPeter Brune Level: intermediate 22b6266c6eSPeter Brune 23db781477SPatrick Sanan .seealso: `SNESSetTrustRegionTolerance()` 24b6266c6eSPeter Brune @*/ 25be95d8f1SBarry Smith PetscErrorCode SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit) 26b6266c6eSPeter Brune { 27be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 28b6266c6eSPeter Brune 29b6266c6eSPeter Brune PetscFunctionBegin; 30b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31b6266c6eSPeter Brune 32b6266c6eSPeter Brune if (abstol != PETSC_DEFAULT) { 3308401ef6SPierre Jolivet PetscCheck(abstol >= 0.0,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol); 34b6266c6eSPeter Brune gs->abstol = abstol; 35b6266c6eSPeter Brune } 36b6266c6eSPeter Brune if (rtol != PETSC_DEFAULT) { 370b121fc5SBarry Smith PetscCheck(rtol >= 0.0 && rtol < 1.0,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol); 38b6266c6eSPeter Brune gs->rtol = rtol; 39b6266c6eSPeter Brune } 40b6266c6eSPeter Brune if (stol != PETSC_DEFAULT) { 4108401ef6SPierre Jolivet PetscCheck(stol >= 0.0,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol); 42b6266c6eSPeter Brune gs->stol = stol; 43b6266c6eSPeter Brune } 44b6266c6eSPeter Brune if (maxit != PETSC_DEFAULT) { 4563a3b9bcSJacob Faibussowitsch PetscCheck(maxit >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %" PetscInt_FMT " must be non-negative",maxit); 46b6266c6eSPeter Brune gs->max_its = maxit; 47b6266c6eSPeter Brune } 48b6266c6eSPeter Brune PetscFunctionReturn(0); 49b6266c6eSPeter Brune } 50b6266c6eSPeter Brune 51b6266c6eSPeter Brune /*@ 52be95d8f1SBarry Smith SNESNGSGetTolerances - Gets various parameters used in convergence tests. 53b6266c6eSPeter Brune 54b6266c6eSPeter Brune Not Collective 55b6266c6eSPeter Brune 56b6266c6eSPeter Brune Input Parameters: 57b6266c6eSPeter Brune + snes - the SNES context 58b6266c6eSPeter Brune . atol - absolute convergence tolerance 59b6266c6eSPeter Brune . rtol - relative convergence tolerance 60b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm 61b6266c6eSPeter Brune of the change in the solution between steps 62b6266c6eSPeter Brune - maxit - maximum number of iterations 63b6266c6eSPeter Brune 64b6266c6eSPeter Brune Notes: 650298fd71SBarry Smith The user can specify NULL for any parameter that is not needed. 66b6266c6eSPeter Brune 67b6266c6eSPeter Brune Level: intermediate 68b6266c6eSPeter Brune 69db781477SPatrick Sanan .seealso: `SNESSetTolerances()` 70b6266c6eSPeter Brune @*/ 71be95d8f1SBarry Smith PetscErrorCode SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit) 72b6266c6eSPeter Brune { 73be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 74b6266c6eSPeter Brune 75b6266c6eSPeter Brune PetscFunctionBegin; 76b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 77b6266c6eSPeter Brune if (atol) *atol = gs->abstol; 78b6266c6eSPeter Brune if (rtol) *rtol = gs->rtol; 79b6266c6eSPeter Brune if (stol) *stol = gs->stol; 80b6266c6eSPeter Brune if (maxit) *maxit = gs->max_its; 81b6266c6eSPeter Brune PetscFunctionReturn(0); 82b6266c6eSPeter Brune } 836cd56b8bSPeter Brune 846cd56b8bSPeter Brune /*@ 85be95d8f1SBarry Smith SNESNGSSetSweeps - Sets the number of sweeps of GS to use. 866cd56b8bSPeter Brune 876cd56b8bSPeter Brune Input Parameters: 886cd56b8bSPeter Brune + snes - the SNES context 896cd56b8bSPeter Brune - sweeps - the number of sweeps of GS to perform. 906cd56b8bSPeter Brune 916cd56b8bSPeter Brune Level: intermediate 926cd56b8bSPeter Brune 93db781477SPatrick Sanan .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()` 946cd56b8bSPeter Brune @*/ 956cd56b8bSPeter Brune 96be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) 976cd56b8bSPeter Brune { 98be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 996cd56b8bSPeter Brune 1006cd56b8bSPeter Brune PetscFunctionBegin; 1016cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1026cd56b8bSPeter Brune gs->sweeps = sweeps; 1036cd56b8bSPeter Brune PetscFunctionReturn(0); 1046cd56b8bSPeter Brune } 1056cd56b8bSPeter Brune 1066cd56b8bSPeter Brune /*@ 107be95d8f1SBarry Smith SNESNGSGetSweeps - Gets the number of sweeps GS will use. 1086cd56b8bSPeter Brune 1096cd56b8bSPeter Brune Input Parameters: 1106cd56b8bSPeter Brune . snes - the SNES context 1116cd56b8bSPeter Brune 1126cd56b8bSPeter Brune Output Parameters: 1136cd56b8bSPeter Brune . sweeps - the number of sweeps of GS to perform. 1146cd56b8bSPeter Brune 1156cd56b8bSPeter Brune Level: intermediate 1166cd56b8bSPeter Brune 117db781477SPatrick Sanan .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()` 1186cd56b8bSPeter Brune @*/ 119be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps) 1206cd56b8bSPeter Brune { 121be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 1226cd56b8bSPeter Brune 1236cd56b8bSPeter Brune PetscFunctionBegin; 1246cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1256cd56b8bSPeter Brune *sweeps = gs->sweeps; 1266cd56b8bSPeter Brune PetscFunctionReturn(0); 1276cd56b8bSPeter Brune } 1286cd56b8bSPeter Brune 129be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes) 1303542a6bcSPeter Brune { 1318a86d3c5SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 1328a86d3c5SBarry Smith 1333542a6bcSPeter Brune PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&gs->coloring)); 1353542a6bcSPeter Brune PetscFunctionReturn(0); 1363542a6bcSPeter Brune } 1373542a6bcSPeter Brune 138be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes) 1393542a6bcSPeter Brune { 1403542a6bcSPeter Brune PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(SNESReset_NGS(snes)); 1429566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 1433542a6bcSPeter Brune PetscFunctionReturn(0); 1443542a6bcSPeter Brune } 1453542a6bcSPeter Brune 146be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes) 1473542a6bcSPeter Brune { 148ef107cf5SPeter Brune PetscErrorCode (*f)(SNES,Vec,Vec,void*); 149ef107cf5SPeter Brune 1503542a6bcSPeter Brune PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(SNESGetNGS(snes,&f,NULL)); 152ef107cf5SPeter Brune if (!f) { 1539566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL)); 154ef107cf5SPeter Brune } 1553542a6bcSPeter Brune PetscFunctionReturn(0); 1563542a6bcSPeter Brune } 1573542a6bcSPeter Brune 158*dbbe0bcdSBarry Smith PetscErrorCode SNESSetFromOptions_NGS(SNES snes,PetscOptionItems *PetscOptionsObject) 1593542a6bcSPeter Brune { 160be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 161b6266c6eSPeter Brune PetscInt sweeps,max_its=PETSC_DEFAULT; 162b6266c6eSPeter Brune PetscReal rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT; 163b6266c6eSPeter Brune PetscBool flg,flg1,flg2,flg3; 1643542a6bcSPeter Brune 1653542a6bcSPeter Brune PetscFunctionBegin; 166d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"SNES GS options"); 1676cd56b8bSPeter Brune /* GS Options */ 1689566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg)); 1691baa6e33SBarry Smith if (flg) PetscCall(SNESNGSSetSweeps(snes,sweeps)); 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg)); 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1)); 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2)); 1739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngs_max_it","Maximum number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3)); 174b6266c6eSPeter Brune if (flg || flg1 || flg2 || flg3) { 1759566063dSJacob Faibussowitsch PetscCall(SNESNGSSetTolerances(snes,atol,rtol,stol,max_its)); 176b6266c6eSPeter Brune } 177ef107cf5SPeter Brune flg = PETSC_FALSE; 1789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL)); 179ef107cf5SPeter Brune if (flg) { 1809566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL)); 1819566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n")); 182ef107cf5SPeter Brune } 1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL)); 1849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngs_secant_mat_coloring","Use the graph coloring of the Jacobian for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg)); 185ef107cf5SPeter Brune 186d0609cedSBarry Smith PetscOptionsHeadEnd(); 1873542a6bcSPeter Brune PetscFunctionReturn(0); 1883542a6bcSPeter Brune } 1893542a6bcSPeter Brune 190be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) 1913542a6bcSPeter Brune { 19298272f38SBarry Smith PetscErrorCode (*f)(SNES,Vec,Vec,void*); 19398272f38SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 19498272f38SBarry Smith PetscBool iascii; 19598272f38SBarry Smith 1963542a6bcSPeter Brune PetscFunctionBegin; 1979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 19898272f38SBarry Smith if (iascii) { 1999566063dSJacob Faibussowitsch PetscCall(DMSNESGetNGS(snes->dm,&f,NULL)); 20098272f38SBarry Smith if (f == SNESComputeNGSDefaultSecant) { 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h)); 20298272f38SBarry Smith } 20398272f38SBarry Smith } 2043542a6bcSPeter Brune PetscFunctionReturn(0); 2053542a6bcSPeter Brune } 2063542a6bcSPeter Brune 207be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes) 2083542a6bcSPeter Brune { 2093542a6bcSPeter Brune Vec F; 2103542a6bcSPeter Brune Vec X; 2113542a6bcSPeter Brune Vec B; 2123542a6bcSPeter Brune PetscInt i; 21306e07b1aSPeter Brune PetscReal fnorm; 214365a6726SPeter Brune SNESNormSchedule normschedule; 2153542a6bcSPeter Brune 2163542a6bcSPeter Brune PetscFunctionBegin; 217c579b300SPatrick Farrell 2180b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 219c579b300SPatrick Farrell 2209566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 2213542a6bcSPeter Brune X = snes->vec_sol; 2223542a6bcSPeter Brune F = snes->vec_func; 2233542a6bcSPeter Brune B = snes->vec_rhs; 2243542a6bcSPeter Brune 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 226534ebe21SPeter Brune snes->iter = 0; 227534ebe21SPeter Brune snes->norm = 0.; 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 229526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 230534ebe21SPeter Brune 2319566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 232365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 2333542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 234e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2359566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2361aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 237e4ed7901SPeter Brune 2389566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 239422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 241e4ed7901SPeter Brune snes->iter = 0; 2423542a6bcSPeter Brune snes->norm = fnorm; 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2449566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 2459566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,snes->norm)); 2463542a6bcSPeter Brune 2473542a6bcSPeter Brune /* test convergence */ 248*dbbe0bcdSBarry Smith PetscUseTypeMethod(snes,converged ,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP); 2493542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 250534ebe21SPeter Brune } else { 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2529566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 253534ebe21SPeter Brune } 254534ebe21SPeter Brune 2553542a6bcSPeter Brune /* Call general purpose update function */ 256*dbbe0bcdSBarry Smith PetscTryTypeMethod(snes,update, snes->iter); 257534ebe21SPeter Brune 2583542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 2599566063dSJacob Faibussowitsch PetscCall(SNESComputeNGS(snes, B, X)); 2604b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 261365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 2629566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2639566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 264422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2653542a6bcSPeter Brune /* Monitor convergence */ 2669566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 2673542a6bcSPeter Brune snes->iter = i+1; 2683542a6bcSPeter Brune snes->norm = fnorm; 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2709566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 2719566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 272ac4a6415SMatthew G. Knepley } 2733542a6bcSPeter Brune /* Test for convergence */ 274*dbbe0bcdSBarry Smith if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes,converged ,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP); 2753542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2763542a6bcSPeter Brune /* Call general purpose update function */ 277*dbbe0bcdSBarry Smith PetscTryTypeMethod(snes,update, snes->iter); 2783542a6bcSPeter Brune } 279365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 280534ebe21SPeter Brune if (i == snes->max_its) { 28163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",snes->max_its)); 282534ebe21SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 283534ebe21SPeter Brune } 2841aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 2853542a6bcSPeter Brune PetscFunctionReturn(0); 2863542a6bcSPeter Brune } 2873542a6bcSPeter Brune 2883542a6bcSPeter Brune /*MC 28998272f38SBarry Smith SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation 29098272f38SBarry Smith using coloring. 2913542a6bcSPeter Brune 2923542a6bcSPeter Brune Level: advanced 2933542a6bcSPeter Brune 29498272f38SBarry Smith Options Database: 29578e24b28SBarry Smith + -snes_ngs_sweeps <n> - Number of sweeps of GS to apply 29678e24b28SBarry Smith . -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration 29778e24b28SBarry Smith . -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration 29878e24b28SBarry Smith . -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration 29978e24b28SBarry Smith . -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply 30098272f38SBarry Smith . -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is 30198272f38SBarry Smith used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring 30298272f38SBarry Smith is available or a Jacobian sparse matrix is provided (from which to get the coloring). 30378e24b28SBarry Smith . -snes_ngs_secant_h <h> - Differencing parameter for secant approximation 30478e24b28SBarry Smith . -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available. 305a5b23f4aSJose E. Roman - -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed 30698272f38SBarry Smith 3073542a6bcSPeter Brune Notes: 308d814e60bSBarry Smith the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetNPC(), it will have 3093542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 3103542a6bcSPeter Brune 31178e24b28SBarry Smith By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with SNESSetNormSchedule() or -snes_norm_schedule 3124f02bc6aSBarry Smith References: 313606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 3144f02bc6aSBarry Smith SIAM Review, 57(4), 2015 3154f02bc6aSBarry Smith 316db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`, 317db781477SPatrick Sanan `SNESSetNormSchedule()` 3183542a6bcSPeter Brune M*/ 3193542a6bcSPeter Brune 320be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) 3213542a6bcSPeter Brune { 322be95d8f1SBarry Smith SNES_NGS *gs; 3233542a6bcSPeter Brune 3243542a6bcSPeter Brune PetscFunctionBegin; 325be95d8f1SBarry Smith snes->ops->destroy = SNESDestroy_NGS; 326be95d8f1SBarry Smith snes->ops->setup = SNESSetUp_NGS; 327be95d8f1SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NGS; 328be95d8f1SBarry Smith snes->ops->view = SNESView_NGS; 329be95d8f1SBarry Smith snes->ops->solve = SNESSolve_NGS; 330be95d8f1SBarry Smith snes->ops->reset = SNESReset_NGS; 3313542a6bcSPeter Brune 3323542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 333efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 3343542a6bcSPeter Brune 3354fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 3364fc747eaSLawrence Mitchell 33788976e71SPeter Brune if (!snes->tolerancesset) { 3380e444f03SPeter Brune snes->max_its = 10000; 3390e444f03SPeter Brune snes->max_funcs = 10000; 34088976e71SPeter Brune } 3410e444f03SPeter Brune 3429566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&gs)); 3436cd56b8bSPeter Brune 3446cd56b8bSPeter Brune gs->sweeps = 1; 345b6266c6eSPeter Brune gs->rtol = 1e-5; 34678e24b28SBarry Smith gs->abstol = PETSC_MACHINE_EPSILON; 34778e24b28SBarry Smith gs->stol = 1000*PETSC_MACHINE_EPSILON; 348b6266c6eSPeter Brune gs->max_its = 50; 34978e24b28SBarry Smith gs->h = PETSC_SQRT_MACHINE_EPSILON; 3506cd56b8bSPeter Brune 3513542a6bcSPeter Brune snes->data = (void*) gs; 3523542a6bcSPeter Brune PetscFunctionReturn(0); 3533542a6bcSPeter Brune } 354