/[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.16 by mlosch, Thu Jan 17 08:51:15 2013 UTC revision 1.30 by mlosch, Thu Jan 28 12:54:12 2016 UTC
# Line 2  C $Header$ Line 2  C $Header$
2  C $Name$  C $Name$
3    
4  #include "SEAICE_OPTIONS.h"  #include "SEAICE_OPTIONS.h"
5    #ifdef ALLOW_AUTODIFF
6    # include "AUTODIFF_OPTIONS.h"
7    #endif
8    
9  C--  File seaice_jfnk.F: seaice jfnk dynamical solver S/R:  C--  File seaice_jfnk.F: seaice jfnk dynamical solver S/R:
10  C--   Contents  C--   Contents
# Line 53  C     myThid :: my Thread Id. number Line 56  C     myThid :: my Thread Id. number
56        INTEGER myIter        INTEGER myIter
57        INTEGER myThid        INTEGER myThid
58    
59  #if ( (defined SEAICE_CGRID) && \  #ifdef SEAICE_ALLOW_JFNK
       (defined SEAICE_ALLOW_JFNK) && \  
       (defined SEAICE_ALLOW_DYNAMICS) )  
60  C     !FUNCTIONS:  C     !FUNCTIONS:
61        LOGICAL  DIFFERENT_MULTIPLE        LOGICAL  DIFFERENT_MULTIPLE
62        EXTERNAL DIFFERENT_MULTIPLE        EXTERNAL DIFFERENT_MULTIPLE
# Line 68  C     loop indices Line 69  C     loop indices
69        INTEGER newtonIter        INTEGER newtonIter
70        INTEGER krylovIter, krylovFails        INTEGER krylovIter, krylovFails
71        INTEGER totalKrylovItersLoc, totalNewtonItersLoc        INTEGER totalKrylovItersLoc, totalNewtonItersLoc
72    C     FGMRES parameters
73    C     im      :: size of Krylov space
74    C     ifgmres :: interation counter
75          INTEGER im
76          PARAMETER ( im = 50 )
77          INTEGER ifgmres
78  C     FGMRES flag that determines amount of output messages of fgmres  C     FGMRES flag that determines amount of output messages of fgmres
79        INTEGER iOutFGMRES        INTEGER iOutFGMRES
80  C     FGMRES flag that indicates what fgmres wants us to do next  C     FGMRES flag that indicates what fgmres wants us to do next
# Line 75  C     FGMRES flag that indicates what fg Line 82  C     FGMRES flag that indicates what fg
82        _RL     JFNKresidual        _RL     JFNKresidual
83        _RL     JFNKresidualKm1        _RL     JFNKresidualKm1
84  C     parameters to compute convergence criterion  C     parameters to compute convergence criterion
85        _RL     phi_e, alp_e, JFNKgamma_lin        _RL     JFNKgamma_lin
86        _RL     FGMRESeps        _RL     FGMRESeps
87        _RL     JFNKtol        _RL     JFNKtol
88  C      C     backward differences extrapolation factors
89          _RL bdfFac, bdfAlpha
90    C
91        _RL     recip_deltaT        _RL     recip_deltaT
92        LOGICAL JFNKconverged, krylovConverged        LOGICAL JFNKconverged, krylovConverged
93        LOGICAL writeNow        LOGICAL writeNow
94        CHARACTER*(MAX_LEN_MBUF) msgBuf        CHARACTER*(MAX_LEN_MBUF) msgBuf
95  C  
96  C     u/vIceRes :: residual of sea-ice momentum equations  C     u/vIceRes :: residual of sea-ice momentum equations
97        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
98        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
99  C     vector version of the residuals  C     extra time level required for backward difference time stepping
100        _RL resTmp (nVec,1,nSx,nSy)        _RL duIcNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
101          _RL dvIcNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
102  C     du/vIce   :: ice velocity increment to be added to u/vIce  C     du/vIce   :: ice velocity increment to be added to u/vIce
103        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
104        _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
105  C     precomputed (= constant per Newton iteration) versions of  C     precomputed (= constant per Newton iteration) versions of
106  C     zeta, eta, and DWATN, press  C     zeta, eta, and DWATN, press
107        _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL zetaPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
108          _RL zetaZPre(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
109        _RL etaPre  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL etaPre  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
110        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
111        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
112    C     work arrays
113          _RL rhs(nVec,nSx,nSy), sol(nVec,nSx,nSy)
114          _RL vv(nVec,im+1,nSx,nSy), w(nVec,im,nSx,nSy)
115          _RL wk1(nVec,nSx,nSy), wk2(nVec,nSx,nSy)
116  CEOP  CEOP
117    
118  C     Initialise  C     Initialise
# Line 117  C     with iOutFgmres=1, seaice_fgmres p Line 132  C     with iOutFgmres=1, seaice_fgmres p
132       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
133       &     iOutFGMRES=1       &     iOutFGMRES=1
134    
135  C      C     backward difference extrapolation factors
136          bdfFac = 0. _d 0
137          IF ( SEAICEuseBDF2 ) THEN
138           IF ( myIter.EQ.nIter0 .AND. SEAICEmomStartBDF.EQ.0 ) THEN
139            bdfFac = 0. _d 0
140           ELSE
141            bdfFac = 0.5 _d 0
142           ENDIF
143          ENDIF
144          bdfAlpha = 1. _d 0 + bdfFac
145    
146        DO bj=myByLo(myThid),myByHi(myThid)        DO bj=myByLo(myThid),myByHi(myThid)
147         DO bi=myBxLo(myThid),myBxHi(myThid)         DO bi=myBxLo(myThid),myBxHi(myThid)
148          DO J=1-Oly,sNy+Oly          DO J=1-OLy,sNy+OLy
149           DO I=1-Olx,sNx+Olx           DO I=1-OLx,sNx+OLx
150            uIceRes(I,J,bi,bj) = 0. _d 0            uIceRes(I,J,bi,bj) = 0. _d 0
151            vIceRes(I,J,bi,bj) = 0. _d 0            vIceRes(I,J,bi,bj) = 0. _d 0
152            duIce  (I,J,bi,bj) = 0. _d 0            duIce  (I,J,bi,bj) = 0. _d 0
153            dvIce  (I,J,bi,bj) = 0. _d 0            dvIce  (I,J,bi,bj) = 0. _d 0
154             ENDDO
155            ENDDO
156    C     cycle ice velocities
157            DO J=1-OLy,sNy+OLy
158             DO I=1-OLx,sNx+OLx
159              duIcNm1(I,J,bi,bj) = uIce(I,J,bi,bj) * bdfAlpha
160         &         + ( uIce(I,J,bi,bj) - uIceNm1(I,J,bi,bj) ) * bdfFac
161              dvIcNm1(I,J,bi,bj) = vIce(I,J,bi,bj) * bdfAlpha
162         &         + ( vIce(I,J,bi,bj) - vIceNm1(I,J,bi,bj) ) * bdfFac
163            uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)            uIceNm1(I,J,bi,bj) = uIce(I,J,bi,bj)
164            vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)            vIceNm1(I,J,bi,bj) = vIce(I,J,bi,bj)
165           ENDDO           ENDDO
166          ENDDO          ENDDO
167    C     As long as IMEX is not properly implemented leave this commented out
168    CML        IF ( .NOT.SEAICEuseIMEX ) THEN
169  C     Compute things that do no change during the Newton iteration:  C     Compute things that do no change during the Newton iteration:
170  C     sea-surface tilt and wind stress:  C     sea-surface tilt and wind stress:
171  C     FORCEX/Y0 - mass*(u/vIceNm1)/deltaT  C     FORCEX/Y0 - mass*(1.5*u/vIceNm1+0.5*(u/vIceNm1-u/vIceNm2))/deltaT
172          DO J=1-Oly,sNy+Oly          DO J=1-OLy,sNy+OLy
173           DO I=1-Olx,sNx+Olx           DO I=1-OLx,sNx+OLx
174            FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)            FORCEX(I,J,bi,bj) = FORCEX0(I,J,bi,bj)
175       &         + seaiceMassU(I,J,bi,bj)*uIceNm1(I,J,bi,bj)*recip_deltaT           &         + seaiceMassU(I,J,bi,bj)*duIcNm1(I,J,bi,bj)*recip_deltaT
176            FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)            FORCEY(I,J,bi,bj) = FORCEY0(I,J,bi,bj)
177       &         + seaiceMassV(I,J,bi,bj)*vIceNm1(I,J,bi,bj)*recip_deltaT           &         + seaiceMassV(I,J,bi,bj)*dvIcNm1(I,J,bi,bj)*recip_deltaT
178           ENDDO           ENDDO
179          ENDDO          ENDDO
180    CML        ENDIF
181         ENDDO         ENDDO
182        ENDDO        ENDDO
183  C     Start nonlinear Newton iteration: outer loop iteration  C     Start nonlinear Newton iteration: outer loop iteration
184        DO WHILE ( newtonIter.LT.SEAICEnewtonIterMax .AND.        DO WHILE ( newtonIter.LT.SEAICEnonLinIterMax .AND.
185       &     .NOT.JFNKconverged )       &     .NOT.JFNKconverged )
186         newtonIter = newtonIter + 1         newtonIter = newtonIter + 1
187  C     Compute initial residual F(u), (includes computation of global  C     Compute initial residual F(u), (includes computation of global
188  C     variables DWATN, zeta, and eta)  C     variables DWATN, zeta, and eta)
189         IF ( newtonIter .EQ. 1 ) CALL SEAICE_JFNK_UPDATE(         IF ( newtonIter .EQ. 1 ) CALL SEAICE_JFNK_UPDATE(
190       I      duIce, dvIce,       I      duIce, dvIce,
191       U      uIce, vIce, JFNKresidual,       U      uIce, vIce, JFNKresidual,
192       O      uIceRes, vIceRes,       O      uIceRes, vIceRes,
193       I      newtonIter, myTime, myIter, myThid )       I      newtonIter, myTime, myIter, myThid )
# Line 158  C     local copies of precomputed coeffi Line 195  C     local copies of precomputed coeffi
195  C     constant for the preconditioner  C     constant for the preconditioner
196         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
197          DO bi=myBxLo(myThid),myBxHi(myThid)          DO bi=myBxLo(myThid),myBxHi(myThid)
198           DO j=1-Oly,sNy+Oly           DO j=1-OLy,sNy+OLy
199            DO i=1-Olx,sNx+Olx            DO i=1-OLx,sNx+OLx
200             zetaPre(I,J,bi,bj) =  zeta(I,J,bi,bj)             zetaPre(I,J,bi,bj) =  zeta(I,J,bi,bj)
201               zetaZPre(I,J,bi,bj)= zetaZ(I,J,bi,bj)
202              etaPre(I,J,bi,bj) =   eta(I,J,bi,bj)              etaPre(I,J,bi,bj) =   eta(I,J,bi,bj)
203             etaZPre(I,J,bi,bj) =  etaZ(I,J,bi,bj)             etaZPre(I,J,bi,bj) =  etaZ(I,J,bi,bj)
204             dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)             dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
# Line 170  C     constant for the preconditioner Line 208  C     constant for the preconditioner
208         ENDDO         ENDDO
209  C     compute convergence criterion for linear preconditioned FGMRES  C     compute convergence criterion for linear preconditioned FGMRES
210         JFNKgamma_lin = JFNKgamma_lin_max         JFNKgamma_lin = JFNKgamma_lin_max
211         IF ( newtonIter.GT.1.AND.newtonIter.LE.100         IF ( newtonIter.GT.1.AND.newtonIter.LE.SEAICE_JFNK_tolIter
212       &      .AND.JFNKresidual.LT.JFNKres_t ) THEN       &      .AND.JFNKresidual.LT.JFNKres_t ) THEN
213  C     Eisenstat, 1996, equ.(2.6)        C     Eisenstat and Walker (1996), eq.(2.6)
214          phi_e = 1. _d 0          JFNKgamma_lin = SEAICE_JFNKphi
215          alp_e = 1. _d 0       &       *( JFNKresidual/JFNKresidualKm1 )**SEAICE_JFNKalpha
         JFNKgamma_lin = phi_e*( JFNKresidual/JFNKresidualKm1 )**alp_e  
