xref: /petsc/src/snes/impls/gs/snesgs.c (revision 88976e7148bb891b38778b9a5d5ae4c2b2cf44bf)
13542a6bcSPeter Brune #include <private/snesimpl.h>             /*I   "petscsnes.h"   I*/
23542a6bcSPeter Brune 
33542a6bcSPeter Brune typedef struct {
43542a6bcSPeter Brune   PetscInt  sweeps;
54b32a720SPeter Brune   PetscBool norms;
63542a6bcSPeter Brune } SNES_GS;
73542a6bcSPeter Brune 
83542a6bcSPeter Brune #undef __FUNCT__
93542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS"
103542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes)
113542a6bcSPeter Brune {
123542a6bcSPeter Brune   PetscFunctionBegin;
133542a6bcSPeter Brune   PetscFunctionReturn(0);
143542a6bcSPeter Brune }
153542a6bcSPeter Brune 
163542a6bcSPeter Brune #undef __FUNCT__
173542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS"
183542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes)
193542a6bcSPeter Brune {
203542a6bcSPeter Brune   PetscErrorCode ierr;
213542a6bcSPeter Brune 
223542a6bcSPeter Brune   PetscFunctionBegin;
233542a6bcSPeter Brune   ierr = SNESReset_GS(snes);CHKERRQ(ierr);
243542a6bcSPeter Brune   ierr = PetscFree(snes->data);
253542a6bcSPeter Brune   PetscFunctionReturn(0);
263542a6bcSPeter Brune }
273542a6bcSPeter Brune 
283542a6bcSPeter Brune #undef __FUNCT__
293542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS"
303542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes)
313542a6bcSPeter Brune {
323542a6bcSPeter Brune   PetscFunctionBegin;
333542a6bcSPeter Brune   PetscFunctionReturn(0);
343542a6bcSPeter Brune }
353542a6bcSPeter Brune 
363542a6bcSPeter Brune #undef __FUNCT__
373542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS"
383542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes)
393542a6bcSPeter Brune {
403542a6bcSPeter Brune   PetscErrorCode ierr;
414b32a720SPeter Brune   SNES_GS        *gs = (SNES_GS *)snes->data;
423542a6bcSPeter Brune 
433542a6bcSPeter Brune   PetscFunctionBegin;
443542a6bcSPeter Brune   ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr);
454b32a720SPeter Brune   ierr = PetscOptionsBool("-snes_gs_norms","Compute norms for monitoring","none",gs->norms,&gs->norms,PETSC_NULL);CHKERRQ(ierr);
464b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
473542a6bcSPeter Brune   PetscFunctionReturn(0);
483542a6bcSPeter Brune }
493542a6bcSPeter Brune 
503542a6bcSPeter Brune #undef __FUNCT__
513542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS"
523542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
533542a6bcSPeter Brune {
543542a6bcSPeter Brune   PetscFunctionBegin;
553542a6bcSPeter Brune   PetscFunctionReturn(0);
563542a6bcSPeter Brune }
573542a6bcSPeter Brune 
583542a6bcSPeter Brune #undef __FUNCT__
593542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS"
603542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes)
613542a6bcSPeter Brune {
623542a6bcSPeter Brune   Vec            F;
633542a6bcSPeter Brune   Vec            X;
643542a6bcSPeter Brune   Vec            B;
653542a6bcSPeter Brune   PetscInt       i;
6606e07b1aSPeter Brune   PetscReal      fnorm;
673542a6bcSPeter Brune   PetscErrorCode ierr;
684b32a720SPeter Brune   SNES_GS        *gs = (SNES_GS *)snes->data;
693542a6bcSPeter Brune 
703542a6bcSPeter Brune   PetscFunctionBegin;
713542a6bcSPeter Brune 
723542a6bcSPeter Brune   X = snes->vec_sol;
733542a6bcSPeter Brune   F = snes->vec_func;
743542a6bcSPeter Brune   B = snes->vec_rhs;
753542a6bcSPeter Brune 
763542a6bcSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
773542a6bcSPeter Brune   snes->iter = 0;
783542a6bcSPeter Brune   snes->norm = 0.;
793542a6bcSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
803542a6bcSPeter Brune 
81526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
823542a6bcSPeter Brune   /* compute the initial function and preconditioned update delX */
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   }
883542a6bcSPeter Brune   /* convergence test */
893542a6bcSPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
903542a6bcSPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
913542a6bcSPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
923542a6bcSPeter Brune   snes->norm = fnorm;
933542a6bcSPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
943542a6bcSPeter Brune   SNESLogConvHistory(snes,fnorm,0);
953542a6bcSPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
963542a6bcSPeter Brune 
973542a6bcSPeter Brune   /* set parameter for default relative tolerance convergence test */
983542a6bcSPeter Brune   snes->ttol = fnorm*snes->rtol;
993542a6bcSPeter Brune   /* test convergence */
1003542a6bcSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1013542a6bcSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1023542a6bcSPeter Brune 
1033542a6bcSPeter Brune   /* Call general purpose update function */
1043542a6bcSPeter Brune   if (snes->ops->update) {
1053542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1063542a6bcSPeter Brune   }
1073542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
1083542a6bcSPeter Brune     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
1094b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
1104b32a720SPeter Brune     if (gs->norms || i == snes->max_its - 1) {
1113542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1123542a6bcSPeter Brune       if (snes->domainerror) {
1133542a6bcSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1143542a6bcSPeter Brune         PetscFunctionReturn(0);
1153542a6bcSPeter Brune       }
1163542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1173542a6bcSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1183542a6bcSPeter Brune       /* Monitor convergence */
1193542a6bcSPeter Brune       ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1203542a6bcSPeter Brune       snes->iter = i+1;
1213542a6bcSPeter Brune       snes->norm = fnorm;
1223542a6bcSPeter Brune       ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1236165829cSJed Brown     }
1243542a6bcSPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1253542a6bcSPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1263542a6bcSPeter Brune     /* Test for convergence */
1273542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1283542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
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     }
1343542a6bcSPeter Brune   }
135526b802eSJed Brown   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
1363542a6bcSPeter Brune   PetscFunctionReturn(0);
1373542a6bcSPeter Brune }
1383542a6bcSPeter Brune 
1393542a6bcSPeter Brune /*MC
1403542a6bcSPeter Brune   SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
1413542a6bcSPeter Brune 
1423542a6bcSPeter Brune    Level: advanced
1433542a6bcSPeter Brune 
1443542a6bcSPeter Brune   Notes:
1453542a6bcSPeter Brune   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
1463542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
1473542a6bcSPeter Brune 
1483542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
1493542a6bcSPeter Brune M*/
1503542a6bcSPeter Brune 
1513542a6bcSPeter Brune EXTERN_C_BEGIN
1523542a6bcSPeter Brune #undef __FUNCT__
1533542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS"
1543542a6bcSPeter Brune PetscErrorCode SNESCreate_GS(SNES snes)
1553542a6bcSPeter Brune {
1563542a6bcSPeter Brune   SNES_GS        *gs;
1573542a6bcSPeter Brune   PetscErrorCode ierr;
1583542a6bcSPeter Brune 
1593542a6bcSPeter Brune   PetscFunctionBegin;
1603542a6bcSPeter Brune   snes->ops->destroy        = SNESDestroy_GS;
1613542a6bcSPeter Brune   snes->ops->setup          = SNESSetUp_GS;
1623542a6bcSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_GS;
1633542a6bcSPeter Brune   snes->ops->view           = SNESView_GS;
1643542a6bcSPeter Brune   snes->ops->solve          = SNESSolve_GS;
1653542a6bcSPeter Brune   snes->ops->reset          = SNESReset_GS;
1663542a6bcSPeter Brune 
1673542a6bcSPeter Brune   snes->usesksp             = PETSC_FALSE;
1683542a6bcSPeter Brune   snes->usespc              = PETSC_FALSE;
1693542a6bcSPeter Brune 
170*88976e71SPeter Brune   if (!snes->tolerancesset) {
1710e444f03SPeter Brune     snes->max_its             = 10000;
1720e444f03SPeter Brune     snes->max_funcs           = 10000;
173*88976e71SPeter Brune   }
1740e444f03SPeter Brune 
1753542a6bcSPeter Brune   ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr);
1763542a6bcSPeter Brune   snes->data = (void*) gs;
1774b32a720SPeter Brune   gs->norms = PETSC_FALSE;
1783542a6bcSPeter Brune   PetscFunctionReturn(0);
1793542a6bcSPeter Brune }
1803542a6bcSPeter Brune EXTERN_C_END
181