/[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.9 by mlosch, Mon Nov 12 09:46:38 2012 UTC revision 1.15 by mlosch, Wed Jan 16 21:20:28 2013 UTC
# Line 3  C $Name$ Line 3  C $Name$
3    
4  #include "SEAICE_OPTIONS.h"  #include "SEAICE_OPTIONS.h"
5    
6    C--  File seaice_jfnk.F: seaice jfnk dynamical solver S/R:
7    C--   Contents
8    C--   o SEAICE_JFNK
9    C--   o SEAICE_JFNK_UPDATE
10    
11  CBOP  CBOP
12  C     !ROUTINE: SEAICE_JFNK  C     !ROUTINE: SEAICE_JFNK
13  C     !INTERFACE:  C     !INTERFACE:
# Line 10  C     !INTERFACE: Line 15  C     !INTERFACE:
15    
16  C     !DESCRIPTION: \bv  C     !DESCRIPTION: \bv
17  C     *==========================================================*  C     *==========================================================*
18  C     | SUBROUTINE SEAICE_JFKF  C     | SUBROUTINE SEAICE_JFNK
19  C     | o Ice dynamics using a Jacobian-free Newton-Krylov solver  C     | o Ice dynamics using a Jacobian-free Newton-Krylov solver
20  C     |   following J.-F. Lemieux et al. Improving the numerical  C     |   following J.-F. Lemieux et al. Improving the numerical
21  C     |   convergence of viscous-plastic sea ice models with the  C     |   convergence of viscous-plastic sea ice models with the
# Line 60  C     i,j,bi,bj :: loop indices Line 65  C     i,j,bi,bj :: loop indices
65  C     loop indices  C     loop indices
66        INTEGER newtonIter        INTEGER newtonIter
67        INTEGER krylovIter, krylovFails        INTEGER krylovIter, krylovFails
68        INTEGER totalKrylovItersLoc        INTEGER totalKrylovItersLoc, totalNewtonItersLoc
69  C     FGMRES flag that determines amount of output messages of fgmres  C     FGMRES flag that determines amount of output messages of fgmres
70        INTEGER iOutFGMRES        INTEGER iOutFGMRES
71  C     FGMRES flag that indicates what fgmres wants us to do next  C     FGMRES flag that indicates what fgmres wants us to do next
72        INTEGER iCode        INTEGER iCode
73        _RL     JFNKresidual, JFNKresidualTile(nSx,nSy)        _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     phi_e, alp_e, JFNKgamma_lin
# Line 80  C Line 85  C
85  C     u/vIceRes :: residual of sea-ice momentum equations  C     u/vIceRes :: residual of sea-ice momentum equations
86        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
87        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
88    C     vector version of the residuals
89          _RL resTmp (nVec,1,nSx,nSy)
90  C     du/vIce   :: ice velocity increment to be added to u/vIce  C     du/vIce   :: ice velocity increment to be added to u/vIce
91        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
92        _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 89  C     zeta, eta, and DWATN, press Line 96  C     zeta, eta, and DWATN, press
96        _RL etaPre  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL etaPre  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
97        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL etaZPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
98        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)        _RL dwatPre (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
       _RL pressPre(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)  
99  CEOP  CEOP
100    
101  C     Initialise  C     Initialise
# Line 104  C     Initialise Line 110  C     Initialise
110        recip_deltaT        = 1. _d 0 / SEAICE_deltaTdyn        recip_deltaT        = 1. _d 0 / SEAICE_deltaTdyn
111    
112        iOutFGMRES=0        iOutFGMRES=0
113  C     iOutFgmres=1 gives a little bit of output  C     with iOutFgmres=1, seaice_fgmres prints the residual at each iteration
114        IF ( debugLevel.GE.debLevA .AND.        IF ( debugLevel.GE.debLevC .AND.
115       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )       &     DIFFERENT_MULTIPLE( SEAICE_monFreq, myTime, deltaTClock ) )
116       &     iOutFGMRES=1       &     iOutFGMRES=1
117    
# Line 141  C     Start nonlinear Newton iteration: Line 147  C     Start nonlinear Newton iteration:
147         newtonIter = newtonIter + 1         newtonIter = newtonIter + 1
148  C     Compute initial residual F(u), (includes computation of global  C     Compute initial residual F(u), (includes computation of global
149  C     variables DWATN, zeta, and eta)  C     variables DWATN, zeta, and eta)
150         CALL SEAICE_CALC_RESIDUAL(  C     Update linear solution vector and return to Newton iteration
151       I      uIce, vIce,  C     Do the linesearch
152       O      uIceRes, vIceRes,         CALL SEAICE_JFNK_UPDATE(
153       I      newtonIter, 0, myTime, myIter, myThid )       I      duIce, dvIce,
154         CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)       U      uIce, vIce, JFNKresidual,
155  C     local copies of precomputed coefficients that are to stay       O      uIceRes, vIceRes,
156  C     constant for the preconditioner       I      newtonIter, myTime, myIter, myThid )
157    C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
158         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
159          DO bi=myBxLo(myThid),myBxHi(myThid)          DO bi=myBxLo(myThid),myBxHi(myThid)
160           DO j=1-Oly,sNy+Oly           DO J=1-Oly,sNy+Oly
161            DO i=1-Olx,sNx+Olx            DO I=1-Olx,sNx+Olx
162              zetaPre(I,J,bi,bj) =  zeta(I,J,bi,bj)             duIce(I,J,bi,bj)= 0. _d 0
163               etaPre(I,J,bi,bj) =   eta(I,J,bi,bj)             dvIce(I,J,bi,bj)= 0. _d 0
             etaZPre(I,J,bi,bj) =  etaZ(I,J,bi,bj)  
             dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)  
            pressPre(I,J,bi,bj) = press(I,J,bi,bj)  
