/[MITgcm]/MITgcm/pkg/seaice/seaice_jfnk.F
ViewVC logotype

Diff of /MITgcm/pkg/seaice/seaice_jfnk.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.20 by jmc, Sat Mar 2 04:35:05 2013 UTC revision 1.23 by mlosch, Thu May 30 14:07:19 2013 UTC
# Line 53  C     myThid :: my Thread Id. number Line 53  C     myThid :: my Thread Id. number
53        INTEGER myIter        INTEGER myIter
54        INTEGER myThid        INTEGER myThid
55    
56  #if ( (defined SEAICE_CGRID) && \  #ifdef SEAICE_ALLOW_JFNK
       (defined SEAICE_ALLOW_JFNK) && \  
       (defined SEAICE_ALLOW_DYNAMICS) )  
57  C     !FUNCTIONS:  C     !FUNCTIONS:
58        LOGICAL  DIFFERENT_MULTIPLE        LOGICAL  DIFFERENT_MULTIPLE
59        EXTERNAL DIFFERENT_MULTIPLE        EXTERNAL DIFFERENT_MULTIPLE
# Line 75  C     FGMRES flag that indicates what fg Line 73  C     FGMRES flag that indicates what fg
73        _RL     JFNKresidual        _RL     JFNKresidual
74        _RL     JFNKresidualKm1        _RL     JFNKresidualKm1
75  C     parameters to compute convergence criterion  C     parameters to compute convergence criterion
76        _RL     phi_e, alp_e, JFNKgamma_lin        _RL     JFNKgamma_lin
77        _RL     FGMRESeps        _RL     FGMRESeps
78        _RL     JFNKtol        _RL     JFNKtol
79    C     Adams-Bashforth extrapolation factors
80          _RL abFac, abAlpha
81    C
82        _RL     recip_deltaT        _RL     recip_deltaT
83        LOGICAL JFNKconverged, krylovConverged        LOGICAL JFNKconverged, krylovConverged
84        LOGICAL writeNow        LOGICAL writeNow
# Line 87  C     parameters to compute convergence Line 87  C     parameters to compute convergence
87  C     u/vIceRes :: residual of sea-ice momentum equations  C     u/vIceRes :: residual of sea-ice momentum equations
88        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
89        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
90    C     extra time level required for Adams-Bashforth-2 time stepping
91          _RL duIcNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
92          _RL dvIcNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
93  C     du/vIce   :: ice velocity increment to be added to u/vIce  C     du/vIce   :: ice velocity increment to be added to u/vIce
94        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
95        _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
# Line 115  C     with iOutFgmres=1, seaice_fgmres p Line 118  C     with iOutFgmres=1, seaice_fgmres p
118       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
119       &     iOutFGMRES=1       &     iOutFGMRES=1
120    
121    C     Adams-Bashforth extrapolation factors
122          abFac = 0. _d 0
123          IF ( SEAICEuseAB2 ) THEN
124           IF ( myIter.EQ.nIter0 .AND. SEAICEmomStartAB.EQ.0 ) THEN
125            abFac = 0. _d 0
126           ELSE
127            abFac = 0.5 _d 0 + SEAICE_abEps
128           ENDIF
129          ENDIF
130          abAlpha = 1. _d 0 + abFac
131    
132        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
133         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
134          DO J=1-OLy,sNy+OLy          DO J=1-OLy,sNy+OLy
# Line 123  C     with iOutFgmres=1, seaice_fgmres p Line 137  C     with iOutFgmres=1, seaice_fgmres p
137            vIceRes(I,J,bi,bj) = 0. _d 0            vIceRes(I,J,bi,bj) = 0. _d 0
138            duIce  (I,J,bi,bj) = 0. _d 0            duIce  (I,J,bi,bj) = 0. _d 0
139            dvIce  (I,J,bi,bj) = 0. _d 0            dvIce  (I,J,bi,bj) = 0. _d 0
140             ENDDO
141            ENDDO
142    C     cycle ice velocities
143            DO J=1-OLy,sNy+OLy
144             DO I=1-OLx,sNx+OLx
145              duIcNm1(I,J,bi,bj) = uIce(I,J,bi,bj) * abAlpha
146         &         + ( uIce(I,J,bi,bj) - uIceNm1(I,J,bi,bj) ) * abFac
147              dvIcNm1(I,J,bi,bj) = vIce(I,J,bi,bj) * abAlpha
148         &         + ( vIce(I,J,bi,bj) - vIceNm1(I,J,bi,bj) ) * abFac
149            uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)            uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
150            vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)            vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
151           ENDDO           ENDDO
152          ENDDO          ENDDO
153            IF ( .NOT.SEAICEuseIMEX ) THEN
154  C     Compute things that do no change during the Newton iteration:  C     Compute things that do no change during the Newton iteration:
155  C     sea-surface tilt and wind stress:  C     sea-surface tilt and wind stress:
156  C     FORCEX/Y0 - mass*(u/vIceNm1)/deltaT  C     FORCEX/Y0 - mass*(abA*u/vIceNm1+abB*(u/vIceNm1-u/vIceNm2))/deltaT
157          DO J=1-OLy,sNy+OLy          DO J=1-OLy,sNy+OLy
158           DO I=1-OLx,sNx+OLx           DO I=1-OLx,sNx+OLx
159            FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)            FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
160       &         + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT       &         + seaiceMassU(I,J,bi,bj)*duIcNm1(I,J,bi,bj)*recip_deltaT
161            FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)            FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
162       &         + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT       &         + seaiceMassV(I,J,bi,bj)*dvIcNm1(I,J,bi,bj)*recip_deltaT
163           ENDDO           ENDDO
164          ENDDO          ENDDO
165            ENDIF
166         ENDDO         ENDDO
167        ENDDO        ENDDO
168  C     Start nonlinear Newton iteration: outer loop iteration  C     Start nonlinear Newton iteration: outer loop iteration
# Line 169  C     compute convergence criterion for Line 194  C     compute convergence criterion for
194         JFNKgamma_lin = JFNKgamma_lin_max         JFNKgamma_lin = JFNKgamma_lin_max
195         IF ( newtonIter.GT.1.AND.newtonIter.LE.SEAICE_JFNK_tolIter         IF ( newtonIter.GT.1.AND.newtonIter.LE.SEAICE_JFNK_tolIter
196       &      .AND.JFNKresidual.LT.JFNKres_t ) THEN       &      .AND.JFNKresidual.LT.JFNKres_t ) THEN
197  C     Eisenstat, 1996, equ.(2.6)  C     Eisenstat and Walker (1996), eq.(2.6)
198          phi_e = 1. _d 0          JFNKgamma_lin = SEAICE_JFNKphi
199          alp_e = 1. _d 0       &       *( JFNKresidual/JFNKresidualKm1 )**SEAICE_JFNKalpha
         JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e  
200          JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)          JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
201          JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)          JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
202         ENDIF         ENDIF
# Line 529  C     iterations, 0.25*du/vIce in the se Line 553  C     iterations, 0.25*du/vIce in the se
553  C     This is the new residual  C     This is the new residual
554        JFNKresidual = resLoc        JFNKresidual = resLoc
555    
556  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */  #endif /* SEAICE_ALLOW_JFNK */
557    
558        RETURN        RETURN
559        END        END

Legend:
Removed from v.1.20  
changed lines
  Added in v.1.23

  ViewVC Help
Powered by ViewVC 1.1.22