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 23b6266c6eSPeter Brune .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) { 33*08401ef6SPierre 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) { 372c71b3e2SJacob Faibussowitsch PetscCheckFalse(rtol < 0.0 || 1.0 <= rtol,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) { 41*08401ef6SPierre 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) { 45*08401ef6SPierre Jolivet PetscCheck(maxit >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D 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 69b6266c6eSPeter Brune .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 93d814e60bSBarry Smith .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 117d814e60bSBarry Smith .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 1584416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes) 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; 1669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(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)); 169b6266c6eSPeter Brune if (flg) { 1709566063dSJacob Faibussowitsch PetscCall(SNESNGSSetSweeps(snes,sweeps)); 171b6266c6eSPeter Brune } 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg)); 1739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1)); 1749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2)); 1759566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngs_max_it","Maximum number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3)); 176b6266c6eSPeter Brune if (flg || flg1 || flg2 || flg3) { 1779566063dSJacob Faibussowitsch PetscCall(SNESNGSSetTolerances(snes,atol,rtol,stol,max_its)); 178b6266c6eSPeter Brune } 179ef107cf5SPeter Brune flg = PETSC_FALSE; 1809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL)); 181ef107cf5SPeter Brune if (flg) { 1829566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL)); 1839566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n")); 184ef107cf5SPeter Brune } 1859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL)); 1869566063dSJacob 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)); 187ef107cf5SPeter Brune 1889566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 1893542a6bcSPeter Brune PetscFunctionReturn(0); 1903542a6bcSPeter Brune } 1913542a6bcSPeter Brune 192be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) 1933542a6bcSPeter Brune { 19498272f38SBarry Smith PetscErrorCode (*f)(SNES,Vec,Vec,void*); 19598272f38SBarry Smith SNES_NGS *gs = (SNES_NGS*)snes->data; 19698272f38SBarry Smith PetscBool iascii; 19798272f38SBarry Smith 1983542a6bcSPeter Brune PetscFunctionBegin; 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 20098272f38SBarry Smith if (iascii) { 2019566063dSJacob Faibussowitsch PetscCall(DMSNESGetNGS(snes->dm,&f,NULL)); 20298272f38SBarry Smith if (f == SNESComputeNGSDefaultSecant) { 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h)); 20498272f38SBarry Smith } 20598272f38SBarry Smith } 2063542a6bcSPeter Brune PetscFunctionReturn(0); 2073542a6bcSPeter Brune } 2083542a6bcSPeter Brune 209be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes) 2103542a6bcSPeter Brune { 2113542a6bcSPeter Brune Vec F; 2123542a6bcSPeter Brune Vec X; 2133542a6bcSPeter Brune Vec B; 2143542a6bcSPeter Brune PetscInt i; 21506e07b1aSPeter Brune PetscReal fnorm; 216365a6726SPeter Brune SNESNormSchedule normschedule; 2173542a6bcSPeter Brune 2183542a6bcSPeter Brune PetscFunctionBegin; 219c579b300SPatrick Farrell 2202c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 221c579b300SPatrick Farrell 2229566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 2233542a6bcSPeter Brune X = snes->vec_sol; 2243542a6bcSPeter Brune F = snes->vec_func; 2253542a6bcSPeter Brune B = snes->vec_rhs; 2263542a6bcSPeter Brune 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 228534ebe21SPeter Brune snes->iter = 0; 229534ebe21SPeter Brune snes->norm = 0.; 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 231526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 232534ebe21SPeter Brune 2339566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 234365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 2353542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 236e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2379566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2381aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 239e4ed7901SPeter Brune 2409566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 241422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 243e4ed7901SPeter Brune snes->iter = 0; 2443542a6bcSPeter Brune snes->norm = fnorm; 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2469566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 2479566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,0,snes->norm)); 2483542a6bcSPeter Brune 2493542a6bcSPeter Brune /* test convergence */ 2509566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 2513542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 252534ebe21SPeter Brune } else { 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2549566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 255534ebe21SPeter Brune } 256534ebe21SPeter Brune 2573542a6bcSPeter Brune /* Call general purpose update function */ 2583542a6bcSPeter Brune if (snes->ops->update) { 2599566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 2603542a6bcSPeter Brune } 261534ebe21SPeter Brune 2623542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 2639566063dSJacob Faibussowitsch PetscCall(SNESComputeNGS(snes, B, X)); 2644b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 265365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 2669566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 2679566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 268422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 2693542a6bcSPeter Brune /* Monitor convergence */ 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 2713542a6bcSPeter Brune snes->iter = i+1; 2723542a6bcSPeter Brune snes->norm = fnorm; 2739566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2749566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 2759566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 276ac4a6415SMatthew G. Knepley } 2773542a6bcSPeter Brune /* Test for convergence */ 2789566063dSJacob Faibussowitsch if (normschedule == SNES_NORM_ALWAYS) PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 2793542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2803542a6bcSPeter Brune /* Call general purpose update function */ 2813542a6bcSPeter Brune if (snes->ops->update) { 2829566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 2833542a6bcSPeter Brune } 2843542a6bcSPeter Brune } 285365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 286534ebe21SPeter Brune if (i == snes->max_its) { 2879566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its)); 288534ebe21SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 289534ebe21SPeter Brune } 2901aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 2913542a6bcSPeter Brune PetscFunctionReturn(0); 2923542a6bcSPeter Brune } 2933542a6bcSPeter Brune 2943542a6bcSPeter Brune /*MC 29598272f38SBarry Smith SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation 29698272f38SBarry Smith using coloring. 2973542a6bcSPeter Brune 2983542a6bcSPeter Brune Level: advanced 2993542a6bcSPeter Brune 30098272f38SBarry Smith Options Database: 30178e24b28SBarry Smith + -snes_ngs_sweeps <n> - Number of sweeps of GS to apply 30278e24b28SBarry Smith . -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration 30378e24b28SBarry Smith . -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration 30478e24b28SBarry Smith . -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration 30578e24b28SBarry Smith . -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply 30698272f38SBarry Smith . -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is 30798272f38SBarry Smith used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring 30898272f38SBarry Smith is available or a Jacobian sparse matrix is provided (from which to get the coloring). 30978e24b28SBarry Smith . -snes_ngs_secant_h <h> - Differencing parameter for secant approximation 31078e24b28SBarry Smith . -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available. 311a5b23f4aSJose E. Roman - -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed 31298272f38SBarry Smith 3133542a6bcSPeter Brune Notes: 314d814e60bSBarry Smith the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetNPC(), it will have 3153542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 3163542a6bcSPeter Brune 31778e24b28SBarry 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 3184f02bc6aSBarry Smith References: 319606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 3204f02bc6aSBarry Smith SIAM Review, 57(4), 2015 3214f02bc6aSBarry Smith 32278e24b28SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances(), 32378e24b28SBarry Smith SNESSetNormSchedule() 3243542a6bcSPeter Brune M*/ 3253542a6bcSPeter Brune 326be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) 3273542a6bcSPeter Brune { 328be95d8f1SBarry Smith SNES_NGS *gs; 3293542a6bcSPeter Brune 3303542a6bcSPeter Brune PetscFunctionBegin; 331be95d8f1SBarry Smith snes->ops->destroy = SNESDestroy_NGS; 332be95d8f1SBarry Smith snes->ops->setup = SNESSetUp_NGS; 333be95d8f1SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NGS; 334be95d8f1SBarry Smith snes->ops->view = SNESView_NGS; 335be95d8f1SBarry Smith snes->ops->solve = SNESSolve_NGS; 336be95d8f1SBarry Smith snes->ops->reset = SNESReset_NGS; 3373542a6bcSPeter Brune 3383542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 339efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 3403542a6bcSPeter Brune 3414fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 3424fc747eaSLawrence Mitchell 34388976e71SPeter Brune if (!snes->tolerancesset) { 3440e444f03SPeter Brune snes->max_its = 10000; 3450e444f03SPeter Brune snes->max_funcs = 10000; 34688976e71SPeter Brune } 3470e444f03SPeter Brune 3489566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&gs)); 3496cd56b8bSPeter Brune 3506cd56b8bSPeter Brune gs->sweeps = 1; 351b6266c6eSPeter Brune gs->rtol = 1e-5; 35278e24b28SBarry Smith gs->abstol = PETSC_MACHINE_EPSILON; 35378e24b28SBarry Smith gs->stol = 1000*PETSC_MACHINE_EPSILON; 354b6266c6eSPeter Brune gs->max_its = 50; 35578e24b28SBarry Smith gs->h = PETSC_SQRT_MACHINE_EPSILON; 3566cd56b8bSPeter Brune 3573542a6bcSPeter Brune snes->data = (void*) gs; 3583542a6bcSPeter Brune PetscFunctionReturn(0); 3593542a6bcSPeter Brune } 360