216          JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)          JFNKgamma_lin = min(JFNKgamma_lin_max, JFNKgamma_lin)
217          JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)          JFNKgamma_lin = max(JFNKgamma_lin_min, JFNKgamma_lin)
218         ENDIF         ENDIF
219  C     save the residual for the next iteration  C     save the residual for the next iteration
220         JFNKresidualKm1 = JFNKresidual         JFNKresidualKm1 = JFNKresidual
221  C  
222  C     The Krylov iteration using FGMRES, the preconditioner is LSOR  C     The Krylov iteration using FGMRES, the preconditioner is LSOR
223  C     for now. The code is adapted from SEAICE_LSR, but heavily stripped  C     for now. The code is adapted from SEAICE_LSR, but heavily stripped
224  C     down.  C     down.
# Line 189  C     krylovIter is mapped into "its" in Line 226  C     krylovIter is mapped into "its" in
226  C     in that routine  C     in that routine
227         krylovIter    = 0         krylovIter    = 0
228         iCode         = 0         iCode         = 0
229  C  
230         JFNKconverged = JFNKresidual.LT.JFNKtol         JFNKconverged = JFNKresidual.LT.JFNKtol
231  C  
232  C     do Krylov loop only if convergence is not reached  C     do Krylov loop only if convergence is not reached
233  C  
234         IF ( .NOT.JFNKconverged ) THEN         IF ( .NOT.JFNKconverged ) THEN
235  C  
236  C     start Krylov iteration (FGMRES)  C     start Krylov iteration (FGMRES)
237  C  
238          krylovConverged = .FALSE.          krylovConverged = .FALSE.
239          FGMRESeps = JFNKgamma_lin * JFNKresidual          FGMRESeps = JFNKgamma_lin * JFNKresidual
240          DO WHILE ( .NOT.krylovConverged )  C     map first guess sol; it is zero because the solution is a correction
241           CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,sol,.TRUE.,myThid)
242    C     map rhs and change its sign because we are solving J*u = -F
243            CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,rhs,.TRUE.,myThid)
244            DO bj=myByLo(myThid),myByHi(myThid)
245             DO bi=myBxLo(myThid),myBxHi(myThid)
246              DO j=1,nVec
247               rhs(j,bi,bj) = - rhs(j,bi,bj)
248              ENDDO
249             ENDDO
250            ENDDO
251            DO WHILE ( .NOT.krylovConverged )
252  C     solution vector sol = du/vIce  C     solution vector sol = du/vIce
253  C     residual vector (rhs) Fu = u/vIceRes  C     residual vector (rhs) Fu = u/vIceRes
254  C     output work vectors wk1, -> input work vector wk2  C     output work vectors wk1, -> input work vector wk2
255    
256    C     map preconditioner results or Jacobian times vector,
257    C     stored in du/vIce to wk2, for iCode=0, wk2 is set to zero,
258    C     because du/vIce = 0
259             CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,wk2,.TRUE.,myThid)
260    C
261             CALL SEAICE_FGMRES (nVec,im,rhs,sol,ifgmres,krylovIter,
262         U        vv,w,wk1,wk2,
263         I        FGMRESeps,SEAICElinearIterMax,iOutFGMRES,
264         U        iCode,
265         I        myThid)
266  C      C    
267           CALL SEAICE_FGMRES_DRIVER(           IF ( iCode .EQ. 0 ) THEN
268       I        uIceRes, vIceRes,  C     map sol(ution) vector to du/vIce
269       U        duIce, dvIce, iCode,            CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,sol,.FALSE.,myThid)
270       I        FGMRESeps, iOutFGMRES,           ELSE
271       I        newtonIter, krylovIter, myTime, myIter, myThid )  C     map work vector to du/vIce to either compute a preconditioner
272    C     solution (wk1=rhs) or a Jacobian times wk1
273              CALL SEAICE_MAP2VEC(nVec,duIce,dvIce,wk1,.FALSE.,myThid)
274             ENDIF
275    C     Fill overlaps in updated fields
276             CALL EXCH_UV_XY_RL( duIce, dvIce,.TRUE.,myThid)
277  C     FGMRES returns iCode either asking for an new preconditioned vector  C     FGMRES returns iCode either asking for an new preconditioned vector
278  C     or product of matrix (Jacobian) times vector. For iCode = 0, terminate  C     or product of matrix (Jacobian) times vector. For iCode = 0, terminate
279  C     iteration  C     iteration
280           IF (iCode.EQ.1) THEN           IF (iCode.EQ.1) THEN
281  C     Call preconditioner  C     Call preconditioner
282            IF ( SOLV_MAX_ITERS .GT. 0 )            IF ( SEAICEpreconLinIter .GT. 0 )
283       &         CALL SEAICE_PRECONDITIONER(       &         CALL SEAICE_PRECONDITIONER(
284       U         duIce, dvIce,       U         duIce, dvIce,
285       I         zetaPre, etaPre, etaZpre, dwatPre,       I         zetaPre, etaPre, etaZpre, zetaZpre, dwatPre,
286       I         newtonIter, krylovIter, myTime, myIter, myThid )       I         newtonIter, krylovIter, myTime, myIter, myThid )
287           ELSEIF (iCode.GE.2) THEN           ELSEIF (iCode.GE.2) THEN
288  C     Compute Jacobian times vector  C     Compute Jacobian times vector
289            CALL SEAICE_JACVEC(            CALL SEAICE_JACVEC(
290       I         uIce, vIce, uIceRes, vIceRes,       I         uIce, vIce, uIceRes, vIceRes,
291       U         duIce, dvIce,         U         duIce, dvIce,
292       I         newtonIter, krylovIter, myTime, myIter, myThid )       I         newtonIter, krylovIter, myTime, myIter, myThid )
293           ENDIF           ENDIF
294           krylovConverged = iCode.EQ.0           krylovConverged = iCode.EQ.0
# Line 234  C     End of Krylov iterate Line 298  C     End of Krylov iterate
298  C     some output diagnostics  C     some output diagnostics
299          IF ( debugLevel.GE.debLevA ) THEN          IF ( debugLevel.GE.debLevA ) THEN
300           _BEGIN_MASTER( myThid )           _BEGIN_MASTER( myThid )
301           totalNewtonItersLoc =           totalNewtonItersLoc =
302       &        SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter       &        SEAICEnonLinIterMax*(myIter-nIter0)+newtonIter
303           WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')           WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
304       &        ' S/R SEAICE_JFNK: Newton iterate / total, ',       &        ' S/R SEAICE_JFNK: Newton iterate / total, ',
305       &        'JFNKgamma_lin, initial norm = ',       &        'JFNKgamma_lin, initial norm = ',
306       &        newtonIter, totalNewtonItersLoc,       &        newtonIter, totalNewtonItersLoc,
# Line 244  C     some output diagnostics Line 308  C     some output diagnostics
308           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
309       &        SQUEEZE_RIGHT, myThid )       &        SQUEEZE_RIGHT, myThid )
310           WRITE(msgBuf,'(3(A,I6))')           WRITE(msgBuf,'(3(A,I6))')
311       &        ' S/R SEAICE_JFNK: Newton iterate / total = ',newtonIter,       &        ' S/R SEAICE_JFNK: Newton iterate / total = ',newtonIter,
312       &        ' / ', totalNewtonItersLoc,       &        ' / ', totalNewtonItersLoc,
313       &        ', Nb. of FGMRES iterations = ', krylovIter       &        ', Nb. of FGMRES iterations = ', krylovIter
314           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
315       &        SQUEEZE_RIGHT, myThid )       &        SQUEEZE_RIGHT, myThid )
316           _END_MASTER( myThid )           _END_MASTER( myThid )
317          ENDIF          ENDIF
318          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN          IF ( krylovIter.EQ.SEAICElinearIterMax ) THEN
319           krylovFails = krylovFails + 1           krylovFails = krylovFails + 1
320          ENDIF          ENDIF
321  C     Set the stopping criterion for the Newton iteration  C     Set the stopping criterion for the Newton iteration and the
322          IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual  C     criterion for the transition from accurate to approximate FGMRES
323            IF ( newtonIter .EQ. 1 ) THEN
324             JFNKtol=SEAICEnonLinTol*JFNKresidual
325             IF ( JFNKres_tFac .NE. UNSET_RL )
326         &        JFNKres_t = JFNKresidual * JFNKres_tFac
327            ENDIF
328  C     Update linear solution vector and return to Newton iteration  C     Update linear solution vector and return to Newton iteration
329  C     Do a linesearch if necessary, and compute a new residual.  C     Do a linesearch if necessary, and compute a new residual.
330  C     Note that it should be possible to do the following operations  C     Note that it should be possible to do the following operations
331  C     at the beginning of the Newton iteration, thereby saving us from  C     at the beginning of the Newton iteration, thereby saving us from
332  C     the extra call of seaice_jfnk_update, but unfortunately that  C     the extra call of seaice_jfnk_update, but unfortunately that
333  C     changes the results, so we leave the stuff here for now.  C     changes the results, so we leave the stuff here for now.
334          CALL SEAICE_JFNK_UPDATE(          CALL SEAICE_JFNK_UPDATE(
335       I       duIce, dvIce,       I       duIce, dvIce,
336       U       uIce, vIce, JFNKresidual,       U       uIce, vIce, JFNKresidual,
337       O       uIceRes, vIceRes,       O       uIceRes, vIceRes,
338       I       newtonIter, myTime, myIter, myThid )       I       newtonIter, myTime, myIter, myThid )
339  C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver  C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
340          DO bj=myByLo(myThid),myByHi(myThid)          DO bj=myByLo(myThid),myByHi(myThid)
341           DO bi=myBxLo(myThid),myBxHi(myThid)           DO bi=myBxLo(myThid),myBxHi(myThid)
342            DO J=1-Oly,sNy+Oly            DO J=1-OLy,sNy+OLy
343             DO I=1-Olx,sNx+Olx             DO I=1-OLx,sNx+OLx
344              duIce(I,J,bi,bj)= 0. _d 0              duIce(I,J,bi,bj)= 0. _d 0
345              dvIce(I,J,bi,bj)= 0. _d 0              dvIce(I,J,bi,bj)= 0. _d 0
346             ENDDO             ENDDO
# Line 281  C     reset du/vIce here instead of sett Line 350  C     reset du/vIce here instead of sett
350         ENDIF         ENDIF
351  C     end of Newton iterate  C     end of Newton iterate
352        ENDDO        ENDDO
353  C  
354  C--   Output diagnostics  C--   Output diagnostics
355  C  
356        IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN        IF ( SEAICE_monFreq .GT. 0. _d 0 ) THEN
357  C     Count iterations  C     Count iterations
358         totalJFNKtimeSteps = totalJFNKtimeSteps + 1         totalJFNKtimeSteps = totalJFNKtimeSteps + 1
# Line 291  C     Count iterations Line 360  C     Count iterations
360         totalKrylovIters   = totalKrylovIters + totalKrylovItersLoc         totalKrylovIters   = totalKrylovIters + totalKrylovItersLoc
361  C     Record failure  C     Record failure
362         totalKrylovFails   = totalKrylovFails + krylovFails         totalKrylovFails   = totalKrylovFails + krylovFails
363         IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN         IF ( newtonIter .EQ. SEAICEnonLinIterMax ) THEN
364          totalNewtonFails = totalNewtonFails + 1          totalNewtonFails = totalNewtonFails + 1
365         ENDIF         ENDIF
366        ENDIF        ENDIF
367  C     Decide whether it is time to dump and reset the counter  C     Decide whether it is time to dump and reset the counter
368        writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,        writeNow = DIFFERENT_MULTIPLE(SEAICE_monFreq,
369       &     myTime+deltaTClock, deltaTClock)       &     myTime+deltaTClock, deltaTClock)
370  #ifdef ALLOW_CAL  #ifdef ALLOW_CAL
371        IF ( useCAL ) THEN        IF ( useCAL ) THEN
372         CALL CAL_TIME2DUMP(         CALL CAL_TIME2DUMP(
373       I      zeroRL, SEAICE_monFreq,  deltaTClock,       I      zeroRL, SEAICE_monFreq,  deltaTClock,
374       U      writeNow,       U      writeNow,
375       I      myTime+deltaTclock, myIter+1, myThid )       I      myTime+deltaTclock, myIter+1, myThid )
# Line 308  C     Decide whether it is time to dump Line 377  C     Decide whether it is time to dump
377  #endif  #endif
378        IF ( writeNow ) THEN        IF ( writeNow ) THEN
379         _BEGIN_MASTER( myThid )         _BEGIN_MASTER( myThid )
380         WRITE(msgBuf,'(A)')         WRITE(msgBuf,'(A)')
381       &' // ======================================================='       &' // ======================================================='
382         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
383       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
384         WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'         WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'
385         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
386       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
387         WRITE(msgBuf,'(A)')         WRITE(msgBuf,'(A)')
388       &' // ======================================================='       &' // ======================================================='
389         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
390       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
391         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
392       &      ' %JFNK_MON: time step              = ', myIter+1       &      ' %JFNK_MON: time step              = ', myIter+1
393         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
394       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
395         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
396       &      ' %JFNK_MON: Nb. of time steps      = ', totalJFNKtimeSteps       &      ' %JFNK_MON: Nb. of time steps      = ', totalJFNKtimeSteps
397         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
398       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
399         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
400       &      ' %JFNK_MON: Nb. of Newton steps    = ', totalNewtonIters       &      ' %JFNK_MON: Nb. of Newton steps    = ', totalNewtonIters
401         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
402       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
403         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
404       &      ' %JFNK_MON: Nb. of Krylov steps    = ', totalKrylovIters       &      ' %JFNK_MON: Nb. of Krylov steps    = ', totalKrylovIters
405         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
406       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
407         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
408       &      ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails       &      ' %JFNK_MON: Nb. of Newton failures = ', totalNewtonFails
409         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
410       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
411         WRITE(msgBuf,'(A,I10)')         WRITE(msgBuf,'(A,I10)')
412       &      ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails       &      ' %JFNK_MON: Nb. of Krylov failures = ', totalKrylovFails
413         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
414       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
415         WRITE(msgBuf,'(A)')         WRITE(msgBuf,'(A)')
416       &' // ======================================================='       &' // ======================================================='
417         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
418       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
419         WRITE(msgBuf,'(A)') ' // End JFNK statistics'         WRITE(msgBuf,'(A)') ' // End JFNK statistics'
420         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
421       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
422         WRITE(msgBuf,'(A)')         WRITE(msgBuf,'(A)')
423       &' // ======================================================='       &' // ======================================================='
424         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
425       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
# Line 365  C     reset and start again Line 434  C     reset and start again
434    
435  C     Print more debugging information  C     Print more debugging information
436        IF ( debugLevel.GE.debLevA ) THEN        IF ( debugLevel.GE.debLevA ) THEN
437         IF ( newtonIter .EQ. SEAICEnewtonIterMax ) THEN         IF ( newtonIter .EQ. SEAICEnonLinIterMax ) THEN
438          _BEGIN_MASTER( myThid )          _BEGIN_MASTER( myThid )
439          WRITE(msgBuf,'(A,I10)')          WRITE(msgBuf,'(A,I10)')
440       &       ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',       &       ' S/R SEAICE_JFNK: JFNK did not converge in timestep ',
441       &       myIter+1       &       myIter+1
442          CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,          CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
# Line 376  C     Print more debugging information Line 445  C     Print more debugging information
445         ENDIF         ENDIF
446         IF ( krylovFails .GT. 0 ) THEN         IF ( krylovFails .GT. 0 ) THEN
447          _BEGIN_MASTER( myThid )          _BEGIN_MASTER( myThid )
448          WRITE(msgBuf,'(A,I4,A,I10)')          WRITE(msgBuf,'(A,I4,A,I10)')
449       &       ' S/R SEAICE_JFNK: FGMRES did not converge ',       &       ' S/R SEAICE_JFNK: FGMRES did not converge ',
450       &       krylovFails, ' times in timestep ', myIter+1       &       krylovFails, ' times in timestep ', myIter+1
451          CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,          CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
# Line 384  C     Print more debugging information Line 453  C     Print more debugging information
453          _END_MASTER( myThid )          _END_MASTER( myThid )
454         ENDIF         ENDIF
455         _BEGIN_MASTER( myThid )         _BEGIN_MASTER( myThid )
456         WRITE(msgBuf,'(A,I6,A,I10)')         WRITE(msgBuf,'(A,I6,A,I10)')
457       &      ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',       &      ' S/R SEAICE_JFNK: Total number FGMRES iterations = ',
458       &      totalKrylovItersLoc, ' in timestep ', myIter+1       &      totalKrylovItersLoc, ' in timestep ', myIter+1
459         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
# Line 400  CBOP Line 469  CBOP
469  C     !ROUTINE: SEAICE_JFNK_UPDATE  C     !ROUTINE: SEAICE_JFNK_UPDATE
470  C     !INTERFACE:  C     !INTERFACE:
471    
472        SUBROUTINE SEAICE_JFNK_UPDATE(        SUBROUTINE SEAICE_JFNK_UPDATE(
473       I     duIce, dvIce,       I     duIce, dvIce,
474       U     uIce, vIce, JFNKresidual,       U     uIce, vIce, JFNKresidual,
475       O     uIceRes, vIceRes,       O     uIceRes, vIceRes,
476       I     newtonIter, myTime, myIter, myThid )       I     newtonIter, myTime, myIter, myThid )
# Line 460  C     i,j,bi,bj :: loop indices Line 529  C     i,j,bi,bj :: loop indices
529        _RL     resLoc, facLS        _RL     resLoc, facLS
530        LOGICAL doLineSearch        LOGICAL doLineSearch
531  C     nVec    :: size of the input vector(s)  C     nVec    :: size of the input vector(s)
532  C     vector version of the residuals  C     resTmp  :: vector version of the residuals
533        INTEGER nVec        INTEGER nVec
534        PARAMETER ( nVec  = 2*sNx*sNy )        PARAMETER ( nVec  = 2*sNx*sNy )
535        _RL resTmp (nVec,1,nSx,nSy)        _RL resTmp (nVec,1,nSx,nSy)
536  C  
537        CHARACTER*(MAX_LEN_MBUF) msgBuf        CHARACTER*(MAX_LEN_MBUF) msgBuf
538  CEOP  CEOP
539    
# Line 474  C     Initialise some local variables Line 543  C     Initialise some local variables
543        facLS = 1. _d 0        facLS = 1. _d 0
544        doLineSearch = .TRUE.        doLineSearch = .TRUE.
545        DO WHILE ( doLineSearch )        DO WHILE ( doLineSearch )
 C     Determine, if we need more iterations  
        doLineSearch = resLoc .GE. JFNKresidual  
 C     Limit the maximum number of iterations arbitrarily to four  
        doLineSearch = doLineSearch .AND. l .LE. 4  
 C     For the first iteration du/vIce = 0 and there will be no  
 C     improvement of the residual possible, so we do only the first  
 C     iteration  
        IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.  
 C     Only start a linesearch after some Newton iterations  
        IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.  
 C     Increment counter  
        l = l + 1  
