1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 23542a6bcSPeter Brune 33542a6bcSPeter Brune typedef struct { 43542a6bcSPeter Brune PetscInt sweeps; 53542a6bcSPeter Brune } SNES_GS; 63542a6bcSPeter Brune 7*6cd56b8bSPeter Brune 8*6cd56b8bSPeter Brune 9*6cd56b8bSPeter Brune #undef __FUNCT__ 10*6cd56b8bSPeter Brune #define __FUNCT__ "SNESGSSetSweeps" 11*6cd56b8bSPeter Brune /*@ 12*6cd56b8bSPeter Brune SNESGSSetSweeps - Sets the number of sweeps of GS to use. 13*6cd56b8bSPeter Brune 14*6cd56b8bSPeter Brune Input Parameters: 15*6cd56b8bSPeter Brune + snes - the SNES context 16*6cd56b8bSPeter Brune - sweeps - the number of sweeps of GS to perform. 17*6cd56b8bSPeter Brune 18*6cd56b8bSPeter Brune Level: intermediate 19*6cd56b8bSPeter Brune 20*6cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 21*6cd56b8bSPeter Brune 22*6cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps() 23*6cd56b8bSPeter Brune @*/ 24*6cd56b8bSPeter Brune 25*6cd56b8bSPeter Brune PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps) 26*6cd56b8bSPeter Brune { 27*6cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 28*6cd56b8bSPeter Brune 29*6cd56b8bSPeter Brune PetscFunctionBegin; 30*6cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31*6cd56b8bSPeter Brune gs->sweeps = sweeps; 32*6cd56b8bSPeter Brune PetscFunctionReturn(0); 33*6cd56b8bSPeter Brune } 34*6cd56b8bSPeter Brune 35*6cd56b8bSPeter Brune #undef __FUNCT__ 36*6cd56b8bSPeter Brune #define __FUNCT__ "SNESGSGetSweeps" 37*6cd56b8bSPeter Brune /*@ 38*6cd56b8bSPeter Brune SNESGSGetSweeps - Gets the number of sweeps GS will use. 39*6cd56b8bSPeter Brune 40*6cd56b8bSPeter Brune Input Parameters: 41*6cd56b8bSPeter Brune . snes - the SNES context 42*6cd56b8bSPeter Brune 43*6cd56b8bSPeter Brune Output Parameters: 44*6cd56b8bSPeter Brune . sweeps - the number of sweeps of GS to perform. 45*6cd56b8bSPeter Brune 46*6cd56b8bSPeter Brune Level: intermediate 47*6cd56b8bSPeter Brune 48*6cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel 49*6cd56b8bSPeter Brune 50*6cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps() 51*6cd56b8bSPeter Brune @*/ 52*6cd56b8bSPeter Brune PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps) 53*6cd56b8bSPeter Brune { 54*6cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 55*6cd56b8bSPeter Brune 56*6cd56b8bSPeter Brune PetscFunctionBegin; 57*6cd56b8bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 58*6cd56b8bSPeter Brune *sweeps = gs->sweeps; 59*6cd56b8bSPeter Brune PetscFunctionReturn(0); 60*6cd56b8bSPeter Brune } 61*6cd56b8bSPeter Brune 62*6cd56b8bSPeter Brune 633542a6bcSPeter Brune #undef __FUNCT__ 64c79d6ae8SPeter Brune #define __FUNCT__ "SNESDefaultApplyGS" 650adebc6cSBarry Smith PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx) 660adebc6cSBarry Smith { 67c79d6ae8SPeter Brune PetscFunctionBegin; 68c79d6ae8SPeter Brune /* see if there's a coloring on the DM */ 691bdf7e34SPeter Brune 70c79d6ae8SPeter Brune PetscFunctionReturn(0); 71c79d6ae8SPeter Brune } 72c79d6ae8SPeter Brune 73c79d6ae8SPeter Brune #undef __FUNCT__ 743542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS" 753542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes) 763542a6bcSPeter Brune { 773542a6bcSPeter Brune PetscFunctionBegin; 783542a6bcSPeter Brune PetscFunctionReturn(0); 793542a6bcSPeter Brune } 803542a6bcSPeter Brune 813542a6bcSPeter Brune #undef __FUNCT__ 823542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS" 833542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes) 843542a6bcSPeter Brune { 853542a6bcSPeter Brune PetscErrorCode ierr; 863542a6bcSPeter Brune 873542a6bcSPeter Brune PetscFunctionBegin; 883542a6bcSPeter Brune ierr = SNESReset_GS(snes);CHKERRQ(ierr); 8922d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 903542a6bcSPeter Brune PetscFunctionReturn(0); 913542a6bcSPeter Brune } 923542a6bcSPeter Brune 933542a6bcSPeter Brune #undef __FUNCT__ 943542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS" 953542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes) 963542a6bcSPeter Brune { 973542a6bcSPeter Brune PetscFunctionBegin; 983542a6bcSPeter Brune PetscFunctionReturn(0); 993542a6bcSPeter Brune } 1003542a6bcSPeter Brune 1013542a6bcSPeter Brune #undef __FUNCT__ 1023542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS" 1033542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes) 1043542a6bcSPeter Brune { 105*6cd56b8bSPeter Brune SNES_GS *gs = (SNES_GS*)snes->data; 1063542a6bcSPeter Brune PetscErrorCode ierr; 1073542a6bcSPeter Brune 1083542a6bcSPeter Brune PetscFunctionBegin; 1093542a6bcSPeter Brune ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr); 110*6cd56b8bSPeter Brune /* GS Options */ 111*6cd56b8bSPeter Brune ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&gs->sweeps,PETSC_NULL);CHKERRQ(ierr); 1124b32a720SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 1133542a6bcSPeter Brune PetscFunctionReturn(0); 1143542a6bcSPeter Brune } 1153542a6bcSPeter Brune 1163542a6bcSPeter Brune #undef __FUNCT__ 1173542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS" 1183542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 1193542a6bcSPeter Brune { 1203542a6bcSPeter Brune PetscFunctionBegin; 1213542a6bcSPeter Brune PetscFunctionReturn(0); 1223542a6bcSPeter Brune } 1233542a6bcSPeter Brune 1243542a6bcSPeter Brune #undef __FUNCT__ 1253542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS" 1263542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes) 1273542a6bcSPeter Brune { 1283542a6bcSPeter Brune Vec F; 1293542a6bcSPeter Brune Vec X; 1303542a6bcSPeter Brune Vec B; 1313542a6bcSPeter Brune PetscInt i; 13206e07b1aSPeter Brune PetscReal fnorm; 1333542a6bcSPeter Brune PetscErrorCode ierr; 134534ebe21SPeter Brune SNESNormType normtype; 1353542a6bcSPeter Brune 1363542a6bcSPeter Brune PetscFunctionBegin; 1373542a6bcSPeter Brune 1383542a6bcSPeter Brune X = snes->vec_sol; 1393542a6bcSPeter Brune F = snes->vec_func; 1403542a6bcSPeter Brune B = snes->vec_rhs; 1413542a6bcSPeter Brune 142534ebe21SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 143534ebe21SPeter Brune snes->iter = 0; 144534ebe21SPeter Brune snes->norm = 0.; 145534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 146526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 147534ebe21SPeter Brune 148534ebe21SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 149fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 1503542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 151e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1523542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1533542a6bcSPeter Brune if (snes->domainerror) { 1543542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1553542a6bcSPeter Brune PetscFunctionReturn(0); 1563542a6bcSPeter Brune } 157e4ed7901SPeter Brune } else { 158e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 159e4ed7901SPeter Brune } 160e4ed7901SPeter Brune 1613542a6bcSPeter Brune /* convergence test */ 162e4ed7901SPeter Brune if (!snes->norm_init_set) { 1633542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1643542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 165e4ed7901SPeter Brune } else { 166e4ed7901SPeter Brune fnorm = snes->norm_init; 167e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 168e4ed7901SPeter Brune } 1693542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 170e4ed7901SPeter Brune snes->iter = 0; 1713542a6bcSPeter Brune snes->norm = fnorm; 1723542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 173e4ed7901SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 174e4ed7901SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 1753542a6bcSPeter Brune 1763542a6bcSPeter Brune /* set parameter for default relative tolerance convergence test */ 1773542a6bcSPeter Brune snes->ttol = fnorm*snes->rtol; 178534ebe21SPeter Brune 1793542a6bcSPeter Brune /* test convergence */ 1803542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1813542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 182534ebe21SPeter Brune } else { 183534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 184534ebe21SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 185534ebe21SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 186534ebe21SPeter Brune } 187534ebe21SPeter Brune 1883542a6bcSPeter Brune 1893542a6bcSPeter Brune /* Call general purpose update function */ 1903542a6bcSPeter Brune if (snes->ops->update) { 1913542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1923542a6bcSPeter Brune } 193534ebe21SPeter Brune 1943542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 1953542a6bcSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 1964b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 197fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 1983542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1993542a6bcSPeter Brune if (snes->domainerror) { 2003542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2013542a6bcSPeter Brune PetscFunctionReturn(0); 2023542a6bcSPeter Brune } 2033542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 2043542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 205534ebe21SPeter Brune } 2063542a6bcSPeter Brune /* Monitor convergence */ 2073542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 2083542a6bcSPeter Brune snes->iter = i+1; 2093542a6bcSPeter Brune snes->norm = fnorm; 2103542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 2113542a6bcSPeter Brune SNESLogConvHistory(snes,snes->norm,0); 2123542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 2133542a6bcSPeter Brune /* Test for convergence */ 214fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 2153542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 2163542a6bcSPeter Brune /* Call general purpose update function */ 2173542a6bcSPeter Brune if (snes->ops->update) { 2183542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 2193542a6bcSPeter Brune } 2203542a6bcSPeter Brune } 221fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 222534ebe21SPeter Brune if (i == snes->max_its) { 223534ebe21SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 224534ebe21SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 225534ebe21SPeter Brune } 226534ebe21SPeter Brune } else { 227526b802eSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 228534ebe21SPeter Brune } 2293542a6bcSPeter Brune PetscFunctionReturn(0); 2303542a6bcSPeter Brune } 2313542a6bcSPeter Brune 2323542a6bcSPeter Brune /*MC 2333542a6bcSPeter Brune SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 2343542a6bcSPeter Brune 2353542a6bcSPeter Brune Level: advanced 2363542a6bcSPeter Brune 2373542a6bcSPeter Brune Notes: 2383542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 2393542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 2403542a6bcSPeter Brune 2413542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 2423542a6bcSPeter Brune M*/ 2433542a6bcSPeter Brune 2443542a6bcSPeter Brune EXTERN_C_BEGIN 2453542a6bcSPeter Brune #undef __FUNCT__ 2463542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS" 2473542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes) 2483542a6bcSPeter Brune { 2493542a6bcSPeter Brune SNES_GS *gs; 2503542a6bcSPeter Brune PetscErrorCode ierr; 2513542a6bcSPeter Brune 2523542a6bcSPeter Brune PetscFunctionBegin; 2533542a6bcSPeter Brune snes->ops->destroy = SNESDestroy_GS; 2543542a6bcSPeter Brune snes->ops->setup = SNESSetUp_GS; 2553542a6bcSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_GS; 2563542a6bcSPeter Brune snes->ops->view = SNESView_GS; 2573542a6bcSPeter Brune snes->ops->solve = SNESSolve_GS; 2583542a6bcSPeter Brune snes->ops->reset = SNESReset_GS; 2593542a6bcSPeter Brune 2603542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 2613542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 2623542a6bcSPeter Brune 26388976e71SPeter Brune if (!snes->tolerancesset) { 2640e444f03SPeter Brune snes->max_its = 10000; 2650e444f03SPeter Brune snes->max_funcs = 10000; 26688976e71SPeter Brune } 2670e444f03SPeter Brune 2683542a6bcSPeter Brune ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 269*6cd56b8bSPeter Brune 270*6cd56b8bSPeter Brune gs->sweeps = 1; 271*6cd56b8bSPeter Brune 2723542a6bcSPeter Brune snes->data = (void*) gs; 273*6cd56b8bSPeter Brune 2743542a6bcSPeter Brune PetscFunctionReturn(0); 2753542a6bcSPeter Brune } 2763542a6bcSPeter Brune EXTERN_C_END 277