New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
icedyn_rhg_evp.F90 in branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/LIM_SRC_3/icedyn_rhg_evp.F90 @ 8637

Last change on this file since 8637 was 8637, checked in by gm, 7 years ago

#1911 (ENHANCE-09): PART I.3 - phasing with updated branch dev_r8183_ICEMODEL revision 8626

File size: 51.3 KB
Line 
1MODULE icedyn_rhg_evp
2   !!======================================================================
3   !!                     ***  MODULE  icedyn_rhg_evp  ***
4   !!   Sea-Ice dynamics : rheology Elasto-Viscous-Plastic
5   !!======================================================================
6   !! History :   -   !  2007-03  (M.A. Morales Maqueda, S. Bouillon) Original code
7   !!            3.0  !  2008-03  (M. Vancoppenolle) LIM3
8   !!             -   !  2008-11  (M. Vancoppenolle, S. Bouillon, Y. Aksenov) add surface tilt in ice rheolohy
9   !!            3.3  !  2009-05  (G.Garric) addition of the evp cas
10   !!            3.4  !  2011-01  (A. Porter)  dynamical allocation
11   !!            3.5  !  2012-08  (R. Benshila)  AGRIF
12   !!            3.6  !  2016-06  (C. Rousset) Rewriting + landfast ice + possibility to use mEVP (Bouillon 2013)
13   !!----------------------------------------------------------------------
14#if defined key_lim3
15   !!----------------------------------------------------------------------
16   !!   'key_lim3'                                       ESIM sea-ice model
17   !!----------------------------------------------------------------------
18   !!   ice_dyn_rhg_evp       : computes ice velocities from EVP rheology
19   !!----------------------------------------------------------------------
20   USE phycst         ! Physical constant
21   USE dom_oce        ! Ocean domain
22   USE sbc_oce , ONLY : ln_ice_embd, nn_fsbc, ssh_m
23   USE sbc_ice , ONLY : utau_ice, vtau_ice, snwice_mass, snwice_mass_b
24   USE ice            ! sea-ice: ice variables
25   USE icedyn_rdgrft  ! sea-ice: ice strength
26   !
27   USE in_out_manager ! I/O manager
28   USE iom            ! I/O manager library
29   USE lib_mpp        ! MPP library
30   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
31   USE lbclnk         ! lateral boundary conditions (or mpp links)
32   USE prtctl         ! Print control
33   !
34#if defined key_agrif
35   USE agrif_lim3_interp
36#endif
37   USE bdy_oce   , ONLY: ln_bdy 
38   USE bdyice 
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   ice_dyn_rhg_evp   ! called by icedyn_rhg.F90
44   PUBLIC   rhg_evp_rst       ! called by icedyn_rhg.F90
45
46   !! * Substitutions
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/ICE 4.0 , NEMO Consortium (2017)
50   !! $Id: icedyn_rhg_evp.F90 8378 2017-07-26 13:55:59Z clem $
51   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
52   !!----------------------------------------------------------------------
53CONTAINS
54
55   SUBROUTINE ice_dyn_rhg_evp( kt, pstress1_i, pstress2_i, pstress12_i, u_ice, v_ice, pshear_i, pdivu_i, pdelta_i )
56      !!-------------------------------------------------------------------
57      !!                 ***  SUBROUTINE ice_dyn_rhg_evp  ***
58      !!                          EVP-C-grid
59      !!
60      !! ** purpose : determines sea ice drift from wind stress, ice-ocean
61      !!  stress and sea-surface slope. Ice-ice interaction is described by
62      !!  a non-linear elasto-viscous-plastic (EVP) law including shear
63      !!  strength and a bulk rheology (Hunke and Dukowicz, 2002).   
64      !!
65      !!  The points in the C-grid look like this, dear reader
66      !!
67      !!                              (ji,jj)
68      !!                                 |
69      !!                                 |
70      !!                      (ji-1,jj)  |  (ji,jj)
71      !!                             ---------   
72      !!                            |         |
73      !!                            | (ji,jj) |------(ji,jj)
74      !!                            |         |
75      !!                             ---------   
76      !!                     (ji-1,jj-1)     (ji,jj-1)
77      !!
78      !! ** Inputs  : - wind forcing (stress), oceanic currents
79      !!                ice total volume (vt_i) per unit area
80      !!                snow total volume (vt_s) per unit area
81      !!
82      !! ** Action  : - compute u_ice, v_ice : the components of the
83      !!                sea-ice velocity vector
84      !!              - compute delta_i, shear_i, divu_i, which are inputs
85      !!                of the ice thickness distribution
86      !!
87      !! ** Steps   : 0) compute mask at F point
88      !!              1) Compute ice snow mass, ice strength
89      !!              2) Compute wind, oceanic stresses, mass terms and
90      !!                 coriolis terms of the momentum equation
91      !!              3) Solve the momentum equation (iterative procedure)
92      !!              4) Recompute delta, shear and divergence
93      !!                 (which are inputs of the ITD) & store stress
94      !!                 for the next time step
95      !!              5) Diagnostics including charge ellipse
96      !!
97      !! ** Notes   : There is the possibility to use mEVP from Bouillon 2013
98      !!              (by uncommenting some lines in part 3 and changing alpha and beta parameters)
99      !!              but this solution appears very unstable (see Kimmritz et al 2016)
100      !!
101      !! References : Hunke and Dukowicz, JPO97
102      !!              Bouillon et al., Ocean Modelling 2009
103      !!              Bouillon et al., Ocean Modelling 2013
104      !!-------------------------------------------------------------------
105      INTEGER, INTENT(in) ::   kt      ! time step
106      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pstress1_i, pstress2_i, pstress12_i
107      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   u_ice, v_ice 
108      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pshear_i, pdivu_i, pdelta_i 
109      !!
110      INTEGER ::   ji, jj       ! dummy loop indices
111      INTEGER ::   jter         ! local integers
112
113      REAL(wp) ::   zrhoco                                                   ! rau0 * rn_cio
114      REAL(wp) ::   zdtevp, z1_dtevp                                         ! time step for subcycling
115      REAL(wp) ::   ecc2, z1_ecc2                                            ! square of yield ellipse eccenticity
116      REAL(wp) ::   zbeta, zalph1, z1_alph1, zalph2, z1_alph2                ! alpha and beta from Bouillon 2009 and 2013
117      REAL(wp) ::   zm1, zm2, zm3, zmassU, zmassV                            ! ice/snow mass
118      REAL(wp) ::   zdelta, zp_delf, zds2, zdt, zdt2, zdiv, zdiv2            ! temporary scalars
119      REAL(wp) ::   zTauO, zTauB, zTauE, zvel                                ! temporary scalars
120
121      REAL(wp) ::   zresm                                                    ! Maximal error on ice velocity
122      REAL(wp) ::   zintb, zintn                                             ! dummy argument
123      REAL(wp) ::   zfac_x, zfac_y
124      REAL(wp) ::   zshear, zdum1, zdum2
125     
126      REAL(wp), DIMENSION(jpi,jpj) ::   z1_e1t0, z1_e2t0                ! scale factors
127      REAL(wp), DIMENSION(jpi,jpj) ::   zp_delt                         ! P/delta at T points
128      !
129      REAL(wp), DIMENSION(jpi,jpj) ::   zaU   , zaV                     ! ice fraction on U/V points
130      REAL(wp), DIMENSION(jpi,jpj) ::   zmU_t, zmV_t                    ! ice/snow mass/dt on U/V points
131      REAL(wp), DIMENSION(jpi,jpj) ::   zmf                             ! coriolis parameter at T points
132      REAL(wp), DIMENSION(jpi,jpj) ::   zTauU_ia , ztauV_ia             ! ice-atm. stress at U-V points
133      REAL(wp), DIMENSION(jpi,jpj) ::   zspgU , zspgV                   ! surface pressure gradient at U/V points
134      REAL(wp), DIMENSION(jpi,jpj) ::   v_oceU, u_oceV, v_iceU, u_iceV  ! ocean/ice u/v component on V/U points                           
135      REAL(wp), DIMENSION(jpi,jpj) ::   zfU   , zfV                     ! internal stresses
136     
137      REAL(wp), DIMENSION(jpi,jpj) ::   zds                             ! shear
138      REAL(wp), DIMENSION(jpi,jpj) ::   zs1, zs2, zs12                  ! stress tensor components
139      REAL(wp), DIMENSION(jpi,jpj) ::   zu_ice, zv_ice, zresr           ! check convergence
140      REAL(wp), DIMENSION(jpi,jpj) ::   zpice                           ! array used for the calculation of ice surface slope:
141                                                                             !   ocean surface (ssh_m) if ice is not embedded
142                                                                             !   ice top surface if ice is embedded   
143      REAL(wp), DIMENSION(jpi,jpj) ::   zCorx, zCory                    ! Coriolis stress array
144      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_oi, ztauy_oi              ! Ocean-to-ice stress array
145
146      REAL(wp), DIMENSION(jpi,jpj) ::   zswitchU, zswitchV              ! dummy arrays
147      REAL(wp), DIMENSION(jpi,jpj) ::   zmaskU, zmaskV                  ! mask for ice presence
148      REAL(wp), DIMENSION(jpi,jpj) ::   zfmask, zwf                     ! mask at F points for the ice
149
150      REAL(wp), PARAMETER          ::   zepsi  = 1.0e-20_wp             ! tolerance parameter
151      REAL(wp), PARAMETER          ::   zmmin  = 1._wp                  ! ice mass (kg/m2) below which ice velocity equals ocean velocity
152      !! --- diags
153      REAL(wp), DIMENSION(jpi,jpj) ::   zswi
154      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zsig1, zsig2, zsig3
155      !! --- SIMIP diags
156      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_sig1      ! Average normal stress in sea ice   
157      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_sig2      ! Maximum shear stress in sea ice
158      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_dssh_dx   ! X-direction sea-surface tilt term (N/m2)
159      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_dssh_dy   ! X-direction sea-surface tilt term (N/m2)
160      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_corstrx   ! X-direction coriolis stress (N/m2)
161      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_corstry   ! Y-direction coriolis stress (N/m2)
162      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_intstrx   ! X-direction internal stress (N/m2)
163      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_intstry   ! Y-direction internal stress (N/m2)
164      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_utau_oi   ! X-direction ocean-ice stress
165      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_vtau_oi   ! Y-direction ocean-ice stress 
166      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s)
167      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s)
168      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s)
169      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s)
170      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xatrp     ! X-component of area transport (m2/s)
171      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_yatrp     ! Y-component of area transport (m2/s)     
172      !!-------------------------------------------------------------------
173
174      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_dyn_rhg_evp: EVP sea-ice rheology'
175
176#if defined key_agrif 
177      CALL agrif_interp_lim3( 'U', 0, nn_nevp )   ! First interpolation of coarse values
178      CALL agrif_interp_lim3( 'V', 0, nn_nevp )
179#endif
180      !
181      !------------------------------------------------------------------------------!
182      ! 0) mask at F points for the ice
183      !------------------------------------------------------------------------------!
184      ! ocean/land mask
185      DO jj = 1, jpjm1
186         DO ji = 1, jpim1      ! NO vector opt.
187            zfmask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
188         END DO
189      END DO
190      CALL lbc_lnk( zfmask, 'F', 1._wp )
191
192      ! Lateral boundary conditions on velocity (modify zfmask)
193      zwf(:,:) = zfmask(:,:)
194      DO jj = 2, jpjm1
195         DO ji = fs_2, fs_jpim1   ! vector opt.
196            IF( zfmask(ji,jj) == 0._wp ) THEN
197               zfmask(ji,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1), zwf(ji-1,jj), zwf(ji,jj-1) ) )
198            ENDIF
199         END DO
200      END DO
201      DO jj = 2, jpjm1
202         IF( zfmask(1,jj) == 0._wp ) THEN
203            zfmask(1  ,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
204         ENDIF
205         IF( zfmask(jpi,jj) == 0._wp ) THEN
206            zfmask(jpi,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
207         ENDIF
208      END DO
209      DO ji = 2, jpim1
210         IF( zfmask(ji,1) == 0._wp ) THEN
211            zfmask(ji,1  ) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
212         ENDIF
213         IF( zfmask(ji,jpj) == 0._wp ) THEN
214            zfmask(ji,jpj) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
215         ENDIF
216      END DO
217      CALL lbc_lnk( zfmask, 'F', 1._wp )
218
219      !------------------------------------------------------------------------------!
220      ! 1) define some variables and initialize arrays
221      !------------------------------------------------------------------------------!
222      zrhoco = rau0 * rn_cio 
223
224      ! ecc2: square of yield ellipse eccenticrity
225      ecc2    = rn_ecc * rn_ecc
226      z1_ecc2 = 1._wp / ecc2
227
228      ! Time step for subcycling
229      zdtevp   = rdt_ice / REAL( nn_nevp )
230      z1_dtevp = 1._wp / zdtevp
231
232      ! alpha parameters (Bouillon 2009)
233      zalph1 = ( 2._wp * rn_relast * rdt_ice ) * z1_dtevp
234      zalph2 = zalph1 * z1_ecc2
235
236      ! alpha and beta parameters (Bouillon 2013)
237      !!zalph1 = 40.
238      !!zalph2 = 40.
239      !!zbeta  = 3000.
240      !!zbeta = REAL( nn_nevp )   ! close to classical EVP of Hunke (2001)
241
242      z1_alph1 = 1._wp / ( zalph1 + 1._wp )
243      z1_alph2 = 1._wp / ( zalph2 + 1._wp )
244
245      ! Initialise stress tensor
246      zs1 (:,:) = pstress1_i (:,:) 
247      zs2 (:,:) = pstress2_i (:,:)
248      zs12(:,:) = pstress12_i(:,:)
249
250      ! Ice strength
251      CALL ice_strength
252
253      ! scale factors
254      DO jj = 2, jpjm1
255         DO ji = fs_2, fs_jpim1
256            z1_e1t0(ji,jj) = 1._wp / ( e1t(ji+1,jj  ) + e1t(ji,jj  ) )
257            z1_e2t0(ji,jj) = 1._wp / ( e2t(ji  ,jj+1) + e2t(ji,jj  ) )
258         END DO
259      END DO
260           
261      !
262      !------------------------------------------------------------------------------!
263      ! 2) Wind / ocean stress, mass terms, coriolis terms
264      !------------------------------------------------------------------------------!
265
266      IF( ln_ice_embd ) THEN             !== embedded sea ice: compute representative ice top surface ==!
267         !                                           
268         ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[n/nn_fsbc], n=0,nn_fsbc-1}
269         !                                               = (1/nn_fsbc)^2 * {SUM[n], n=0,nn_fsbc-1}
270         zintn = REAL( nn_fsbc - 1 ) / REAL( nn_fsbc ) * 0.5_wp     
271         !
272         ! average interpolation coeff as used in dynspg = (1/nn_fsbc) * {SUM[1-n/nn_fsbc], n=0,nn_fsbc-1}
273         !                                               = (1/nn_fsbc)^2 * (nn_fsbc^2 - {SUM[n], n=0,nn_fsbc-1})
274         zintb = REAL( nn_fsbc + 1 ) / REAL( nn_fsbc ) * 0.5_wp
275         !
276         zpice(:,:) = ssh_m(:,:) + ( zintn * snwice_mass(:,:) + zintb * snwice_mass_b(:,:) ) * r1_rau0
277         !
278      ELSE                                    !== non-embedded sea ice: use ocean surface for slope calculation ==!
279         zpice(:,:) = ssh_m(:,:)
280      ENDIF
281
282      DO jj = 2, jpjm1
283         DO ji = fs_2, fs_jpim1
284
285            ! ice fraction at U-V points
286            zaU(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e1e2t(ji,jj) + at_i(ji+1,jj) * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1)
287            zaV(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e1e2t(ji,jj) + at_i(ji,jj+1) * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1)
288
289            ! Ice/snow mass at U-V points
290            zm1 = ( rhosn * vt_s(ji  ,jj  ) + rhoic * vt_i(ji  ,jj  ) )
291            zm2 = ( rhosn * vt_s(ji+1,jj  ) + rhoic * vt_i(ji+1,jj  ) )
292            zm3 = ( rhosn * vt_s(ji  ,jj+1) + rhoic * vt_i(ji  ,jj+1) )
293            zmassU = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm2 * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1)
294            zmassV = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm3 * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1)
295
296            ! Ocean currents at U-V points
297            v_oceU(ji,jj)   = 0.5_wp * ( ( v_oce(ji  ,jj) + v_oce(ji  ,jj-1) ) * e1t(ji+1,jj)    &
298               &                       + ( v_oce(ji+1,jj) + v_oce(ji+1,jj-1) ) * e1t(ji  ,jj) ) * z1_e1t0(ji,jj) * umask(ji,jj,1)
299           
300            u_oceV(ji,jj)   = 0.5_wp * ( ( u_oce(ji,jj  ) + u_oce(ji-1,jj  ) ) * e2t(ji,jj+1)    &
301               &                       + ( u_oce(ji,jj+1) + u_oce(ji-1,jj+1) ) * e2t(ji,jj  ) ) * z1_e2t0(ji,jj) * vmask(ji,jj,1)
302
303            ! Coriolis at T points (m*f)
304            zmf(ji,jj)      = zm1 * ff_t(ji,jj)
305
306            ! m/dt
307            zmU_t(ji,jj)    = zmassU * z1_dtevp
308            zmV_t(ji,jj)    = zmassV * z1_dtevp
309
310            ! Drag ice-atm.
311            zTauU_ia(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj)
312            zTauV_ia(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj)
313
314            ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points
315            zspgU(ji,jj)    = - zmassU * grav * ( zpice(ji+1,jj) - zpice(ji,jj) ) * r1_e1u(ji,jj)
316            zspgV(ji,jj)    = - zmassV * grav * ( zpice(ji,jj+1) - zpice(ji,jj) ) * r1_e2v(ji,jj)
317
318            ! masks
319            zmaskU(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) )  ! 0 if no ice
320            zmaskV(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) )  ! 0 if no ice
321
322            ! switches
323            zswitchU(ji,jj) = MAX( 0._wp, SIGN( 1._wp, zmassU - zmmin ) ) ! 0 if ice mass < zmmin
324            zswitchV(ji,jj) = MAX( 0._wp, SIGN( 1._wp, zmassV - zmmin ) ) ! 0 if ice mass < zmmin
325
326         END DO
327      END DO
328      CALL lbc_lnk( zmf, 'T', 1. )
329      !
330      !------------------------------------------------------------------------------!
331      ! 3) Solution of the momentum equation, iterative procedure
332      !------------------------------------------------------------------------------!
333      !
334      !                                               !----------------------!
335      DO jter = 1 , nn_nevp                           !    loop over jter    !
336         !                                            !----------------------!       
337         IF(ln_ctl) THEN   ! Convergence test
338            DO jj = 1, jpjm1
339               zu_ice(:,jj) = u_ice(:,jj) ! velocity at previous time step
340               zv_ice(:,jj) = v_ice(:,jj)
341            END DO
342         ENDIF
343
344         ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- !
345         DO jj = 1, jpjm1         ! loops start at 1 since there is no boundary condition (lbc_lnk) at i=1 and j=1 for F points
346            DO ji = 1, jpim1
347
348               ! shear at F points
349               zds(ji,jj) = ( ( u_ice(ji,jj+1) * r1_e1u(ji,jj+1) - u_ice(ji,jj) * r1_e1u(ji,jj) ) * e1f(ji,jj) * e1f(ji,jj)   &
350                  &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
351                  &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
352
353            END DO
354         END DO
355         CALL lbc_lnk( zds, 'F', 1. )
356
357         DO jj = 2, jpj    ! loop to jpi,jpj to avoid making a communication for zs1,zs2,zs12
358            DO ji = 2, jpi ! no vector loop
359
360               ! shear**2 at T points (doc eq. A16)
361               zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
362                  &   + zds(ji,jj-1) * zds(ji,jj-1) * e1e2f(ji,jj-1) + zds(ji-1,jj-1) * zds(ji-1,jj-1) * e1e2f(ji-1,jj-1)  &
363                  &   ) * 0.25_wp * r1_e1e2t(ji,jj)
364             
365               ! divergence at T points
366               zdiv  = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
367                  &    + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
368                  &    ) * r1_e1e2t(ji,jj)
369               zdiv2 = zdiv * zdiv
370               
371               ! tension at T points
372               zdt  = ( ( u_ice(ji,jj) * r1_e2u(ji,jj) - u_ice(ji-1,jj) * r1_e2u(ji-1,jj) ) * e2t(ji,jj) * e2t(ji,jj)   &
373                  &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
374                  &   ) * r1_e1e2t(ji,jj)
375               zdt2 = zdt * zdt
376               
377               ! delta at T points
378               zdelta = SQRT( zdiv2 + ( zdt2 + zds2 ) * z1_ecc2 ) 
379
380               ! P/delta at T points
381               zp_delt(ji,jj) = strength(ji,jj) / ( zdelta + rn_creepl )
382               
383               ! stress at T points
384               zs1(ji,jj) = ( zs1(ji,jj) * zalph1 + zp_delt(ji,jj) * ( zdiv - zdelta ) ) * z1_alph1
385               zs2(ji,jj) = ( zs2(ji,jj) * zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 ) ) * z1_alph2
386             
387            END DO
388         END DO
389         CALL lbc_lnk( zp_delt, 'T', 1. )
390
391         DO jj = 1, jpjm1
392            DO ji = 1, jpim1
393
394               ! P/delta at F points
395               zp_delf = 0.25_wp * ( zp_delt(ji,jj) + zp_delt(ji+1,jj) + zp_delt(ji,jj+1) + zp_delt(ji+1,jj+1) )
396               
397               ! stress at F points
398               zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 ) * 0.5_wp ) * z1_alph2
399
400            END DO
401         END DO
402
403         ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- !
404         DO jj = 2, jpjm1
405            DO ji = fs_2, fs_jpim1               
406               !                   !--- U points
407               zfU(ji,jj) = 0.5_wp * ( ( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj)                                             &
408                  &                  + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj)    &
409                  &                    ) * r1_e2u(ji,jj)                                                                      &
410                  &                  + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1)  &
411                  &                    ) * 2._wp * r1_e1u(ji,jj)                                                              &
412                  &                  ) * r1_e1e2u(ji,jj)
413                  !
414                  !                !--- V points
415               zfV(ji,jj) = 0.5_wp * ( ( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj)                                             &
416                  &                  - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj)    &
417                  &                    ) * r1_e1v(ji,jj)                                                                      &
418                  &                  + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj)  &
419                  &                    ) * 2._wp * r1_e2v(ji,jj)                                                              &
420                  &                  ) * r1_e1e2v(ji,jj)
421                  !
422                  !                !--- u_ice at V point
423               u_iceV(ji,jj) = 0.5_wp * ( ( u_ice(ji,jj  ) + u_ice(ji-1,jj  ) ) * e2t(ji,jj+1)     &
424                  &                     + ( u_ice(ji,jj+1) + u_ice(ji-1,jj+1) ) * e2t(ji,jj  ) ) * z1_e2t0(ji,jj) * vmask(ji,jj,1)
425                  !
426                  !                !--- v_ice at U point
427               v_iceU(ji,jj) = 0.5_wp * ( ( v_ice(ji  ,jj) + v_ice(ji  ,jj-1) ) * e1t(ji+1,jj)     &
428                  &                     + ( v_ice(ji+1,jj) + v_ice(ji+1,jj-1) ) * e1t(ji  ,jj) ) * z1_e1t0(ji,jj) * umask(ji,jj,1)
429               !
430            END DO
431         END DO
432         !
433         ! --- Computation of ice velocity --- !
434         !  Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta are chosen wisely and large nn_nevp
435         !  Bouillon et al. 2009 (eq 34-35) => stable
436         IF( MOD(jter,2) == 0 ) THEN ! even iterations
437            !
438            DO jj = 2, jpjm1
439               DO ji = fs_2, fs_jpim1
440                     !                 !--- tau_io/(v_oce - v_ice)
441                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
442                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
443                     !                 !--- Ocean-to-Ice stress
444                  ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
445                     !
446                     !                 !--- tau_bottom/v_ice
447                  zvel  = MAX( zepsi, SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) ) )
448                  zTauB = - tau_icebfr(ji,jj) / zvel
449                     !
450                     !                 !--- Coriolis at V-points (energy conserving formulation)
451                  zCory(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
452                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
453                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
454                     !
455                     !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
456                  zTauE = zfV(ji,jj) + zTauV_ia(ji,jj) + zCory(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
457                     !
458                     !                 !--- landfast switch => 0 = static friction ; 1 = sliding friction
459                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, ztauE - tau_icebfr(ji,jj) ) - SIGN( 1._wp, zTauE ) ) )
460                     !
461                     !                 !--- ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
462                  v_ice(ji,jj) = ( (           rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                             &  ! previous velocity
463                     &                                     + zTauE + zTauO * v_ice(ji,jj)                             &  ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
464                     &                                     ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )             &  ! m/dt + tau_io(only ice part) + landfast
465                     &              + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )  &  ! static friction => slow decrease to v=0
466                     &              ) * zswitchV(ji,jj) + v_oce(ji,jj) * ( 1._wp - zswitchV(ji,jj) )                  &  ! v_ice = v_oce if mass < zmmin
467                     &            ) * zmaskV(ji,jj)
468                     !
469                  ! Bouillon 2013
470                  !!v_ice(ji,jj) = ( zmV_t(ji,jj) * ( zbeta * v_ice(ji,jj) + v_ice_b(ji,jj) )                  &
471                  !!   &           + zfV(ji,jj) + zCory(ji,jj) + zTauV_ia(ji,jj) + zTauO * v_oce(ji,jj) + zspgV(ji,jj)  &
472                  !!   &           ) / MAX( zmV_t(ji,jj) * ( zbeta + 1._wp ) + zTauO - zTauB, zepsi ) * zswitchV(ji,jj)
473                  !
474               END DO
475            END DO
476            CALL lbc_lnk( v_ice, 'V', -1. )
477            !
478#if defined key_agrif
479!!            CALL agrif_interp_lim3( 'V', jter, nn_nevp )
480            CALL agrif_interp_lim3( 'V' )
481#endif
482            IF( ln_bdy ) CALL bdy_ice_dyn( 'V' )
483            !
484            DO jj = 2, jpjm1
485               DO ji = fs_2, fs_jpim1
486                               
487                  ! tau_io/(u_oce - u_ice)
488                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
489                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
490
491                  ! Ocean-to-Ice stress
492                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
493
494                  ! tau_bottom/u_ice
495                  zvel  = MAX( zepsi, SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) ) )
496                  zTauB = - tau_icebfr(ji,jj) / zvel
497
498                  ! Coriolis at U-points (energy conserving formulation)
499                  zCorx(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
500                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
501                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
502                 
503                  ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
504                  zTauE = zfU(ji,jj) + zTauU_ia(ji,jj) + zCorx(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
505
506                  ! landfast switch => 0 = static friction ; 1 = sliding friction
507                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, ztauE - tau_icebfr(ji,jj) ) - SIGN( 1._wp, zTauE ) ) )
508
509                  ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
510                  u_ice(ji,jj) = ( (           rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                             &  ! previous velocity
511                     &                                     + zTauE + zTauO * u_ice(ji,jj)                             &  ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
512                     &                                     ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )             &  ! m/dt + tau_io(only ice part) + landfast
513                     &              + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )  &  ! static friction => slow decrease to v=0
514                     &              ) * zswitchU(ji,jj) + u_oce(ji,jj) * ( 1._wp - zswitchU(ji,jj) )                  &  ! v_ice = v_oce if mass < zmmin
515                     &            ) * zmaskU(ji,jj)
516                  ! Bouillon 2013
517                  !!u_ice(ji,jj) = ( zmU_t(ji,jj) * ( zbeta * u_ice(ji,jj) + u_ice_b(ji,jj) )                  &
518                  !!   &           + zfU(ji,jj) + zCorx(ji,jj) + zTauU_ia(ji,jj) + zTauO * u_oce(ji,jj) + zspgU(ji,jj)  &
519                  !!   &           ) / MAX( zmU_t(ji,jj) * ( zbeta + 1._wp ) + zTauO - zTauB, zepsi ) * zswitchU(ji,jj)
520               END DO
521            END DO
522            CALL lbc_lnk( u_ice, 'U', -1. )
523            !
524#if defined key_agrif
525!!            CALL agrif_interp_lim3( 'U', jter, nn_nevp )
526            CALL agrif_interp_lim3( 'U' )
527#endif
528            IF( ln_bdy ) CALL bdy_ice_dyn( 'U' )
529            !
530         ELSE ! odd iterations
531            !
532            DO jj = 2, jpjm1
533               DO ji = fs_2, fs_jpim1
534
535                  ! tau_io/(u_oce - u_ice)
536                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
537                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
538
539                  ! Ocean-to-Ice stress
540                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
541
542                  ! tau_bottom/u_ice
543                  zvel  = MAX( zepsi, SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) ) )
544                  zTauB = - tau_icebfr(ji,jj) / zvel
545
546                  ! Coriolis at U-points (energy conserving formulation)
547                  zCorx(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
548                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
549                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
550
551                  ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
552                  zTauE = zfU(ji,jj) + zTauU_ia(ji,jj) + zCorx(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
553
554                  ! landfast switch => 0 = static friction ; 1 = sliding friction
555                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, ztauE - tau_icebfr(ji,jj) ) - SIGN( 1._wp, zTauE ) ) )
556
557                  ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
558                  u_ice(ji,jj) = ( (           rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                             &  ! previous velocity
559                     &                                     + zTauE + zTauO * u_ice(ji,jj)                             &  ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
560                     &                                     ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )             &  ! m/dt + tau_io(only ice part) + landfast
561                     &              + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )  &  ! static friction => slow decrease to v=0
562                     &              ) * zswitchU(ji,jj) + u_oce(ji,jj) * ( 1._wp - zswitchU(ji,jj) )                  &  ! v_ice = v_oce if mass < zmmin
563                     &            ) * zmaskU(ji,jj)
564                  ! Bouillon 2013
565                  !!u_ice(ji,jj) = ( zmU_t(ji,jj) * ( zbeta * u_ice(ji,jj) + u_ice_b(ji,jj) )                  &
566                  !!   &           + zfU(ji,jj) + zCorx(ji,jj) + zTauU_ia(ji,jj) + zTauO * u_oce(ji,jj) + zspgU(ji,jj)  &
567                  !!   &           ) / MAX( zmU_t(ji,jj) * ( zbeta + 1._wp ) + zTauO - zTauB, zepsi ) * zswitchU(ji,jj)
568               END DO
569            END DO
570            CALL lbc_lnk( u_ice, 'U', -1. )
571            !
572#if defined key_agrif
573!!            CALL agrif_interp_lim3( 'U', jter, nn_nevp )
574            CALL agrif_interp_lim3( 'U' )
575#endif
576            IF( ln_bdy ) CALL bdy_ice_dyn( 'U' )
577            !
578            DO jj = 2, jpjm1
579               DO ji = fs_2, fs_jpim1
580         
581                  ! tau_io/(v_oce - v_ice)
582                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
583                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
584
585                  ! Ocean-to-Ice stress
586                  ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
587
588                  ! tau_bottom/v_ice
589                  zvel  = MAX( zepsi, SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) ) )
590                  ztauB = - tau_icebfr(ji,jj) / zvel
591                 
592                  ! Coriolis at V-points (energy conserving formulation)
593                  zCory(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
594                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
595                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
596
597                  ! Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
598                  zTauE = zfV(ji,jj) + zTauV_ia(ji,jj) + zCory(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
599
600                  ! landfast switch => 0 = static friction (tau_icebfr > zTauE); 1 = sliding friction
601                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zTauE - tau_icebfr(ji,jj) ) - SIGN( 1._wp, zTauE ) ) )
602                 
603                  ! ice velocity using implicit formulation (cf Madec doc & Bouillon 2009)
604                  v_ice(ji,jj) = ( (           rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                             &  ! previous velocity
605                     &                                     + zTauE + zTauO * v_ice(ji,jj)                             &  ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
606                     &                                     ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )             &  ! m/dt + tau_io(only ice part) + landfast
607                     &              + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )  &  ! static friction => slow decrease to v=0
608                     &              ) * zswitchV(ji,jj) + v_oce(ji,jj) * ( 1._wp - zswitchV(ji,jj) )                  &  ! v_ice = v_oce if mass < zmmin
609                     &            ) * zmaskV(ji,jj)
610                  ! Bouillon 2013
611                  !!v_ice(ji,jj) = ( zmV_t(ji,jj) * ( zbeta * v_ice(ji,jj) + v_ice_b(ji,jj) )                  &
612                  !!   &           + zfV(ji,jj) + zCory(ji,jj) + zTauV_ia(ji,jj) + zTauO * v_oce(ji,jj) + zspgV(ji,jj)  &
613                  !!   &           ) / MAX( zmV_t(ji,jj) * ( zbeta + 1._wp ) + zTauO - zTauB, zepsi ) * zswitchV(ji,jj)
614               END DO
615            END DO
616            CALL lbc_lnk( v_ice, 'V', -1. )
617            !
618#if defined key_agrif
619!!            CALL agrif_interp_lim3( 'V', jter, nn_nevp )
620            CALL agrif_interp_lim3( 'V' )
621#endif
622            IF( ln_bdy ) CALL bdy_ice_dyn( 'V' )
623            !
624         ENDIF
625         
626         IF(ln_ctl) THEN   ! Convergence test
627            DO jj = 2 , jpjm1
628               zresr(:,jj) = MAX( ABS( u_ice(:,jj) - zu_ice(:,jj) ), ABS( v_ice(:,jj) - zv_ice(:,jj) ) )
629            END DO
630            zresm = MAXVAL( zresr( 1:jpi, 2:jpjm1 ) )
631            IF( lk_mpp )   CALL mpp_max( zresm )   ! max over the global domain
632         ENDIF
633         !
634         !                                                ! ==================== !
635      END DO                                              !  end loop over jter  !
636      !                                                   ! ==================== !
637      !
638      !------------------------------------------------------------------------------!
639      ! 4) Recompute delta, shear and div (inputs for mechanical redistribution)
640      !------------------------------------------------------------------------------!
641      DO jj = 1, jpjm1
642         DO ji = 1, jpim1
643
644            ! shear at F points
645            zds(ji,jj) = ( ( u_ice(ji,jj+1) * r1_e1u(ji,jj+1) - u_ice(ji,jj) * r1_e1u(ji,jj) ) * e1f(ji,jj) * e1f(ji,jj)   &
646               &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
647               &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
648
649         END DO
650      END DO           
651     
652      DO jj = 2, jpjm1
653         DO ji = 2, jpim1 ! no vector loop
654           
655            ! tension**2 at T points
656            zdt  = ( ( u_ice(ji,jj) * r1_e2u(ji,jj) - u_ice(ji-1,jj) * r1_e2u(ji-1,jj) ) * e2t(ji,jj) * e2t(ji,jj)   &
657               &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
658               &   ) * r1_e1e2t(ji,jj)
659            zdt2 = zdt * zdt
660           
661            ! shear**2 at T points (doc eq. A16)
662            zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
663               &   + zds(ji,jj-1) * zds(ji,jj-1) * e1e2f(ji,jj-1) + zds(ji-1,jj-1) * zds(ji-1,jj-1) * e1e2f(ji-1,jj-1)  &
664               &   ) * 0.25_wp * r1_e1e2t(ji,jj)
665           
666            ! shear at T points
667            pshear_i(ji,jj) = SQRT( zdt2 + zds2 )
668
669            ! divergence at T points
670            pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
671               &             + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
672               &             ) * r1_e1e2t(ji,jj)
673           
674            ! delta at T points
675            zdelta         = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) 
676            rswitch        = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zdelta ) ) ! 0 if delta=0
677            pdelta_i(ji,jj) = zdelta + rn_creepl * rswitch
678
679         END DO
680      END DO
681      CALL lbc_lnk_multi( pshear_i, 'T', 1., pdivu_i, 'T', 1., pdelta_i, 'T', 1. )
682     
683      ! --- Store the stress tensor for the next time step --- !
684      CALL lbc_lnk_multi( zs1, 'T', 1., zs2, 'T', 1., zs12, 'F', 1. )
685      pstress1_i (:,:) = zs1 (:,:)
686      pstress2_i (:,:) = zs2 (:,:)
687      pstress12_i(:,:) = zs12(:,:)
688      !
689
690      !------------------------------------------------------------------------------!
691      ! 5) diagnostics
692      !------------------------------------------------------------------------------!
693      DO jj = 1, jpj
694         DO ji = 1, jpi
695            zswi(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) ! 1 if ice, 0 if no ice
696         END DO
697      END DO
698
699      ! --- divergence, shear and strength --- !
700      IF( iom_use('idive')  )   CALL iom_put( "idive"  , pdivu_i(:,:)  * zswi(:,:) )   ! divergence
701      IF( iom_use('ishear') )   CALL iom_put( "ishear" , pshear_i(:,:) * zswi(:,:) )   ! shear
702      IF( iom_use('icestr') )   CALL iom_put( "icestr" , strength(:,:) * zswi(:,:) )   ! Ice strength
703
704      ! --- charge ellipse --- !
705      IF( iom_use('isig1') .OR. iom_use('isig2') .OR. iom_use('isig3') ) THEN
706         !
707         ALLOCATE( zsig1(jpi,jpj) , zsig2(jpi,jpj) , zsig3(jpi,jpj) )
708         !         
709         DO jj = 2, jpjm1
710            DO ji = 2, jpim1
711               zdum1 = ( zswi(ji-1,jj) * pstress12_i(ji-1,jj) + zswi(ji  ,jj-1) * pstress12_i(ji  ,jj-1) +  &  ! stress12_i at T-point
712                  &      zswi(ji  ,jj) * pstress12_i(ji  ,jj) + zswi(ji-1,jj-1) * pstress12_i(ji-1,jj-1) )  &
713                  &    / MAX( 1._wp, zswi(ji-1,jj) + zswi(ji,jj-1) + zswi(ji,jj) + zswi(ji-1,jj-1) )
714
715               zshear = SQRT( pstress2_i(ji,jj) * pstress2_i(ji,jj) + 4._wp * zdum1 * zdum1 ) ! shear stress 
716
717               zdum2 = zswi(ji,jj) / MAX( 1._wp, strength(ji,jj) )
718
719!!               zsig1(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) + zshear ) ! principal stress (y-direction, see Hunke & Dukowicz 2002)
720!!               zsig2(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) - zshear ) ! principal stress (x-direction, see Hunke & Dukowicz 2002)
721!!               zsig3(ji,jj) = zdum2**2 * ( ( pstress1_i(ji,jj) + strength(ji,jj) )**2 + ( rn_ecc * zshear )**2 ) ! quadratic relation linking compressive stress to shear stress
722!!                                                                                                               ! (scheme converges if this value is ~1, see Bouillon et al 2009 (eq. 11))
723               zsig1(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) )          ! compressive stress, see Bouillon et al. 2015
724               zsig2(ji,jj) = 0.5_wp * zdum2 * ( zshear )                     ! shear stress
725               zsig3(ji,jj) = zdum2**2 * ( ( pstress1_i(ji,jj) + strength(ji,jj) )**2 + ( rn_ecc * zshear )**2 )
726            END DO
727         END DO
728         CALL lbc_lnk_multi( zsig1, 'T', 1., zsig2, 'T', 1., zsig3, 'T', 1. )
729         !
730         IF( iom_use('isig1') )   CALL iom_put( "isig1" , zsig1 )
731         IF( iom_use('isig2') )   CALL iom_put( "isig2" , zsig2 )
732         IF( iom_use('isig3') )   CALL iom_put( "isig3" , zsig3 )
733         !
734         DEALLOCATE( zsig1 , zsig2 , zsig3 )
735      ENDIF
736     
737      ! --- SIMIP --- !
738      IF ( iom_use( 'normstr'  ) .OR. iom_use( 'sheastr'  ) .OR. iom_use( 'dssh_dx'  ) .OR. iom_use( 'dssh_dy'  ) .OR. &
739         & iom_use( 'corstrx'  ) .OR. iom_use( 'corstry'  ) .OR. iom_use( 'intstrx'  ) .OR. iom_use( 'intstry'  ) .OR. &
740         & iom_use( 'utau_oi'  ) .OR. iom_use( 'vtau_oi'  ) .OR. iom_use( 'xmtrpice' ) .OR. iom_use( 'ymtrpice' ) .OR. &
741         & iom_use( 'xmtrpsnw' ) .OR. iom_use( 'ymtrpsnw' ) .OR. iom_use( 'xatrp'    ) .OR. iom_use( 'yatrp'    ) ) THEN
742
743         ALLOCATE( zdiag_sig1     (jpi,jpj) , zdiag_sig2     (jpi,jpj) , zdiag_dssh_dx  (jpi,jpj) , zdiag_dssh_dy  (jpi,jpj) ,  &
744            &      zdiag_corstrx  (jpi,jpj) , zdiag_corstry  (jpi,jpj) , zdiag_intstrx  (jpi,jpj) , zdiag_intstry  (jpi,jpj) ,  &
745            &      zdiag_utau_oi  (jpi,jpj) , zdiag_vtau_oi  (jpi,jpj) , zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) ,  &
746            &      zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp    (jpi,jpj) , zdiag_yatrp    (jpi,jpj) )
747         
748         DO jj = 2, jpjm1
749            DO ji = 2, jpim1
750               rswitch  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) ! 1 if ice, 0 if no ice
751               
752               ! Stress tensor invariants (normal and shear stress N/m)
753               zdiag_sig1(ji,jj) = ( zs1(ji,jj) + zs2(ji,jj) ) * rswitch                                 ! normal stress
754               zdiag_sig2(ji,jj) = SQRT( ( zs1(ji,jj) - zs2(ji,jj) )**2 + 4*zs12(ji,jj)**2 ) * rswitch   ! shear stress
755               
756               ! Stress terms of the momentum equation (N/m2)
757               zdiag_dssh_dx(ji,jj) = zspgU(ji,jj) * rswitch     ! sea surface slope stress term
758               zdiag_dssh_dy(ji,jj) = zspgV(ji,jj) * rswitch
759               
760               zdiag_corstrx(ji,jj) = zCorx(ji,jj) * rswitch     ! Coriolis stress term
761               zdiag_corstry(ji,jj) = zCory(ji,jj) * rswitch
762               
763               zdiag_intstrx(ji,jj) = zfU(ji,jj)   * rswitch     ! internal stress term
764               zdiag_intstry(ji,jj) = zfV(ji,jj)   * rswitch
765               
766               zdiag_utau_oi(ji,jj) = ztaux_oi(ji,jj) * rswitch  ! oceanic stress
767               zdiag_vtau_oi(ji,jj) = ztauy_oi(ji,jj) * rswitch
768               
769               ! 2D ice mass, snow mass, area transport arrays (X, Y)
770               zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * rswitch
771               zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * rswitch
772               
773               zdiag_xmtrp_ice(ji,jj) = rhoic * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component
774               zdiag_ymtrp_ice(ji,jj) = rhoic * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) !        ''           Y-   ''
775               
776               zdiag_xmtrp_snw(ji,jj) = rhosn * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
777               zdiag_ymtrp_snw(ji,jj) = rhosn * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) !          ''          Y-   ''
778               
779               zdiag_xatrp(ji,jj)     = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) )         ! area transport,      X-component
780               zdiag_yatrp(ji,jj)     = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) )         !        ''            Y-   ''
781               
782            END DO
783         END DO
784         
785         CALL lbc_lnk_multi(   zdiag_sig1   , 'T',  1., zdiag_sig2   , 'T',  1.,   &
786            &                  zdiag_dssh_dx, 'U', -1., zdiag_dssh_dy, 'V', -1.,   &
787            &                  zdiag_corstrx, 'U', -1., zdiag_corstry, 'V', -1.,   & 
788            &                  zdiag_intstrx, 'U', -1., zdiag_intstry, 'V', -1.    )
789         
790         CALL lbc_lnk_multi(   zdiag_utau_oi, 'U', -1., zdiag_vtau_oi, 'V', -1.    )
791         
792         CALL lbc_lnk_multi(   zdiag_xmtrp_ice, 'U', -1., zdiag_xmtrp_snw, 'U', -1.,   &
793            &                  zdiag_xatrp    , 'U', -1., zdiag_ymtrp_ice, 'V', -1.,   &
794            &                  zdiag_ymtrp_snw, 'V', -1., zdiag_yatrp    , 'V', -1.    )
795         
796         IF ( iom_use( 'normstr'  ) )   CALL iom_put( 'normstr'  ,  zdiag_sig1(:,:)      )   ! Normal stress
797         IF ( iom_use( 'sheastr'  ) )   CALL iom_put( 'sheastr'  ,  zdiag_sig2(:,:)      )   ! Shear stress
798         IF ( iom_use( 'dssh_dx'  ) )   CALL iom_put( 'dssh_dx'  ,  zdiag_dssh_dx(:,:)   )   ! Sea-surface tilt term in force balance (x)
799         IF ( iom_use( 'dssh_dy'  ) )   CALL iom_put( 'dssh_dy'  ,  zdiag_dssh_dy(:,:)   )   ! Sea-surface tilt term in force balance (y)
800         IF ( iom_use( 'corstrx'  ) )   CALL iom_put( 'corstrx'  ,  zdiag_corstrx(:,:)   )   ! Coriolis force term in force balance (x)
801         IF ( iom_use( 'corstry'  ) )   CALL iom_put( 'corstry'  ,  zdiag_corstry(:,:)   )   ! Coriolis force term in force balance (y)
802         IF ( iom_use( 'intstrx'  ) )   CALL iom_put( 'intstrx'  ,  zdiag_intstrx(:,:)   )   ! Internal force term in force balance (x)
803         IF ( iom_use( 'intstry'  ) )   CALL iom_put( 'intstry'  ,  zdiag_intstry(:,:)   )   ! Internal force term in force balance (y)
804         IF ( iom_use( 'utau_oi'  ) )   CALL iom_put( 'utau_oi'  ,  zdiag_utau_oi(:,:)   )   ! Ocean stress term in force balance (x)
805         IF ( iom_use( 'vtau_oi'  ) )   CALL iom_put( 'vtau_oi'  ,  zdiag_vtau_oi(:,:)   )   ! Ocean stress term in force balance (y)
806         IF ( iom_use( 'xmtrpice' ) )   CALL iom_put( 'xmtrpice' ,  zdiag_xmtrp_ice(:,:) )   ! X-component of sea-ice mass transport (kg/s)
807         IF ( iom_use( 'ymtrpice' ) )   CALL iom_put( 'ymtrpice' ,  zdiag_ymtrp_ice(:,:) )   ! Y-component of sea-ice mass transport
808         IF ( iom_use( 'xmtrpsnw' ) )   CALL iom_put( 'xmtrpsnw' ,  zdiag_xmtrp_snw(:,:) )   ! X-component of snow mass transport (kg/s)
809         IF ( iom_use( 'ymtrpsnw' ) )   CALL iom_put( 'ymtrpsnw' ,  zdiag_ymtrp_snw(:,:) )   ! Y-component of snow mass transport
810         IF ( iom_use( 'xatrp'    ) )   CALL iom_put( 'xatrp'    ,  zdiag_xatrp(:,:)     )   ! X-component of ice area transport
811         IF ( iom_use( 'yatrp'    ) )   CALL iom_put( 'yatrp'    ,  zdiag_yatrp(:,:)     )   ! Y-component of ice area transport
812
813         DEALLOCATE( zdiag_sig1      , zdiag_sig2      , zdiag_dssh_dx   , zdiag_dssh_dy   ,  &
814            &        zdiag_corstrx   , zdiag_corstry   , zdiag_intstrx   , zdiag_intstry   ,  &
815            &        zdiag_utau_oi   , zdiag_vtau_oi   , zdiag_xmtrp_ice , zdiag_ymtrp_ice ,  &
816            &        zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp     , zdiag_yatrp     )
817
818      ENDIF
819      !
820   END SUBROUTINE ice_dyn_rhg_evp
821
822   SUBROUTINE rhg_evp_rst( cdrw, kt )
823      !!---------------------------------------------------------------------
824      !!                   ***  ROUTINE rhg_evp_rst  ***
825      !!                     
826      !! ** Purpose :   Read or write RHG file in restart file
827      !!
828      !! ** Method  :   use of IOM library
829      !!----------------------------------------------------------------------
830      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
831      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
832      !
833      INTEGER  ::   iter            ! local integer
834      INTEGER  ::   id1, id2, id3   ! local integers
835      !!----------------------------------------------------------------------
836      !
837      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
838         !                                   ! ---------------
839         IF( ln_rstart ) THEN                   !* Read the restart file
840            !
841            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
842            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
843            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
844            !
845            IF( MIN( id1, id2, id3 ) > 0 ) THEN      ! fields exist
846               CALL iom_get( numrir, jpdom_autoglo, 'stress1_i' , stress1_i  )
847               CALL iom_get( numrir, jpdom_autoglo, 'stress2_i' , stress2_i  )
848               CALL iom_get( numrir, jpdom_autoglo, 'stress12_i', stress12_i )
849            ELSE                                     ! start rheology from rest
850               IF(lwp) WRITE(numout,*) '   ==>>   previous run without rheology, set stresses to 0'
851               stress1_i (:,:) = 0._wp
852               stress2_i (:,:) = 0._wp
853               stress12_i(:,:) = 0._wp
854            ENDIF
855         ELSE                                   !* Start from rest
856            IF(lwp) WRITE(numout,*) '   ==>>   start from rest: set stresses to 0'
857            stress1_i (:,:) = 0._wp
858            stress2_i (:,:) = 0._wp
859            stress12_i(:,:) = 0._wp
860         ENDIF
861         !
862      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
863         !                                   ! -------------------
864         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
865         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
866         !
867         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i  )
868         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i  )
869         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
870         !
871      ENDIF
872      !
873   END SUBROUTINE rhg_evp_rst
874
875#else
876   !!----------------------------------------------------------------------
877   !!   Default option         Empty module           NO ESIM sea-ice model
878   !!----------------------------------------------------------------------
879#endif
880
881   !!==============================================================================
882END MODULE icedyn_rhg_evp
Note: See TracBrowser for help on using the repository browser.