1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 23542a6bcSPeter Brune 33542a6bcSPeter Brune typedef struct { 4b6266c6eSPeter Brune PetscInt sweeps; /* number of sweeps through the local subdomain before neighbor communication */ 5b6266c6eSPeter Brune PetscInt max_its; /* maximum iterations of the inner pointblock solver */ 6b6266c6eSPeter Brune PetscReal rtol; /* relative tolerance of the inner pointblock solver */ 7b6266c6eSPeter Brune PetscReal abstol; /* absolute tolerance of the inner pointblock solver */ 8b6266c6eSPeter Brune PetscReal stol; /* step tolerance of the inner pointblock solver */ 93542a6bcSPeter Brune } SNES_GS; 103542a6bcSPeter Brune 116cd56b8bSPeter Brune 12b6266c6eSPeter Brune #undef __FUNCT__ 13b6266c6eSPeter Brune #define __FUNCT__ "SNESGSSetTolerances" 14b6266c6eSPeter Brune /*@ 15b6266c6eSPeter Brune SNESGSSetTolerances - Sets various parameters used in convergence tests. 16b6266c6eSPeter Brune 17b6266c6eSPeter Brune Logically Collective on SNES 18b6266c6eSPeter Brune 19b6266c6eSPeter Brune Input Parameters: 20b6266c6eSPeter Brune + snes - the SNES context 21b6266c6eSPeter Brune . abstol - absolute convergence tolerance 22b6266c6eSPeter Brune . rtol - relative convergence tolerance 23b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x || 24b6266c6eSPeter Brune - maxit - maximum number of iterations 25b6266c6eSPeter Brune 26b6266c6eSPeter Brune 27b6266c6eSPeter Brune Options Database Keys: 28b6266c6eSPeter Brune + -snes_gs_atol <abstol> - Sets abstol 29b6266c6eSPeter Brune . -snes_gs_rtol <rtol> - Sets rtol 30b6266c6eSPeter Brune . -snes_gs_stol <stol> - Sets stol 31b6266c6eSPeter Brune - -snes_max_it <maxit> - Sets maxit 32b6266c6eSPeter Brune 33b6266c6eSPeter Brune Level: intermediate 34b6266c6eSPeter Brune 35b6266c6eSPeter Brune .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances 36b6266c6eSPeter Brune 37b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance() 38b6266c6eSPeter Brune @*/ 39b6266c6eSPeter Brune PetscErrorCode SNESGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit) 40b6266c6eSPeter Brune { 41b6266c6eSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 42b6266c6eSPeter Brune 43b6266c6eSPeter Brune PetscFunctionBegin; 44b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 45b6266c6eSPeter Brune 46b6266c6eSPeter Brune if (abstol != PETSC_DEFAULT) { 47*ce94432eSBarry Smith if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol); 48b6266c6eSPeter Brune gs->abstol = abstol; 49b6266c6eSPeter Brune } 50b6266c6eSPeter Brune if (rtol != PETSC_DEFAULT) { 51*ce94432eSBarry 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",rtol); 52b6266c6eSPeter Brune gs->rtol = rtol; 53b6266c6eSPeter Brune } 54b6266c6eSPeter Brune if (stol != PETSC_DEFAULT) { 55*ce94432eSBarry Smith if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol); 56b6266c6eSPeter Brune gs->stol = stol; 57b6266c6eSPeter Brune } 58b6266c6eSPeter Brune if (maxit != PETSC_DEFAULT) { 59*ce94432eSBarry Smith if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit); 60b6266c6eSPeter Brune gs->max_its = maxit; 61b6266c6eSPeter Brune } 62b6266c6eSPeter Brune PetscFunctionReturn(0); 63b6266c6eSPeter Brune } 64b6266c6eSPeter Brune 65b6266c6eSPeter Brune #undef __FUNCT__ 66b6266c6eSPeter Brune #define __FUNCT__ "SNESGSGetTolerances" 67b6266c6eSPeter Brune /*@ 68b6266c6eSPeter Brune SNESGSGetTolerances - Gets various parameters used in convergence tests. 69b6266c6eSPeter Brune 70b6266c6eSPeter Brune Not Collective 71b6266c6eSPeter Brune 72b6266c6eSPeter Brune Input Parameters: 73b6266c6eSPeter Brune + snes - the SNES context 74b6266c6eSPeter Brune . atol - absolute convergence tolerance 75b6266c6eSPeter Brune . rtol - relative convergence tolerance 76b6266c6eSPeter Brune . stol - convergence tolerance in terms of the norm 77b6266c6eSPeter Brune of the change in the solution between steps 78b6266c6eSPeter Brune - maxit - maximum number of iterations 79b6266c6eSPeter Brune 80b6266c6eSPeter Brune Notes: 810298fd71SBarry Smith The user can specify NULL for any parameter that is not needed. 82b6266c6eSPeter Brune 83b6266c6eSPeter Brune Level: intermediate 84b6266c6eSPeter Brune 85b6266c6eSPeter Brune .keywords: SNES, nonlinear, get, convergence, tolerances 86b6266c6eSPeter Brune 87b6266c6eSPeter Brune .seealso: SNESSetTolerances() 88b6266c6eSPeter Brune @*/ 89b6266c6eSPeter Brune PetscErrorCode SNESGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit) 90b6266c6eSPeter Brune { 91b6266c6eSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 92b6266c6eSPeter Brune 93b6266c6eSPeter Brune PetscFunctionBegin; 94b6266c6eSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 95b6266c6eSPeter Brune if (atol) *atol = gs->abstol; 96b6266c6eSPeter Brune if (rtol) *rtol = gs->rtol; 97b6266c6eSPeter Brune if (stol) *stol = gs->stol; 98b6266c6eSPeter Brune if (maxit) *maxit = gs->max_its; 99b6266c6eSPeter Brune PetscFunctionReturn(0); 100b6266c6eSPeter Brune } 1016cd56b8bSPeter Brune 1026cd56b8bSPeter Brune #undef __FUNCT__ 1036cd56b8bSPeter Brune #define __FUNCT__ "SNESGSSetSweeps" 1046cd56b8bSPeter Brune /*@ 1056cd56b8bSPeter Brune SNESGSSetSweeps - Sets the number of sweeps of GS to use. 1066cd56b8bSPeter Brune 1076cd56b8bSPeter Brune Input Parameters: 1086cd56b8bSPeter Brune + snes - the SNES context 1096cd56b8bSPeter Brune - sweeps - the number of sweeps of GS to perform. 1106cd56b8bSPeter Brune 1116cd56b8bSPeter Brune Level: intermediate 1126cd56b8bSPeter Brune 1136cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 1146cd56b8bSPeter Brune 1156cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps() 1166cd56b8bSPeter Brune @*/ 1176cd56b8bSPeter Brune 1186cd56b8bSPeter Brune PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps) 1196cd56b8bSPeter Brune { 1206cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 1216cd56b8bSPeter Brune 1226cd56b8bSPeter Brune PetscFunctionBegin; 1236cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1246cd56b8bSPeter Brune gs->sweeps = sweeps; 1256cd56b8bSPeter Brune PetscFunctionReturn(0); 1266cd56b8bSPeter Brune } 1276cd56b8bSPeter Brune 1286cd56b8bSPeter Brune #undef __FUNCT__ 1296cd56b8bSPeter Brune #define __FUNCT__ "SNESGSGetSweeps" 1306cd56b8bSPeter Brune /*@ 1316cd56b8bSPeter Brune SNESGSGetSweeps - Gets the number of sweeps GS will use. 1326cd56b8bSPeter Brune 1336cd56b8bSPeter Brune Input Parameters: 1346cd56b8bSPeter Brune . snes - the SNES context 1356cd56b8bSPeter Brune 1366cd56b8bSPeter Brune Output Parameters: 1376cd56b8bSPeter Brune . sweeps - the number of sweeps of GS to perform. 1386cd56b8bSPeter Brune 1396cd56b8bSPeter Brune Level: intermediate 1406cd56b8bSPeter Brune 1416cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 1426cd56b8bSPeter Brune 1436cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps() 1446cd56b8bSPeter Brune @*/ 1456cd56b8bSPeter Brune PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps) 1466cd56b8bSPeter Brune { 1476cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 1486cd56b8bSPeter Brune 1496cd56b8bSPeter Brune PetscFunctionBegin; 1506cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1516cd56b8bSPeter Brune *sweeps = gs->sweeps; 1526cd56b8bSPeter Brune PetscFunctionReturn(0); 1536cd56b8bSPeter Brune } 1546cd56b8bSPeter Brune 1556cd56b8bSPeter Brune 1563542a6bcSPeter Brune #undef __FUNCT__ 157c79d6ae8SPeter Brune #define __FUNCT__ "SNESDefaultApplyGS" 1580adebc6cSBarry Smith PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx) 1590adebc6cSBarry Smith { 160c79d6ae8SPeter Brune PetscFunctionBegin; 161c79d6ae8SPeter Brune /* see if there's a coloring on the DM */ 162c79d6ae8SPeter Brune PetscFunctionReturn(0); 163c79d6ae8SPeter Brune } 164c79d6ae8SPeter Brune 165c79d6ae8SPeter Brune #undef __FUNCT__ 1663542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS" 1673542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes) 1683542a6bcSPeter Brune { 1693542a6bcSPeter Brune PetscFunctionBegin; 1703542a6bcSPeter Brune PetscFunctionReturn(0); 1713542a6bcSPeter Brune } 1723542a6bcSPeter Brune 1733542a6bcSPeter Brune #undef __FUNCT__ 1743542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS" 1753542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes) 1763542a6bcSPeter Brune { 1773542a6bcSPeter Brune PetscErrorCode ierr; 1783542a6bcSPeter Brune 1793542a6bcSPeter Brune PetscFunctionBegin; 1803542a6bcSPeter Brune ierr = SNESReset_GS(snes);CHKERRQ(ierr); 18122d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 1823542a6bcSPeter Brune PetscFunctionReturn(0); 1833542a6bcSPeter Brune } 1843542a6bcSPeter Brune 1853542a6bcSPeter Brune #undef __FUNCT__ 1863542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS" 1873542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes) 1883542a6bcSPeter Brune { 1893542a6bcSPeter Brune PetscFunctionBegin; 1903542a6bcSPeter Brune PetscFunctionReturn(0); 1913542a6bcSPeter Brune } 1923542a6bcSPeter Brune 1933542a6bcSPeter Brune #undef __FUNCT__ 1943542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS" 1953542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes) 1963542a6bcSPeter Brune { 1976cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 1983542a6bcSPeter Brune PetscErrorCode ierr; 199b6266c6eSPeter Brune PetscInt sweeps,max_its=PETSC_DEFAULT; 200b6266c6eSPeter Brune PetscReal rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT; 201b6266c6eSPeter Brune PetscBool flg,flg1,flg2,flg3; 2023542a6bcSPeter Brune 2033542a6bcSPeter Brune PetscFunctionBegin; 2043542a6bcSPeter Brune ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr); 2056cd56b8bSPeter Brune /* GS Options */ 206b6266c6eSPeter Brune ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr); 207b6266c6eSPeter Brune if (flg) { 208b6266c6eSPeter Brune ierr = SNESGSSetSweeps(snes,sweeps);CHKERRQ(ierr); 209b6266c6eSPeter Brune } 210b6266c6eSPeter Brune ierr = PetscOptionsReal("-snes_gs_atol","Number of sweeps of GS to apply","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr); 211b6266c6eSPeter Brune ierr = PetscOptionsReal("-snes_gs_rtol","Number of sweeps of GS to apply","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr); 212b6266c6eSPeter Brune ierr = PetscOptionsReal("-snes_gs_stol","Number of sweeps of GS to apply","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr); 213b6266c6eSPeter Brune ierr = PetscOptionsInt("-snes_gs_max_it","Number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr); 214b6266c6eSPeter Brune if (flg || flg1 || flg2 || flg3) { 215b6266c6eSPeter Brune ierr = SNESGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr); 216b6266c6eSPeter Brune } 2174b32a720SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2183542a6bcSPeter Brune PetscFunctionReturn(0); 2193542a6bcSPeter Brune } 2203542a6bcSPeter Brune 2213542a6bcSPeter Brune #undef __FUNCT__ 2223542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS" 2233542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 2243542a6bcSPeter Brune { 2253542a6bcSPeter Brune PetscFunctionBegin; 2263542a6bcSPeter Brune PetscFunctionReturn(0); 2273542a6bcSPeter Brune } 2283542a6bcSPeter Brune 2293542a6bcSPeter Brune #undef __FUNCT__ 2303542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS" 2313542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes) 2323542a6bcSPeter Brune { 2333542a6bcSPeter Brune Vec F; 2343542a6bcSPeter Brune Vec X; 2353542a6bcSPeter Brune Vec B; 2363542a6bcSPeter Brune PetscInt i; 23706e07b1aSPeter Brune PetscReal fnorm; 2383542a6bcSPeter Brune PetscErrorCode ierr; 239534ebe21SPeter Brune SNESNormType normtype; 2403542a6bcSPeter Brune 2413542a6bcSPeter Brune PetscFunctionBegin; 2423542a6bcSPeter Brune X = snes->vec_sol; 2433542a6bcSPeter Brune F = snes->vec_func; 2443542a6bcSPeter Brune B = snes->vec_rhs; 2453542a6bcSPeter Brune 246534ebe21SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 247534ebe21SPeter Brune snes->iter = 0; 248534ebe21SPeter Brune snes->norm = 0.; 249534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 250526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 251534ebe21SPeter Brune 252534ebe21SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 253fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 2543542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 255e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 2563542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2573542a6bcSPeter Brune if (snes->domainerror) { 2583542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2593542a6bcSPeter Brune PetscFunctionReturn(0); 2603542a6bcSPeter Brune } 2611aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 262e4ed7901SPeter Brune 2633542a6bcSPeter Brune /* convergence test */ 264e4ed7901SPeter Brune if (!snes->norm_init_set) { 2653542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 2663542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 267e4ed7901SPeter Brune } else { 268e4ed7901SPeter Brune fnorm = snes->norm_init; 269e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 270e4ed7901SPeter Brune } 2713542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 272e4ed7901SPeter Brune snes->iter = 0; 2733542a6bcSPeter Brune snes->norm = fnorm; 2743542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 275a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 276e4ed7901SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 2773542a6bcSPeter Brune 2783542a6bcSPeter Brune /* set parameter for default relative tolerance convergence test */ 2793542a6bcSPeter Brune snes->ttol = fnorm*snes->rtol; 280534ebe21SPeter 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 { 285534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 286a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 287534ebe21SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 288534ebe21SPeter Brune } 289534ebe21SPeter Brune 2903542a6bcSPeter Brune 2913542a6bcSPeter Brune /* Call general purpose update function */ 2923542a6bcSPeter Brune if (snes->ops->update) { 2933542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 2943542a6bcSPeter Brune } 295534ebe21SPeter Brune 2963542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 2973542a6bcSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 2984b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 299fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 3003542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 3013542a6bcSPeter Brune if (snes->domainerror) { 3023542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 3033542a6bcSPeter Brune PetscFunctionReturn(0); 3043542a6bcSPeter Brune } 3053542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 3063542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 307534ebe21SPeter Brune } 3083542a6bcSPeter Brune /* Monitor convergence */ 3093542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 3103542a6bcSPeter Brune snes->iter = i+1; 3113542a6bcSPeter Brune snes->norm = fnorm; 3123542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 313a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 3143542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 3153542a6bcSPeter Brune /* Test for convergence */ 316fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 3173542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 3183542a6bcSPeter Brune /* Call general purpose update function */ 3193542a6bcSPeter Brune if (snes->ops->update) { 3203542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 3213542a6bcSPeter Brune } 3223542a6bcSPeter Brune } 323fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 324534ebe21SPeter Brune if (i == snes->max_its) { 325534ebe21SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 326534ebe21SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 327534ebe21SPeter Brune } 3281aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 3293542a6bcSPeter Brune PetscFunctionReturn(0); 3303542a6bcSPeter Brune } 3313542a6bcSPeter Brune 3323542a6bcSPeter Brune /*MC 3333542a6bcSPeter Brune SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 3343542a6bcSPeter Brune 3353542a6bcSPeter Brune Level: advanced 3363542a6bcSPeter Brune 3373542a6bcSPeter Brune Notes: 3383542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 3393542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 3403542a6bcSPeter Brune 3413542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 3423542a6bcSPeter Brune M*/ 3433542a6bcSPeter Brune 3443542a6bcSPeter Brune EXTERN_C_BEGIN 3453542a6bcSPeter Brune #undef __FUNCT__ 3463542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS" 3473542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes) 3483542a6bcSPeter Brune { 3493542a6bcSPeter Brune SNES_GS *gs; 3503542a6bcSPeter Brune PetscErrorCode ierr; 3513542a6bcSPeter Brune 3523542a6bcSPeter Brune PetscFunctionBegin; 3533542a6bcSPeter Brune snes->ops->destroy = SNESDestroy_GS; 3543542a6bcSPeter Brune snes->ops->setup = SNESSetUp_GS; 3553542a6bcSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_GS; 3563542a6bcSPeter Brune snes->ops->view = SNESView_GS; 3573542a6bcSPeter Brune snes->ops->solve = SNESSolve_GS; 3583542a6bcSPeter Brune snes->ops->reset = SNESReset_GS; 3593542a6bcSPeter Brune 3603542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 3613542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 3623542a6bcSPeter Brune 36388976e71SPeter Brune if (!snes->tolerancesset) { 3640e444f03SPeter Brune snes->max_its = 10000; 3650e444f03SPeter Brune snes->max_funcs = 10000; 36688976e71SPeter Brune } 3670e444f03SPeter Brune 3683542a6bcSPeter Brune ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 3696cd56b8bSPeter Brune 3706cd56b8bSPeter Brune gs->sweeps = 1; 371b6266c6eSPeter Brune gs->rtol = 1e-5; 372b6266c6eSPeter Brune gs->abstol = 1e-15; 373b6266c6eSPeter Brune gs->stol = 1e-12; 374b6266c6eSPeter Brune gs->max_its = 50; 3756cd56b8bSPeter Brune 3763542a6bcSPeter Brune snes->data = (void*) gs; 3773542a6bcSPeter Brune PetscFunctionReturn(0); 3783542a6bcSPeter Brune } 3793542a6bcSPeter Brune EXTERN_C_END 380