164            ENDDO            ENDDO
165           ENDDO           ENDDO
166          ENDDO          ENDDO
167         ENDDO         ENDDO
168  C      CMLC     Do it again, Sam
169    CML       CALL SEAICE_CALC_RESIDUAL(
170    CML     I      uIce, vIce,
171    CML     O      uIceRes, vIceRes,
172    CML     I      newtonIter, 0, myTime, myIter, myThid )
173    CMLC     probably not necessary, will be removed later:
174    CML       CALL EXCH_UV_XY_RL( uIceRes, vIceRes,.TRUE.,myThid)
175    CMLC     Important: Compute the norm of the residual using the same scalar
176    CMLC     product that SEAICE_FGMRES does
177    CML       CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
178    CML       CALL SEAICE_SCALPROD(
179    CML     &      nVec,1,1,1,resTmp,resTmp,JFNKresidual,myThid)
180    CML       JFNKresidual = SQRT(JFNKresidual)
181    C     local copies of precomputed coefficients that are to stay
182    C     constant for the preconditioner
183         DO bj=myByLo(myThid),myByHi(myThid)         DO bj=myByLo(myThid),myByHi(myThid)
184          DO bi=myBxLo(myThid),myBxHi(myThid)          DO bi=myBxLo(myThid),myBxHi(myThid)
185           JFNKresidualTile(bi,bj) = 0. _d 0           DO j=1-Oly,sNy+Oly
186           DO J=1,sNy            DO i=1-Olx,sNx+Olx
187            DO I=1,sNx             zetaPre(I,J,bi,bj) =  zeta(I,J,bi,bj)
188  #ifdef CG2D_SINGLECPU_SUM              etaPre(I,J,bi,bj) =   eta(I,J,bi,bj)
189             JFNKlocalBuf(I,J,bi,bj) =             etaZPre(I,J,bi,bj) =  etaZ(I,J,bi,bj)
190  #else             dwatPre(I,J,bi,bj) = DWATN(I,J,bi,bj)
            JFNKresidualTile(bi,bj) = JFNKresidualTile(bi,bj) +  
 #endif  
      &          uIceRes(I,J,bi,bj)*uIceRes(I,J,bi,bj) +  
      &          vIceRes(I,J,bi,bj)*vIceRes(I,J,bi,bj)  
191            ENDDO            ENDDO
192           ENDDO           ENDDO
193          ENDDO          ENDDO
194         ENDDO         ENDDO
        JFNKresidual = 0. _d 0  
 #ifdef CG2D_SINGLECPU_SUM  
        CALL GLOBAL_SUM_SINGLECPU_RL(  
      &         JFNKlocalBuf,JFNKresidual, 0, 0, myThid)  
 #else  
        CALL GLOBAL_SUM_TILE_RL( JFNKresidualTile,JFNKresidual,myThid )  
 #endif  
        JFNKresidual = SQRT(JFNKresidual)  