546  C     Create update  C     Create update
547         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
548          DO bi=myBxLo(myThid),myBxHi(myThid)          DO bi=myBxLo(myThid),myBxHi(myThid)
549           DO J=1-Oly,sNy+Oly           DO J=1-OLy,sNy+OLy
550            DO I=1-Olx,sNx+Olx            DO I=1-OLx,sNx+OLx
551             uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)             uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)
552             vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)             vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)
553            ENDDO            ENDDO
# Line 499  C     Create update Line 556  C     Create update
556         ENDDO         ENDDO
557  C     Compute current residual F(u), (includes re-computation of global  C     Compute current residual F(u), (includes re-computation of global
558  C     variables DWATN, zeta, and eta, i.e. they are different after this)  C     variables DWATN, zeta, and eta, i.e. they are different after this)
559         CALL SEAICE_CALC_RESIDUAL(         CALL SEAICE_CALC_RESIDUAL(
560       I      uIce, vIce,       I      uIce, vIce,
561       O      uIceRes, vIceRes,       O      uIceRes, vIceRes,
562       I      newtonIter, 0, myTime, myIter, myThid )       I      newtonIter, 0, myTime, myIter, myThid )
563  C     Important: Compute the norm of the residual using the same scalar  C     Important: Compute the norm of the residual using the same scalar
564  C     product that SEAICE_FGMRES does  C     product that SEAICE_FGMRES does
565         CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)         CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
566         CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)         CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)
567         resLoc = SQRT(resLoc)         resLoc = SQRT(resLoc)
568    C     Determine, if we need more iterations
569           doLineSearch = resLoc .GE. JFNKresidual
570    C     Limit the maximum number of iterations arbitrarily to four
571           doLineSearch = doLineSearch .AND. l .LT. 4
572    C     For the first iteration du/vIce = 0 and there will be no
573    C     improvement of the residual possible, so we do only the first
574    C     iteration
575           IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.
576    C     Only start a linesearch after some Newton iterations
577           IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.
578    C     Increment counter
579           l = l + 1
580  C     some output diagnostics  C     some output diagnostics
581         IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN         IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN
582          _BEGIN_MASTER( myThid )          _BEGIN_MASTER( myThid )
583          WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')          WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')
584       &       ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',       &       ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',
585       &       'facLS, JFNKresidual, resLoc = ',       &       'facLS, JFNKresidual, resLoc = ',
586       &        newtonIter, l, facLS, JFNKresidual, resLoc       &        newtonIter, l, facLS, JFNKresidual, resLoc
# Line 527  C     iterations, 0.25*du/vIce in the se Line 596  C     iterations, 0.25*du/vIce in the se
596  C     This is the new residual  C     This is the new residual
597        JFNKresidual = resLoc        JFNKresidual = resLoc
598    
599  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */  #endif /* SEAICE_ALLOW_JFNK */
600    
601        RETURN        RETURN
602        END        END

Legend:
Removed from v.1.16  
changed lines
  Added in v.1.30

  ViewVC Help
Powered by ViewVC 1.1.22