13542a6bcSPeter Brune #include <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); 433542a6bcSPeter Brune PetscFunctionReturn(0); 443542a6bcSPeter Brune } 453542a6bcSPeter Brune 463542a6bcSPeter Brune #undef __FUNCT__ 473542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS" 483542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer) 493542a6bcSPeter Brune { 503542a6bcSPeter Brune PetscFunctionBegin; 513542a6bcSPeter Brune PetscFunctionReturn(0); 523542a6bcSPeter Brune } 533542a6bcSPeter Brune 543542a6bcSPeter Brune #undef __FUNCT__ 553542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS" 563542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes) 573542a6bcSPeter Brune { 583542a6bcSPeter Brune Vec F; 593542a6bcSPeter Brune Vec X; 603542a6bcSPeter Brune Vec B; 613542a6bcSPeter Brune PetscInt i; 62*06e07b1aSPeter Brune PetscReal fnorm; 633542a6bcSPeter Brune PetscErrorCode ierr; 643542a6bcSPeter Brune 653542a6bcSPeter Brune PetscFunctionBegin; 663542a6bcSPeter Brune 673542a6bcSPeter Brune X = snes->vec_sol; 683542a6bcSPeter Brune F = snes->vec_func; 693542a6bcSPeter Brune B = snes->vec_rhs; 703542a6bcSPeter Brune 713542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 723542a6bcSPeter Brune snes->iter = 0; 733542a6bcSPeter Brune snes->norm = 0.; 743542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 753542a6bcSPeter Brune 763542a6bcSPeter Brune snes->reason = SNES_CONVERGED_ITS; 773542a6bcSPeter Brune /* compute the initial function and preconditioned update delX */ 783542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 793542a6bcSPeter Brune if (snes->domainerror) { 803542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 813542a6bcSPeter Brune PetscFunctionReturn(0); 823542a6bcSPeter Brune } 833542a6bcSPeter Brune /* convergence test */ 843542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 853542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 863542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 873542a6bcSPeter Brune snes->norm = fnorm; 883542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 893542a6bcSPeter Brune SNESLogConvHistory(snes,fnorm,0); 903542a6bcSPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 913542a6bcSPeter Brune 923542a6bcSPeter Brune /* set parameter for default relative tolerance convergence test */ 933542a6bcSPeter Brune snes->ttol = fnorm*snes->rtol; 943542a6bcSPeter Brune /* test convergence */ 953542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 963542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 973542a6bcSPeter Brune 983542a6bcSPeter Brune /* Call general purpose update function */ 993542a6bcSPeter Brune if (snes->ops->update) { 1003542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1013542a6bcSPeter Brune } 1023542a6bcSPeter Brune for(i = 0; i < snes->max_its; i++) { 1033542a6bcSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 1043542a6bcSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1053542a6bcSPeter Brune if (snes->domainerror) { 1063542a6bcSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1073542a6bcSPeter Brune PetscFunctionReturn(0); 1083542a6bcSPeter Brune } 1093542a6bcSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1103542a6bcSPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1113542a6bcSPeter Brune /* Monitor convergence */ 1123542a6bcSPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1133542a6bcSPeter Brune snes->iter = i+1; 1143542a6bcSPeter Brune snes->norm = fnorm; 1153542a6bcSPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1163542a6bcSPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1173542a6bcSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1183542a6bcSPeter Brune 1193542a6bcSPeter Brune /* Test for convergence */ 1203542a6bcSPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1213542a6bcSPeter Brune if (snes->reason) PetscFunctionReturn(0); 1223542a6bcSPeter Brune 1233542a6bcSPeter Brune /* Call general purpose update function */ 1243542a6bcSPeter Brune if (snes->ops->update) { 1253542a6bcSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1263542a6bcSPeter Brune } 1273542a6bcSPeter Brune } 1283542a6bcSPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);CHKERRQ(ierr); 1293542a6bcSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1303542a6bcSPeter Brune PetscFunctionReturn(0); 1313542a6bcSPeter Brune } 1323542a6bcSPeter Brune 1333542a6bcSPeter Brune /*MC 1343542a6bcSPeter Brune SNESGS - Just calls the user-provided solution routine provided with SNESSetGS() 1353542a6bcSPeter Brune 1363542a6bcSPeter Brune Level: advanced 1373542a6bcSPeter Brune 1383542a6bcSPeter Brune Notes: 1393542a6bcSPeter Brune the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have 1403542a6bcSPeter Brune its parent's Gauss-Seidel routine associated with it. 1413542a6bcSPeter Brune 1423542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types) 1433542a6bcSPeter Brune M*/ 1443542a6bcSPeter Brune 1453542a6bcSPeter Brune EXTERN_C_BEGIN 1463542a6bcSPeter Brune #undef __FUNCT__ 1473542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS" 1483542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes) 1493542a6bcSPeter Brune { 1503542a6bcSPeter Brune SNES_GS *gs; 1513542a6bcSPeter Brune PetscErrorCode ierr; 1523542a6bcSPeter Brune 1533542a6bcSPeter Brune PetscFunctionBegin; 1543542a6bcSPeter Brune snes->ops->destroy = SNESDestroy_GS; 1553542a6bcSPeter Brune snes->ops->setup = SNESSetUp_GS; 1563542a6bcSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_GS; 1573542a6bcSPeter Brune snes->ops->view = SNESView_GS; 1583542a6bcSPeter Brune snes->ops->solve = SNESSolve_GS; 1593542a6bcSPeter Brune snes->ops->reset = SNESReset_GS; 1603542a6bcSPeter Brune 1613542a6bcSPeter Brune snes->usesksp = PETSC_FALSE; 1623542a6bcSPeter Brune snes->usespc = PETSC_FALSE; 1633542a6bcSPeter Brune 1643542a6bcSPeter Brune ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr); 1653542a6bcSPeter Brune snes->data = (void*) gs; 1663542a6bcSPeter Brune PetscFunctionReturn(0); 1673542a6bcSPeter Brune } 1683542a6bcSPeter Brune EXTERN_C_END 169