1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h" I*/ 26cd56b8bSPeter Brune 3b6266c6eSPeter Brune /*@ 4f6dfbefdSBarry Smith SNESNGSSetTolerances - Sets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG` 5b6266c6eSPeter Brune 6c3339decSBarry Smith Logically Collective 7b6266c6eSPeter Brune 8b6266c6eSPeter Brune Input Parameters: 9f6dfbefdSBarry Smith + 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 23f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetTrustRegionTolerance()` 24b6266c6eSPeter Brune @*/ 25d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit) 26d71ae5a4SJacob Faibussowitsch { 27be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 28b6266c6eSPeter Brune 29b6266c6eSPeter Brune PetscFunctionBegin; 30b6266c6eSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 31b6266c6eSPeter Brune 3213bcc0bdSJacob Faibussowitsch if (abstol != (PetscReal)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 } 3613bcc0bdSJacob Faibussowitsch if (rtol != (PetscReal)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 } 4013bcc0bdSJacob Faibussowitsch if (stol != (PetscReal)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 } 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49b6266c6eSPeter Brune } 50b6266c6eSPeter Brune 51b6266c6eSPeter Brune /*@ 52f6dfbefdSBarry Smith SNESNGSGetTolerances - Gets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG` 53b6266c6eSPeter Brune 54b6266c6eSPeter Brune Not Collective 55b6266c6eSPeter Brune 56b6266c6eSPeter Brune Input Parameters: 57f6dfbefdSBarry Smith + 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 642fe279fdSBarry Smith Level: intermediate 652fe279fdSBarry Smith 66f6dfbefdSBarry Smith Note: 670298fd71SBarry Smith The user can specify NULL for any parameter that is not needed. 68b6266c6eSPeter Brune 69f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetTolerances()` 70b6266c6eSPeter Brune @*/ 71d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit) 72d71ae5a4SJacob Faibussowitsch { 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; 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82b6266c6eSPeter Brune } 836cd56b8bSPeter Brune 846cd56b8bSPeter Brune /*@ 85f6dfbefdSBarry Smith SNESNGSSetSweeps - Sets the number of sweeps of nonlinear GS to use in `SNESNCG` 866cd56b8bSPeter Brune 876cd56b8bSPeter Brune Input Parameters: 88f6dfbefdSBarry Smith + snes - the `SNES` context 89f6dfbefdSBarry Smith - sweeps - the number of sweeps of nonlinear GS to perform. 90f6dfbefdSBarry Smith 91f6dfbefdSBarry Smith Options Database Key: 92f6dfbefdSBarry Smith . -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply 936cd56b8bSPeter Brune 946cd56b8bSPeter Brune Level: intermediate 956cd56b8bSPeter Brune 96f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()` 976cd56b8bSPeter Brune @*/ 986cd56b8bSPeter Brune 99d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) 100d71ae5a4SJacob Faibussowitsch { 101be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 1026cd56b8bSPeter Brune 1036cd56b8bSPeter Brune PetscFunctionBegin; 1046cd56b8bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1056cd56b8bSPeter Brune gs->sweeps = sweeps; 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1076cd56b8bSPeter Brune } 1086cd56b8bSPeter Brune 1096cd56b8bSPeter Brune /*@ 110f6dfbefdSBarry Smith SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG` 1116cd56b8bSPeter Brune 1122fe279fdSBarry Smith Input Parameter: 113f6dfbefdSBarry Smith . snes - the `SNES` context 1146cd56b8bSPeter Brune 1152fe279fdSBarry Smith Output Parameter: 116f6dfbefdSBarry Smith . sweeps - the number of sweeps of nonlinear GS to perform. 1176cd56b8bSPeter Brune 1186cd56b8bSPeter Brune Level: intermediate 1196cd56b8bSPeter Brune 120f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()` 1216cd56b8bSPeter Brune @*/ 122d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps) 123d71ae5a4SJacob Faibussowitsch { 124be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 1256cd56b8bSPeter Brune 1266cd56b8bSPeter Brune PetscFunctionBegin; 1276cd56b8bSPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1286cd56b8bSPeter Brune *sweeps = gs->sweeps; 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1306cd56b8bSPeter Brune } 1316cd56b8bSPeter Brune 132d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NGS(SNES snes) 133d71ae5a4SJacob Faibussowitsch { 1348a86d3c5SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 1358a86d3c5SBarry Smith 1363542a6bcSPeter Brune PetscFunctionBegin; 1379566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&gs->coloring)); 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1393542a6bcSPeter Brune } 1403542a6bcSPeter Brune 141d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NGS(SNES snes) 142d71ae5a4SJacob Faibussowitsch { 1433542a6bcSPeter Brune PetscFunctionBegin; 1449566063dSJacob Faibussowitsch PetscCall(SNESReset_NGS(snes)); 1459566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1473542a6bcSPeter Brune } 1483542a6bcSPeter Brune 149d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NGS(SNES snes) 150d71ae5a4SJacob Faibussowitsch { 151ef107cf5SPeter Brune PetscErrorCode (*f)(SNES, Vec, Vec, void *); 152ef107cf5SPeter Brune 1533542a6bcSPeter Brune PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(SNESGetNGS(snes, &f, NULL)); 15548a46eb9SPierre Jolivet if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1573542a6bcSPeter Brune } 1583542a6bcSPeter Brune 159d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject) 160d71ae5a4SJacob Faibussowitsch { 161be95d8f1SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 162b6266c6eSPeter Brune PetscInt sweeps, max_its = PETSC_DEFAULT; 163b6266c6eSPeter Brune PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT; 164b6266c6eSPeter Brune PetscBool flg, flg1, flg2, flg3; 1653542a6bcSPeter Brune 1663542a6bcSPeter Brune PetscFunctionBegin; 167d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options"); 1686cd56b8bSPeter Brune /* GS Options */ 1699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg)); 1701baa6e33SBarry Smith if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps)); 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg)); 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1)); 1739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2)); 1749566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3)); 17548a46eb9SPierre Jolivet if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its)); 176ef107cf5SPeter Brune flg = PETSC_FALSE; 1779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL)); 178ef107cf5SPeter Brune if (flg) { 1799566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL)); 1809566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n")); 181ef107cf5SPeter Brune } 1829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL)); 1839566063dSJacob 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)); 184ef107cf5SPeter Brune 185d0609cedSBarry Smith PetscOptionsHeadEnd(); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1873542a6bcSPeter Brune } 1883542a6bcSPeter Brune 189d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) 190d71ae5a4SJacob Faibussowitsch { 19198272f38SBarry Smith PetscErrorCode (*f)(SNES, Vec, Vec, void *); 19298272f38SBarry Smith SNES_NGS *gs = (SNES_NGS *)snes->data; 19398272f38SBarry Smith PetscBool iascii; 19498272f38SBarry Smith 1953542a6bcSPeter Brune PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 19798272f38SBarry Smith if (iascii) { 1989566063dSJacob Faibussowitsch PetscCall(DMSNESGetNGS(snes->dm, &f, NULL)); 19948a46eb9SPierre Jolivet if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, " Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h)); 20098272f38SBarry Smith } 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2023542a6bcSPeter Brune } 2033542a6bcSPeter Brune 204d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NGS(SNES snes) 205d71ae5a4SJacob Faibussowitsch { 2063542a6bcSPeter Brune Vec F; 2073542a6bcSPeter Brune Vec X; 2083542a6bcSPeter Brune Vec B; 2093542a6bcSPeter Brune PetscInt i; 21006e07b1aSPeter Brune PetscReal fnorm; 211365a6726SPeter Brune SNESNormSchedule normschedule; 2123542a6bcSPeter Brune 2133542a6bcSPeter Brune PetscFunctionBegin; 214c579b300SPatrick Farrell 2150b121fc5SBarry 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); 216c579b300SPatrick Farrell 2179566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 2183542a6bcSPeter Brune X = snes->vec_sol; 2193542a6bcSPeter Brune F = snes->vec_func; 2203542a6bcSPeter Brune B = snes->vec_rhs; 2213542a6bcSPeter Brune 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 223534ebe21SPeter Brune snes->iter = 0; 224534ebe21SPeter Brune snes->norm = 0.; 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 226526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 227534ebe21SPeter Brune 2289566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 229365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 2303542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 231e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2329566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2331aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 234e4ed7901SPeter Brune 2359566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 236422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 238e4ed7901SPeter Brune snes->iter = 0; 2393542a6bcSPeter Brune snes->norm = fnorm; 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2419566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 2423542a6bcSPeter Brune 2433542a6bcSPeter Brune /* test convergence */ 244*d76a863bSStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 2452d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, snes->norm)); 2463ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 247534ebe21SPeter Brune } else { 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2499566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 250534ebe21SPeter Brune } 251534ebe21SPeter Brune 2523542a6bcSPeter Brune /* Call general purpose update function */ 253dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 254534ebe21SPeter Brune 2553542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 2569566063dSJacob Faibussowitsch PetscCall(SNESComputeNGS(snes, B, X)); 2574b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 258365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 2599566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 2609566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 261422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2622d157150SStefano Zampini } 2633542a6bcSPeter Brune /* Monitor convergence */ 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 2653542a6bcSPeter Brune snes->iter = i + 1; 2663542a6bcSPeter Brune snes->norm = fnorm; 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2682d157150SStefano Zampini PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 2693542a6bcSPeter Brune /* Test for convergence */ 2702d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 2712d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 2723ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 2733542a6bcSPeter Brune /* Call general purpose update function */ 274dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 2753542a6bcSPeter Brune } 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2773542a6bcSPeter Brune } 2783542a6bcSPeter Brune 2793542a6bcSPeter Brune /*MC 280f6dfbefdSBarry Smith SNESNGS - Either calls the user-provided solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation 28198272f38SBarry Smith using coloring. 2823542a6bcSPeter Brune 2833542a6bcSPeter Brune Level: advanced 2843542a6bcSPeter Brune 285f6dfbefdSBarry Smith Options Database Keys: 286f6dfbefdSBarry Smith + -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply 287f6dfbefdSBarry Smith . -snes_ngs_atol <atol> - Absolute residual tolerance for nonlinear GS iteration 288f6dfbefdSBarry Smith . -snes_ngs_rtol <rtol> - Relative residual tolerance for nonlinear GS iteration 289f6dfbefdSBarry Smith . -snes_ngs_stol <stol> - Absolute update tolerance for nonlinear GS iteration 290f6dfbefdSBarry Smith . -snes_ngs_max_it <maxit> - Maximum number of sweeps of nonlinea GS to apply 29198272f38SBarry Smith . -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is 292f6dfbefdSBarry Smith used by default if no user provided Gauss-Seidel routine is available. Requires either that a `DM` that can compute a coloring 29398272f38SBarry Smith is available or a Jacobian sparse matrix is provided (from which to get the coloring). 29478e24b28SBarry Smith . -snes_ngs_secant_h <h> - Differencing parameter for secant approximation 29578e24b28SBarry Smith . -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available. 296a5b23f4aSJose E. Roman - -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed 29798272f38SBarry Smith 2983542a6bcSPeter Brune Notes: 299f6dfbefdSBarry Smith the Gauss-Seidel smoother is inherited through composition. If a solver has been created with `SNESGetNPC()`, it will have 3003542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 3013542a6bcSPeter Brune 302f6dfbefdSBarry Smith By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()` 303f6dfbefdSBarry Smith or -snes_norm_schedule none 304f6dfbefdSBarry Smith 3054f02bc6aSBarry Smith References: 306606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 3074f02bc6aSBarry Smith SIAM Review, 57(4), 2015 3084f02bc6aSBarry Smith 309f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`, 310db781477SPatrick Sanan `SNESSetNormSchedule()` 3113542a6bcSPeter Brune M*/ 3123542a6bcSPeter Brune 313d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) 314d71ae5a4SJacob Faibussowitsch { 315be95d8f1SBarry Smith SNES_NGS *gs; 3163542a6bcSPeter Brune 3173542a6bcSPeter Brune PetscFunctionBegin; 318be95d8f1SBarry Smith snes->ops->destroy = SNESDestroy_NGS; 319be95d8f1SBarry Smith snes->ops->setup = SNESSetUp_NGS; 320be95d8f1SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_NGS; 321be95d8f1SBarry Smith snes->ops->view = SNESView_NGS; 322be95d8f1SBarry Smith snes->ops->solve = SNESSolve_NGS; 323be95d8f1SBarry Smith snes->ops->reset = SNESReset_NGS; 3243542a6bcSPeter Brune 3253542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 326efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 3273542a6bcSPeter Brune 3284fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 3294fc747eaSLawrence Mitchell 33088976e71SPeter Brune if (!snes->tolerancesset) { 3310e444f03SPeter Brune snes->max_its = 10000; 3320e444f03SPeter Brune snes->max_funcs = 10000; 33388976e71SPeter Brune } 3340e444f03SPeter Brune 3354dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&gs)); 3366cd56b8bSPeter Brune 3376cd56b8bSPeter Brune gs->sweeps = 1; 338b6266c6eSPeter Brune gs->rtol = 1e-5; 33978e24b28SBarry Smith gs->abstol = PETSC_MACHINE_EPSILON; 34078e24b28SBarry Smith gs->stol = 1000 * PETSC_MACHINE_EPSILON; 341b6266c6eSPeter Brune gs->max_its = 50; 34278e24b28SBarry Smith gs->h = PETSC_SQRT_MACHINE_EPSILON; 3436cd56b8bSPeter Brune 3443542a6bcSPeter Brune snes->data = (void *)gs; 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3463542a6bcSPeter Brune } 347