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 NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap/src/ICE – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap/src/ICE/icedyn_rhg_evp.F90 @ 11053

Last change on this file since 11053 was 11053, checked in by davestorkey, 5 years ago

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps_rewrite_time_filterswap : Merge in latest changes from main branch and finish conversion of "h" variables. NB. This version still doesn't work!

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