195  C     compute convergence criterion for linear preconditioned FGMRES  C     compute convergence criterion for linear preconditioned FGMRES
196         JFNKgamma_lin = JFNKgamma_lin_max         JFNKgamma_lin = JFNKgamma_lin_max
197         IF ( newtonIter.GT.1.AND.newtonIter.LE.100         IF ( newtonIter.GT.1.AND.newtonIter.LE.100
# Line 207  C     krylovIter is mapped into "its" in Line 213  C     krylovIter is mapped into "its" in
213  C     in that routine  C     in that routine
214         krylovIter    = 0         krylovIter    = 0
215         iCode         = 0         iCode         = 0
        IF ( debugLevel.GE.debLevA ) THEN    
         _BEGIN_MASTER( myThid )  
         WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')  
      &       ' S/R SEAICE_JFNK: newtonIter,',  
      &       ' total newtonIter, JFNKgamma_lin, initial norm = ',  
      &       newtonIter,SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,  
      &       JFNKgamma_lin, JFNKresidual  
         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,  
      &       SQUEEZE_RIGHT, myThid )  
         _END_MASTER( myThid )  
        ENDIF  
216  C  C
217         JFNKconverged = JFNKresidual.LT.JFNKtol         JFNKconverged = JFNKresidual.LT.JFNKtol
218  C  C
# Line 247  C     Call preconditioner Line 242  C     Call preconditioner
242            IF ( SOLV_MAX_ITERS .GT. 0 )            IF ( SOLV_MAX_ITERS .GT. 0 )
243       &         CALL SEAICE_PRECONDITIONER(       &         CALL SEAICE_PRECONDITIONER(
244       U         duIce, dvIce,       U         duIce, dvIce,
245       I         zetaPre, etaPre, etaZpre, dwatPre, pressPre,       I         zetaPre, etaPre, etaZpre, dwatPre,
246       I         newtonIter, krylovIter, myTime, myIter, myThid )       I         newtonIter, krylovIter, myTime, myIter, myThid )
247           ELSEIF (iCode.GE.2) THEN           ELSEIF (iCode.GE.2) THEN
248  C     Compute Jacobian times vector  C     Compute Jacobian times vector
# Line 263  C     End of Krylov iterate Line 258  C     End of Krylov iterate
258  C     some output diagnostics  C     some output diagnostics
259          IF ( debugLevel.GE.debLevA ) THEN          IF ( debugLevel.GE.debLevA ) THEN
260           _BEGIN_MASTER( myThid )           _BEGIN_MASTER( myThid )
261             totalNewtonItersLoc =
262         &        SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter
263             WRITE(msgBuf,'(2A,2(1XI6),2E12.5)')
264         &        ' S/R SEAICE_JFNK: Newton iterate / total, ',
265         &        'JFNKgamma_lin, initial norm = ',
266         &        newtonIter, totalNewtonItersLoc,
267         &        JFNKgamma_lin,JFNKresidual
268             CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
269         &        SQUEEZE_RIGHT, myThid )
270           WRITE(msgBuf,'(3(A,I6))')           WRITE(msgBuf,'(3(A,I6))')
271       &        ' S/R SEAICE_JFNK: Newton iterate / total = ', newtonIter,       &        ' S/R SEAICE_JFNK: Newton iterate / total = ',newtonIter,
272       &        ' / ', SEAICEnewtonIterMax*(myIter-nIter0)+newtonIter,       &        ' / ', totalNewtonItersLoc,
273       &        ', Nb. of FGMRES iterations = ', krylovIter       &        ', Nb. of FGMRES iterations = ', krylovIter
274           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,           CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
275       &        SQUEEZE_RIGHT, myThid )       &        SQUEEZE_RIGHT, myThid )
# Line 274  C     some output diagnostics Line 278  C     some output diagnostics
278          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN          IF ( krylovIter.EQ.SEAICEkrylovIterMax ) THEN
279           krylovFails = krylovFails + 1           krylovFails = krylovFails + 1
280          ENDIF          ENDIF
281    C     Set the stopping criterion for the Newton iteration
282            IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual
283  C     Update linear solution vector and return to Newton iteration  C     Update linear solution vector and return to Newton iteration
284    C     Do a linesearch if necessary, and compute a new residual.
285    C     Note that it should be possible to do the following operations
286    C     at the beginning of the Newton iteration, thereby saving us from
287    C     the extra call of seaice_jfnk_update, but unfortunately that
288    C     changes the results, so we leave the stuff here for now.
289            CALL SEAICE_JFNK_UPDATE(
290         I       duIce, dvIce,
291         U       uIce, vIce, JFNKresidual,
292         O       uIceRes, vIceRes,
293         I       newtonIter, myTime, myIter, myThid )
294    C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver
295          DO bj=myByLo(myThid),myByHi(myThid)          DO bj=myByLo(myThid),myByHi(myThid)
296           DO bi=myBxLo(myThid),myBxHi(myThid)           DO bi=myBxLo(myThid),myBxHi(myThid)
297            DO J=1-Oly,sNy+Oly            DO J=1-Oly,sNy+Oly
298             DO I=1-Olx,sNx+Olx             DO I=1-Olx,sNx+Olx
             uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+duIce(I,J,bi,bj)  
             vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+dvIce(I,J,bi,bj)  
 C     reset du/vIce here instead of setting sol = 0 in seaice_fgmres_driver  
