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__ 83542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS" 93542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes) 103542a6bcSPeter Brune { 113542a6bcSPeter Brune PetscFunctionBegin; 123542a6bcSPeter Brune PetscFunctionReturn(0); 133542a6bcSPeter Brune } 143542a6bcSPeter Brune 153542a6bcSPeter Brune #undef __FUNCT__ 163542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS" 173542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes) 183542a6bcSPeter Brune { 193542a6bcSPeter Brune PetscErrorCode ierr; 203542a6bcSPeter Brune 213542a6bcSPeter Brune PetscFunctionBegin; 223542a6bcSPeter Brune ierr = SNESReset_GS(snes);CHKERRQ(ierr); 233542a6bcSPeter Brune ierr = PetscFree(snes->data); 243542a6bcSPeter Brune PetscFunctionReturn(0); 253542a6bcSPeter Brune } 263542a6bcSPeter Brune 273542a6bcSPeter Brune #undef __FUNCT__ 283542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS" 293542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes) 303542a6bcSPeter Brune { 313542a6bcSPeter Brune PetscFunctionBegin; 323542a6bcSPeter Brune PetscFunctionReturn(0); 333542a6bcSPeter Brune } 343542a6bcSPeter Brune 353542a6bcSPeter Brune #undef __FUNCT__ 363542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS" 373542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes) 383542a6bcSPeter Brune { 393542a6bcSPeter Brune PetscErrorCode ierr; 403542a6bcSPeter Brune 413542a6bcSPeter Brune PetscFunctionBegin; 423542a6bcSPeter Brune ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr); 434b32a720SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 443542a6bcSPeter Brune PetscFunctionReturn(0); 453542a6bcSPeter Brune } 463542a6bcSPeter Brune 473542a6bcSPeter Brune #undef __FUNCT__ 483542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS" 493542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 503542a6bcSPeter Brune { 513542a6bcSPeter Brune PetscFunctionBegin; 523542a6bcSPeter Brune PetscFunctionReturn(0); 533542a6bcSPeter Brune } 543542a6bcSPeter Brune 553542a6bcSPeter Brune #undef __FUNCT__ 563542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS" 573542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes) 583542a6bcSPeter Brune { 593542a6bcSPeter Brune Vec F; 603542a6bcSPeter Brune Vec X; 613542a6bcSPeter Brune Vec B; 623542a6bcSPeter Brune PetscInt i; 6306e07b1aSPeter Brune PetscReal fnorm; 643542a6bcSPeter Brune PetscErrorCode ierr; 65*534ebe21SPeter Brune SNESNormType normtype; 663542a6bcSPeter Brune 673542a6bcSPeter Brune PetscFunctionBegin; 683542a6bcSPeter Brune 693542a6bcSPeter Brune X = snes->vec_sol; 703542a6bcSPeter Brune F = snes->vec_func; 713542a6bcSPeter Brune B = snes->vec_rhs; 723542a6bcSPeter Brune 73*534ebe21SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 74*534ebe21SPeter Brune snes->iter = 0; 75*534ebe21SPeter Brune snes->norm = 0.; 76*534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 77526b802eSJed Brown snes->reason = SNES_CONVERGED_ITERATING; 78*534ebe21SPeter Brune 79*534ebe21SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 80fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 813542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 82e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 833542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 843542a6bcSPeter Brune if (snes->domainerror) { 853542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 863542a6bcSPeter Brune PetscFunctionReturn(0); 873542a6bcSPeter Brune } 88e4ed7901SPeter Brune } else { 89e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 90e4ed7901SPeter Brune } 91e4ed7901SPeter Brune 923542a6bcSPeter Brune /* convergence test */ 93e4ed7901SPeter Brune if (!snes->norm_init_set) { 943542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 953542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 96e4ed7901SPeter Brune } else { 97e4ed7901SPeter Brune fnorm = snes->norm_init; 98e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 99e4ed7901SPeter Brune } 1003542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 101e4ed7901SPeter Brune snes->iter = 0; 1023542a6bcSPeter Brune snes->norm = fnorm; 1033542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 104e4ed7901SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 105e4ed7901SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 1063542a6bcSPeter Brune 1073542a6bcSPeter Brune /* set parameter for default relative tolerance convergence test */ 1083542a6bcSPeter Brune snes->ttol = fnorm*snes->rtol; 109*534ebe21SPeter Brune 1103542a6bcSPeter Brune /* test convergence */ 1113542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1123542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 113*534ebe21SPeter Brune } else { 114*534ebe21SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 115*534ebe21SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 116*534ebe21SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 117*534ebe21SPeter Brune } 118*534ebe21SPeter Brune 1193542a6bcSPeter Brune 1203542a6bcSPeter Brune /* Call general purpose update function */ 1213542a6bcSPeter Brune if (snes->ops->update) { 1223542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1233542a6bcSPeter Brune } 124*534ebe21SPeter Brune 1253542a6bcSPeter Brune for (i = 0; i < snes->max_its; i++) { 1263542a6bcSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 1274b32a720SPeter Brune /* only compute norms if requested or about to exit due to maximum iterations */ 128fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 1293542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1303542a6bcSPeter Brune if (snes->domainerror) { 1313542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1323542a6bcSPeter Brune PetscFunctionReturn(0); 1333542a6bcSPeter Brune } 1343542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1353542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 136*534ebe21SPeter Brune } 1373542a6bcSPeter Brune /* Monitor convergence */ 1383542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1393542a6bcSPeter Brune snes->iter = i+1; 1403542a6bcSPeter Brune snes->norm = fnorm; 1413542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1423542a6bcSPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1433542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1443542a6bcSPeter Brune /* Test for convergence */ 145fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1463542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1473542a6bcSPeter Brune /* Call general purpose update function */ 1483542a6bcSPeter Brune if (snes->ops->update) { 1493542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1503542a6bcSPeter Brune } 1513542a6bcSPeter Brune } 152fdacfa88SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 153*534ebe21SPeter Brune if (i == snes->max_its) { 154*534ebe21SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 155*534ebe21SPeter Brune if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 156*534ebe21SPeter Brune } 157*534ebe21SPeter Brune } else { 158526b802eSJed Brown if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */ 159*534ebe21SPeter Brune } 1603542a6bcSPeter Brune PetscFunctionReturn(0); 1613542a6bcSPeter Brune } 1623542a6bcSPeter Brune 1633542a6bcSPeter Brune /*MC 1643542a6bcSPeter Brune SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 1653542a6bcSPeter Brune 1663542a6bcSPeter Brune Level: advanced 1673542a6bcSPeter Brune 1683542a6bcSPeter Brune Notes: 1693542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 1703542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 1713542a6bcSPeter Brune 1723542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 1733542a6bcSPeter Brune M*/ 1743542a6bcSPeter Brune 1753542a6bcSPeter Brune EXTERN_C_BEGIN 1763542a6bcSPeter Brune #undef __FUNCT__ 1773542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS" 1783542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes) 1793542a6bcSPeter Brune { 1803542a6bcSPeter Brune SNES_GS *gs; 1813542a6bcSPeter Brune PetscErrorCode ierr; 1823542a6bcSPeter Brune 1833542a6bcSPeter Brune PetscFunctionBegin; 1843542a6bcSPeter Brune snes->ops->destroy = SNESDestroy_GS; 1853542a6bcSPeter Brune snes->ops->setup = SNESSetUp_GS; 1863542a6bcSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_GS; 1873542a6bcSPeter Brune snes->ops->view = SNESView_GS; 1883542a6bcSPeter Brune snes->ops->solve = SNESSolve_GS; 1893542a6bcSPeter Brune snes->ops->reset = SNESReset_GS; 1903542a6bcSPeter Brune 1913542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 1923542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 1933542a6bcSPeter Brune 19488976e71SPeter Brune if (!snes->tolerancesset) { 1950e444f03SPeter Brune snes->max_its = 10000; 1960e444f03SPeter Brune snes->max_funcs = 10000; 19788976e71SPeter Brune } 1980e444f03SPeter Brune 1993542a6bcSPeter Brune ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 2003542a6bcSPeter Brune snes->data = (void*) gs; 2013542a6bcSPeter Brune PetscFunctionReturn(0); 2023542a6bcSPeter Brune } 2033542a6bcSPeter Brune EXTERN_C_END 204