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 73542a6bcSPeter Brune #undef __FUNCT__ 8*c79d6ae8SPeter Brune #define __FUNCT__ "SNESDefaultApplyGS" 9*c79d6ae8SPeter Brune PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx) { 10*c79d6ae8SPeter Brune DM dm; 11*c79d6ae8SPeter Brune PetscErrorCode ierr; 12*c79d6ae8SPeter Brune PetscFunctionBegin; 13*c79d6ae8SPeter Brune /* see if there's a coloring on the DM */ 14*c79d6ae8SPeter Brune PetscFunctionReturn(0); 15*c79d6ae8SPeter Brune } 16*c79d6ae8SPeter Brune 17*c79d6ae8SPeter Brune #undef __FUNCT__ 183542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS" 193542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes) 203542a6bcSPeter Brune { 213542a6bcSPeter Brune PetscFunctionBegin; 223542a6bcSPeter Brune PetscFunctionReturn(0); 233542a6bcSPeter Brune } 243542a6bcSPeter Brune 253542a6bcSPeter Brune #undef __FUNCT__ 263542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS" 273542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes) 283542a6bcSPeter Brune { 293542a6bcSPeter Brune PetscErrorCode ierr; 303542a6bcSPeter Brune 313542a6bcSPeter Brune PetscFunctionBegin; 323542a6bcSPeter Brune ierr = SNESReset_GS(snes);CHKERRQ(ierr); 333542a6bcSPeter Brune ierr = PetscFree(snes->data); 343542a6bcSPeter Brune PetscFunctionReturn(0); 353542a6bcSPeter Brune } 363542a6bcSPeter Brune 373542a6bcSPeter Brune #undef __FUNCT__ 383542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS" 393542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes) 403542a6bcSPeter Brune { 413542a6bcSPeter Brune PetscFunctionBegin; 423542a6bcSPeter Brune PetscFunctionReturn(0); 433542a6bcSPeter Brune } 443542a6bcSPeter Brune 453542a6bcSPeter Brune #undef __FUNCT__ 463542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS" 473542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes) 483542a6bcSPeter Brune { 493542a6bcSPeter Brune PetscErrorCode ierr; 503542a6bcSPeter Brune 513542a6bcSPeter Brune PetscFunctionBegin; 523542a6bcSPeter Brune ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr); 534b32a720SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 543542a6bcSPeter Brune PetscFunctionReturn(0); 553542a6bcSPeter Brune } 563542a6bcSPeter Brune 573542a6bcSPeter Brune #undef __FUNCT__ 583542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS" 593542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 603542a6bcSPeter Brune { 613542a6bcSPeter Brune PetscFunctionBegin; 623542a6bcSPeter Brune PetscFunctionReturn(0); 633542a6bcSPeter Brune } 643542a6bcSPeter Brune 653542a6bcSPeter Brune #undef __FUNCT__ 663542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS" 673542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes) 683542a6bcSPeter Brune { 693542a6bcSPeter Brune Vec F; 703542a6bcSPeter Brune Vec X; 713542a6bcSPeter Brune Vec B; 723542a6bcSPeter Brune PetscInt i; 7306e07b1aSPeter Brune PetscReal fnorm; 743542a6bcSPeter Brune PetscErrorCode ierr; 75534ebe21SPeter Brune SNESNormType normtype; 763542a6bcSPeter Brune 773542a6bcSPeter Brune PetscFunctionBegin; 783542a6bcSPeter Brune 793542a6bcSPeter Brune X = snes->vec_sol; 803542a6bcSPeter Brune F = snes->vec_func; 813542a6bcSPeter Brune B = snes->vec_rhs; 823542a6bcSPeter Brune 83534ebe21SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 84534ebe21SPeter Brune snes->iter = 0; 85534ebe21SPeter Brune snes->norm = 0.; 86534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 87526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 88534ebe21SPeter Brune 89534ebe21SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 90fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 913542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 92e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 933542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 943542a6bcSPeter Brune if (snes->domainerror) { 953542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 963542a6bcSPeter Brune PetscFunctionReturn(0); 973542a6bcSPeter Brune } 98e4ed7901SPeter Brune } else { 99e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 100e4ed7901SPeter Brune } 101e4ed7901SPeter Brune 1023542a6bcSPeter Brune /* convergence test */ 103e4ed7901SPeter Brune if (!snes->norm_init_set) { 1043542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1053542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 106e4ed7901SPeter Brune } else { 107e4ed7901SPeter Brune fnorm = snes->norm_init; 108e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 109e4ed7901SPeter Brune } 1103542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 111e4ed7901SPeter Brune snes->iter = 0; 1123542a6bcSPeter Brune snes->norm = fnorm; 1133542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 114e4ed7901SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 115e4ed7901SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 1163542a6bcSPeter Brune 1173542a6bcSPeter Brune /* set parameter for default relative tolerance convergence test */ 1183542a6bcSPeter Brune snes->ttol = fnorm*snes->rtol; 119534ebe21SPeter Brune 1203542a6bcSPeter Brune /* test convergence */ 1213542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1223542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 123534ebe21SPeter Brune } else { 124534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 125534ebe21SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 126534ebe21SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 127534ebe21SPeter Brune } 128534ebe21SPeter Brune 1293542a6bcSPeter Brune 1303542a6bcSPeter Brune /* Call general purpose update function */ 1313542a6bcSPeter Brune if (snes->ops->update) { 1323542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1333542a6bcSPeter Brune } 134534ebe21SPeter Brune 1353542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 1363542a6bcSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 1374b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 138fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 1393542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1403542a6bcSPeter Brune if (snes->domainerror) { 1413542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1423542a6bcSPeter Brune PetscFunctionReturn(0); 1433542a6bcSPeter Brune } 1443542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1453542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 146534ebe21SPeter Brune } 1473542a6bcSPeter Brune /* Monitor convergence */ 1483542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1493542a6bcSPeter Brune snes->iter = i+1; 1503542a6bcSPeter Brune snes->norm = fnorm; 1513542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1523542a6bcSPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1533542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1543542a6bcSPeter Brune /* Test for convergence */ 155fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1563542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1573542a6bcSPeter Brune /* Call general purpose update function */ 1583542a6bcSPeter Brune if (snes->ops->update) { 1593542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1603542a6bcSPeter Brune } 1613542a6bcSPeter Brune } 162fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 163534ebe21SPeter Brune if (i == snes->max_its) { 164534ebe21SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 165534ebe21SPeter Brune if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 166534ebe21SPeter Brune } 167534ebe21SPeter Brune } else { 168526b802eSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 169534ebe21SPeter Brune } 1703542a6bcSPeter Brune PetscFunctionReturn(0); 1713542a6bcSPeter Brune } 1723542a6bcSPeter Brune 1733542a6bcSPeter Brune /*MC 1743542a6bcSPeter Brune SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 1753542a6bcSPeter Brune 1763542a6bcSPeter Brune Level: advanced 1773542a6bcSPeter Brune 1783542a6bcSPeter Brune Notes: 1793542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 1803542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 1813542a6bcSPeter Brune 1823542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 1833542a6bcSPeter Brune M*/ 1843542a6bcSPeter Brune 1853542a6bcSPeter Brune EXTERN_C_BEGIN 1863542a6bcSPeter Brune #undef __FUNCT__ 1873542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS" 1883542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes) 1893542a6bcSPeter Brune { 1903542a6bcSPeter Brune SNES_GS *gs; 1913542a6bcSPeter Brune PetscErrorCode ierr; 1923542a6bcSPeter Brune 1933542a6bcSPeter Brune PetscFunctionBegin; 1943542a6bcSPeter Brune snes->ops->destroy = SNESDestroy_GS; 1953542a6bcSPeter Brune snes->ops->setup = SNESSetUp_GS; 1963542a6bcSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_GS; 1973542a6bcSPeter Brune snes->ops->view = SNESView_GS; 1983542a6bcSPeter Brune snes->ops->solve = SNESSolve_GS; 1993542a6bcSPeter Brune snes->ops->reset = SNESReset_GS; 2003542a6bcSPeter Brune 2013542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 2023542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 2033542a6bcSPeter Brune 20488976e71SPeter Brune if (!snes->tolerancesset) { 2050e444f03SPeter Brune snes->max_its = 10000; 2060e444f03SPeter Brune snes->max_funcs = 10000; 20788976e71SPeter Brune } 2080e444f03SPeter Brune 2093542a6bcSPeter Brune ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 2103542a6bcSPeter Brune snes->data = (void*) gs; 2113542a6bcSPeter Brune PetscFunctionReturn(0); 2123542a6bcSPeter Brune } 2133542a6bcSPeter Brune EXTERN_C_END 214