299              duIce(I,J,bi,bj)= 0. _d 0              duIce(I,J,bi,bj)= 0. _d 0
300              dvIce(I,J,bi,bj)= 0. _d 0              dvIce(I,J,bi,bj)= 0. _d 0
301             ENDDO             ENDDO
302            ENDDO            ENDDO
303           ENDDO           ENDDO
304          ENDDO          ENDDO
 C     Set the stopping criterion for the Newton iteration  
         IF ( newtonIter .EQ. 1 ) JFNKtol=JFNKgamma_nonlin*JFNKresidual  
305         ENDIF         ENDIF
306  C     end of Newton iterate  C     end of Newton iterate
307        ENDDO        ENDDO
# Line 359  C     Decide whether it is time to dump Line 371  C     Decide whether it is time to dump
371       &' // ======================================================='       &' // ======================================================='
372         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
373       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
374         WRITE(msgBuf,'(A)') ' // Begin JFNK statistics'         WRITE(msgBuf,'(A)') ' // End JFNK statistics'
375         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,         CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
376       &      SQUEEZE_RIGHT, myThid )       &      SQUEEZE_RIGHT, myThid )
377         WRITE(msgBuf,'(A)')         WRITE(msgBuf,'(A)')
# Line 404  C     Print more debugging information Line 416  C     Print more debugging information
416         _END_MASTER( myThid )         _END_MASTER( myThid )
417        ENDIF        ENDIF
418    
419          RETURN
420          END
421    
422    CBOP
423    C     !ROUTINE: SEAICE_JFNK_UPDATE
424    C     !INTERFACE:
425    
426          SUBROUTINE SEAICE_JFNK_UPDATE(
427         I     duIce, dvIce,
428         U     uIce, vIce, JFNKresidual,
429         O     uIceRes, vIceRes,
430         I     newtonIter, myTime, myIter, myThid )
431    
432    C     !DESCRIPTION: \bv
433    C     *==========================================================*
434    C     | SUBROUTINE SEAICE_JFNK_UPDATE
435    C     | o Update velocities with incremental solutions of FGMRES
436    C     | o compute residual of updated solutions and do
437    C     | o linesearch:
438    C     |   reduce update until residual is smaller than previous
439    C     |   one (input)
440    C     *==========================================================*
441    C     | written by Martin Losch, Jan 2013
442    C     *==========================================================*
443    C     \ev
444    
445    C     !USES:
446          IMPLICIT NONE
447    
448    C     === Global variables ===
449    #include "SIZE.h"
450    #include "EEPARAMS.h"
451    #include "PARAMS.h"
452    #include "SEAICE_SIZE.h"
453    #include "SEAICE_PARAMS.h"
454    
455    C     !INPUT/OUTPUT PARAMETERS:
456    C     === Routine arguments ===
457    C     myTime :: Simulation time
458    C     myIter :: Simulation timestep number
459    C     myThid :: my Thread Id. number
460    C     newtonIter :: current iterate of Newton iteration
461          _RL     myTime
462          INTEGER myIter
463          INTEGER myThid
464          INTEGER newtonIter
465    C     JFNKresidual :: Residual at the beginning of the FGMRES iteration,
466    C                     changes with newtonIter (updated)
467          _RL     JFNKresidual
468    C     du/vIce   :: ice velocity increment to be added to u/vIce (input)
469          _RL duIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
470          _RL dvIce  (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
471    C     u/vIce    :: ice velocity increment to be added to u/vIce (updated)
472          _RL uIce   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
473          _RL vIce   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
474    C     u/vIceRes :: residual of sea-ice momentum equations (output)
475          _RL uIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
476          _RL vIceRes(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
477    
478    C     Local variables:
479    C     i,j,bi,bj :: loop indices
480          INTEGER i,j,bi,bj
481          INTEGER l
482          _RL     resLoc, facLS
483          LOGICAL doLineSearch
484    C     nVec    :: size of the input vector(s)
485    C     vector version of the residuals
486          INTEGER nVec
487          PARAMETER ( nVec  = 2*sNx*sNy )
488          _RL resTmp (nVec,1,nSx,nSy)
489    C
490          CHARACTER*(MAX_LEN_MBUF) msgBuf
491    CEOP
492    
493    C     Initialise some local variables
494          l = 0
495          resLoc = JFNKresidual
496          facLS = 1. _d 0
497          doLineSearch = .TRUE.
498          DO WHILE ( doLineSearch )
499    C     Determine, if we need more iterations
500           doLineSearch = resLoc .GE. JFNKresidual
501           doLineSearch = doLineSearch .AND. l .LE. 4
502    C     For the first iteration du/vIce = 0 and there will be no
503    C     improvement of the residual possible, so we do only the first
504    C     iteration
505           IF ( newtonIter .EQ. 1 ) doLineSearch = .FALSE.
506    C     Only start a linesearch after some Newton iterations
507           IF ( newtonIter .LE. SEAICE_JFNK_lsIter ) doLineSearch = .FALSE.
508    C     Increment counter
509           l = l + 1
510    C     Create update
511           DO bj=myByLo(myThid),myByHi(myThid)
512            DO bi=myBxLo(myThid),myBxHi(myThid)
513             DO J=1-Oly,sNy+Oly
514              DO I=1-Olx,sNx+Olx
515               uIce(I,J,bi,bj) = uIce(I,J,bi,bj)+facLS*duIce(I,J,bi,bj)
516               vIce(I,J,bi,bj) = vIce(I,J,bi,bj)+facLS*dvIce(I,J,bi,bj)
517              ENDDO
518             ENDDO
519            ENDDO
520           ENDDO
521    C     Compute current residual F(u), (includes re-computation of global
522    C     variables DWATN, zeta, and eta, i.e. they are different after this)
523           CALL SEAICE_CALC_RESIDUAL(
524         I      uIce, vIce,
525         O      uIceRes, vIceRes,
526         I      newtonIter, 0, myTime, myIter, myThid )
527    C     Important: Compute the norm of the residual using the same scalar
528    C     product that SEAICE_FGMRES does
529           CALL SEAICE_MAP2VEC(nVec,uIceRes,vIceRes,resTmp,.TRUE.,myThid)
530           CALL SEAICE_SCALPROD(nVec,1,1,1,resTmp,resTmp,resLoc,myThid)
531           resLoc = SQRT(resLoc)
532    C     some output diagnostics
533           IF ( debugLevel.GE.debLevA .AND. doLineSearch ) THEN
534            _BEGIN_MASTER( myThid )
535            WRITE(msgBuf,'(2A,2(1XI6),3E12.5)')
536         &       ' S/R SEAICE_JFNK_UPDATE: Newton iter, LSiter, ',
537         &       'facLS, JFNKresidual, resLoc = ',
538         &        newtonIter, l, facLS, JFNKresidual, resLoc
539            CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
540         &       SQUEEZE_RIGHT, myThid )
541            _END_MASTER( myThid )
542           ENDIF
543    C     Get ready for the next iteration: after adding du/vIce in the first
544    C     iteration, we substract 0.5*du/vIce from u/vIce in the next
545    C     iterations, 0.25*du/vIce in the second, etc.
546           facLS = - 0.5 _d 0 * ABS(facLS)
547          ENDDO
548    C     This is the new residual
549          JFNKresidual = resLoc
550    
551  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */  #endif /* SEAICE_ALLOW_DYNAMICS and SEAICE_CGRID and SEAICE_ALLOW_JFNK */
552    
553        RETURN        RETURN

Legend:
Removed from v.1.9  
changed lines
  Added in v.1.15

  ViewVC Help
Powered by ViewVC 1.1.22