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_eap.F90 in NEMO/branches/2019/dev_r11842_SI3-10_EAP/src/ICE – NEMO

source: NEMO/branches/2019/dev_r11842_SI3-10_EAP/src/ICE/icedyn_rhg_eap.F90 @ 13662

Last change on this file since 13662 was 13662, checked in by clem, 4 years ago

update to almost r4.0.4

File size: 106.2 KB
Line 
1MODULE icedyn_rhg_eap
2   !!======================================================================
3   !!                     ***  MODULE  icedyn_rhg_eap  ***
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   !!                 !  2019     (S. Rynders)  change into eap rheology from
16   !!                                           CICE code (Tsamados, Heorton)
17   !!----------------------------------------------------------------------
18#if defined key_si3
19   !!----------------------------------------------------------------------
20   !!   'key_si3'                                       SI3 sea-ice model
21   !!----------------------------------------------------------------------
22   !!   ice_dyn_rhg_eap : computes ice velocities from EVP rheology
23   !!   rhg_eap_rst     : read/write EVP fields in ice restart
24   !!----------------------------------------------------------------------
25   USE phycst         ! Physical constant
26   USE dom_oce        ! Ocean domain
27   USE sbc_oce , ONLY : ln_ice_embd, nn_fsbc, ssh_m
28   USE sbc_ice , ONLY : utau_ice, vtau_ice, snwice_mass, snwice_mass_b
29   USE ice            ! sea-ice: ice variables
30   USE icevar         ! ice_var_sshdyn
31   USE icedyn_rdgrft  ! sea-ice: ice strength
32   USE bdy_oce , ONLY : ln_bdy 
33   USE bdyice 
34#if defined key_agrif
35   USE agrif_ice_interp
36#endif
37   !
38   USE in_out_manager ! I/O manager
39   USE iom            ! I/O manager library
40   USE lib_mpp        ! MPP library
41   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
42   USE lbclnk         ! lateral boundary conditions (or mpp links)
43   USE prtctl         ! Print control
44
45   USE netcdf         ! NetCDF library for convergence test
46   IMPLICIT NONE
47   PRIVATE
48
49   PUBLIC   ice_dyn_rhg_eap   ! called by icedyn_rhg.F90
50   PUBLIC   rhg_eap_rst       ! called by icedyn_rhg.F90
51
52   REAL(wp), PARAMETER ::   pphi = 3.141592653589793_wp/12._wp    ! diamond shaped floe smaller angle (default phi = 30 deg)
53
54   ! look-up table for calculating structure tensor
55   INTEGER, PARAMETER ::   nx_yield = 41
56   INTEGER, PARAMETER ::   ny_yield = 41
57   INTEGER, PARAMETER ::   na_yield = 21
58
59   REAL(wp), DIMENSION(nx_yield, ny_yield, na_yield) ::   s11r, s12r, s22r, s11s, s12s, s22s
60
61   !! * Substitutions
62#  include "vectopt_loop_substitute.h90"
63
64   !! for convergence tests
65   INTEGER ::   ncvgid   ! netcdf file id
66   INTEGER ::   nvarid   ! netcdf variable id
67   REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zmsk00, zmsk15
68   !!----------------------------------------------------------------------
69   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
70   !! $Id: icedyn_rhg_eap.F90 11536 2019-09-11 13:54:18Z smasson $
71   !! Software governed by the CeCILL license (see ./LICENSE)
72   !!----------------------------------------------------------------------
73CONTAINS
74
75   SUBROUTINE ice_dyn_rhg_eap( kt, pstress1_i, pstress2_i, pstress12_i, pshear_i, pdivu_i, pdelta_i, paniso_11, paniso_12, prdg_conv )
76      !!-------------------------------------------------------------------
77      !!                 ***  SUBROUTINE ice_dyn_rhg_eap  ***
78      !!                             EAP-C-grid
79      !!
80      !! ** purpose : determines sea ice drift from wind stress, ice-ocean
81      !!  stress and sea-surface slope. Ice-ice interaction is described by
82      !!  a non-linear elasto-anisotropic-plastic (EAP) law including shear
83      !!  strength and a bulk rheology . 
84      !!
85      !!  The points in the C-grid look like this, dear reader
86      !!
87      !!                              (ji,jj)
88      !!                                 |
89      !!                                 |
90      !!                      (ji-1,jj)  |  (ji,jj)
91      !!                             ---------   
92      !!                            |         |
93      !!                            | (ji,jj) |------(ji,jj)
94      !!                            |         |
95      !!                             ---------   
96      !!                     (ji-1,jj-1)     (ji,jj-1)
97      !!
98      !! ** Inputs  : - wind forcing (stress), oceanic currents
99      !!                ice total volume (vt_i) per unit area
100      !!                snow total volume (vt_s) per unit area
101      !!
102      !! ** Action  : - compute u_ice, v_ice : the components of the
103      !!                sea-ice velocity vector
104      !!              - compute delta_i, shear_i, divu_i, which are inputs
105      !!                of the ice thickness distribution
106      !!
107      !! ** Steps   : 0) compute mask at F point
108      !!              1) Compute ice snow mass, ice strength
109      !!              2) Compute wind, oceanic stresses, mass terms and
110      !!                 coriolis terms of the momentum equation
111      !!              3) Solve the momentum equation (iterative procedure)
112      !!              4) Recompute delta, shear and divergence
113      !!                 (which are inputs of the ITD) & store stress
114      !!                 for the next time step
115      !!              5) Diagnostics including charge ellipse
116      !!
117      !! ** Notes   : There is the possibility to use aEVP from the nice work of Kimmritz et al. (2016 & 2017)
118      !!              by setting up ln_aEVP=T (i.e. changing alpha and beta parameters).
119      !!              This is an upgraded version of mEVP from Bouillon et al. 2013
120      !!              (i.e. more stable and better convergence)
121      !!
122      !! References : Hunke and Dukowicz, JPO97
123      !!              Bouillon et al., Ocean Modelling 2009
124      !!              Bouillon et al., Ocean Modelling 2013
125      !!              Kimmritz et al., Ocean Modelling 2016 & 2017
126      !!-------------------------------------------------------------------
127      INTEGER                 , INTENT(in   ) ::   kt                                    ! time step
128      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pstress1_i, pstress2_i, pstress12_i   !
129      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pshear_i  , pdivu_i   , pdelta_i      !
130      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   paniso_11 , paniso_12                 ! structure tensor components
131      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   prdg_conv                             ! for ridging
132      !!
133      INTEGER ::   ji, jj       ! dummy loop indices
134      INTEGER ::   jter         ! local integers
135      !
136      REAL(wp) ::   zrhoco                                              ! rau0 * rn_cio
137      REAL(wp) ::   zdtevp, z1_dtevp                                    ! time step for subcycling
138      REAL(wp) ::   ecc2, z1_ecc2                                       ! square of yield ellipse eccenticity
139      REAL(wp) ::   zalph1, z1_alph1, zalph2, z1_alph2                  ! alpha coef from Bouillon 2009 or Kimmritz 2017
140      REAl(wp) ::   zbetau, zbetav
141      REAL(wp) ::   zm1, zm2, zm3, zmassU, zmassV, zvU, zvV             ! ice/snow mass and volume
142      REAL(wp) ::   zds2, zdt, zdt2, zdiv, zdiv2, zdsT                  ! temporary scalars
143      REAL(wp) ::   zTauO, zTauB, zRHS, zvel                            ! temporary scalars
144      REAL(wp) ::   zkt                                                 ! isotropic tensile strength for landfast ice
145      REAL(wp) ::   zvCr                                                ! critical ice volume above which ice is landfast
146      !
147      REAL(wp) ::   zintb, zintn                                        ! dummy argument
148      REAL(wp) ::   zfac_x, zfac_y
149      REAL(wp) ::   zshear, zdum1, zdum2
150      REAL(wp) ::   zstressptmp, zstressmtmp, zstress12tmpF             ! anisotropic stress tensor components
151      REAL(wp) ::   zalphar, zalphas                                    ! for mechanical redistribution
152      REAL(wp) ::   zmresult11, zmresult12, z1dtevpkth, zp5kth, z1_dtevp_A  ! for structure tensor evolution
153      !
154      REAL(wp), DIMENSION(jpi,jpj) ::   zstress12tmp                    ! anisotropic stress tensor component for regridding
155      REAL(wp), DIMENSION(jpi,jpj) ::   zyield11, zyield22, zyield12    ! yield surface tensor for history
156      REAL(wp), DIMENSION(jpi,jpj) ::   zdelta, zp_delt                 ! delta and P/delta at T points
157      REAL(wp), DIMENSION(jpi,jpj) ::   zten_i                          ! tension
158      REAL(wp), DIMENSION(jpi,jpj) ::   zbeta                           ! beta coef from Kimmritz 2017
159      !
160      REAL(wp), DIMENSION(jpi,jpj) ::   zdt_m                           ! (dt / ice-snow_mass) on T points
161      REAL(wp), DIMENSION(jpi,jpj) ::   zaU  , zaV                      ! ice fraction on U/V points
162      REAL(wp), DIMENSION(jpi,jpj) ::   zmU_t, zmV_t                    ! (ice-snow_mass / dt) on U/V points
163      REAL(wp), DIMENSION(jpi,jpj) ::   zmf                             ! coriolis parameter at T points
164      REAL(wp), DIMENSION(jpi,jpj) ::   v_oceU, u_oceV, v_iceU, u_iceV  ! ocean/ice u/v component on V/U points                           
165      !
166      REAL(wp), DIMENSION(jpi,jpj) ::   zds                             ! shear
167      REAL(wp), DIMENSION(jpi,jpj) ::   zs1, zs2, zs12                  ! stress tensor components
168      REAL(wp), DIMENSION(jpi,jpj) ::   zsshdyn                         ! array used for the calculation of ice surface slope:
169      !                                                                 !    ocean surface (ssh_m) if ice is not embedded
170      !                                                                 !    ice bottom surface if ice is embedded   
171      REAL(wp), DIMENSION(jpi,jpj) ::   zfU  , zfV                      ! internal stresses
172      REAL(wp), DIMENSION(jpi,jpj) ::   zspgU, zspgV                    ! surface pressure gradient at U/V points
173      REAL(wp), DIMENSION(jpi,jpj) ::   zCorU, zCorV                    ! Coriolis stress array
174      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_ai, ztauy_ai              ! ice-atm. stress at U-V points
175      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_oi, ztauy_oi              ! ice-ocean stress at U-V points
176      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_bi, ztauy_bi              ! ice-OceanBottom stress at U-V points (landfast)
177      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_base, ztauy_base          ! ice-bottom stress at U-V points (landfast)
178      !
179      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk01x, zmsk01y                ! dummy arrays
180      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk00x, zmsk00y                ! mask for ice presence
181      REAL(wp), DIMENSION(jpi,jpj) ::   zfmask                          ! mask at F points for the ice
182
183      REAL(wp), PARAMETER          ::   zepsi  = 1.0e-20_wp             ! tolerance parameter
184      REAL(wp), PARAMETER          ::   zmmin  = 1._wp                  ! ice mass (kg/m2)  below which ice velocity becomes very small
185      REAL(wp), PARAMETER          ::   zamin  = 0.001_wp               ! ice concentration below which ice velocity becomes very small
186      !! --- check convergence
187      REAL(wp), DIMENSION(jpi,jpj) ::   zu_ice, zv_ice
188      !! --- diags
189      REAL(wp) ::   zsig1, zsig2, zsig12, zfac, z1_strength
190      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zsig_I, zsig_II, zsig1_p, zsig2_p         
191      !! --- SIMIP diags
192      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s)
193      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s)
194      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s)
195      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s)
196      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xatrp     ! X-component of area transport (m2/s)
197      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_yatrp     ! Y-component of area transport (m2/s)     
198      !!-------------------------------------------------------------------
199
200      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_dyn_rhg_eap: EAP sea-ice rheology'
201      !
202      ! for diagnostics and convergence tests
203      ALLOCATE( zmsk00(jpi,jpj), zmsk15(jpi,jpj) )
204      DO jj = 1, jpj
205         DO ji = 1, jpi
206            zmsk00(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06  ) ) ! 1 if ice    , 0 if no ice
207            zmsk15(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - 0.15_wp ) ) ! 1 if 15% ice, 0 if less
208         END DO
209      END DO
210      !
211!!gm for Clem:  OPTIMIZATION:  I think zfmask can be computed one for all at the initialization....
212      !------------------------------------------------------------------------------!
213      ! 0) mask at F points for the ice
214      !------------------------------------------------------------------------------!
215      ! ocean/land mask
216      DO jj = 1, jpjm1
217         DO ji = 1, jpim1      ! NO vector opt.
218            zfmask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
219         END DO
220      END DO
221      CALL lbc_lnk( 'icedyn_rhg_eap', zfmask, 'F', 1._wp )
222
223      ! Lateral boundary conditions on velocity (modify zfmask)
224      DO jj = 2, jpjm1
225         DO ji = fs_2, fs_jpim1   ! vector opt.
226            IF( zfmask(ji,jj) == 0._wp ) THEN
227               zfmask(ji,jj) = rn_ishlat * MIN( 1._wp , MAX( umask(ji,jj,1), umask(ji,jj+1,1), &
228                  &                                          vmask(ji,jj,1), vmask(ji+1,jj,1) ) )
229            ENDIF
230         END DO
231      END DO
232      DO jj = 2, jpjm1
233         IF( zfmask(1,jj) == 0._wp ) THEN
234            zfmask(1  ,jj) = rn_ishlat * MIN( 1._wp , MAX( vmask(2,jj,1), umask(1,jj+1,1), umask(1,jj,1) ) )
235         ENDIF
236         IF( zfmask(jpi,jj) == 0._wp ) THEN
237            zfmask(jpi,jj) = rn_ishlat * MIN( 1._wp , MAX( umask(jpi,jj+1,1), vmask(jpim1,jj,1), umask(jpi,jj-1,1) ) )
238         ENDIF
239      END DO
240      DO ji = 2, jpim1
241         IF( zfmask(ji,1) == 0._wp ) THEN
242            zfmask(ji, 1 ) = rn_ishlat * MIN( 1._wp , MAX( vmask(ji+1,1,1), umask(ji,2,1), vmask(ji,1,1) ) )
243         ENDIF
244         IF( zfmask(ji,jpj) == 0._wp ) THEN
245            zfmask(ji,jpj) = rn_ishlat * MIN( 1._wp , MAX( vmask(ji+1,jpj,1), vmask(ji-1,jpj,1), umask(ji,jpjm1,1) ) )
246         ENDIF
247      END DO
248      CALL lbc_lnk( 'icedyn_rhg_eap', zfmask, 'F', 1._wp )
249
250      !------------------------------------------------------------------------------!
251      ! 1) define some variables and initialize arrays
252      !------------------------------------------------------------------------------!
253      zrhoco = rau0 * rn_cio 
254
255      ! ecc2: square of yield ellipse eccenticrity
256      ecc2    = rn_ecc * rn_ecc
257      z1_ecc2 = 1._wp / ecc2
258
259      ! alpha parameters (Bouillon 2009)
260      IF( .NOT. ln_aEVP ) THEN
261         zdtevp   = rdt_ice / REAL( nn_nevp )
262         zalph1 =   2._wp * rn_relast * REAL( nn_nevp )
263         zalph2 = zalph1 * z1_ecc2
264
265         z1_alph1 = 1._wp / ( zalph1 + 1._wp )
266         z1_alph2 = 1._wp / ( zalph2 + 1._wp )
267      ELSE
268         zdtevp   = rdt_ice
269         ! zalpha parameters set later on adaptatively
270      ENDIF
271      z1_dtevp = 1._wp / zdtevp
272         
273      ! Initialise stress tensor
274      zs1 (:,:) = pstress1_i (:,:) 
275      zs2 (:,:) = pstress2_i (:,:)
276      zs12(:,:) = pstress12_i(:,:)
277
278      ! constants for structure tensor
279      z1_dtevp_A = z1_dtevp/10.0_wp
280      z1dtevpkth = 1._wp / (z1_dtevp_A + 0.00002_wp)
281      zp5kth = 0.5_wp * 0.00002_wp
282
283      ! Ice strength
284      CALL ice_strength
285
286      ! landfast param from Lemieux(2016): add isotropic tensile strength (following Konig Beatty and Holland, 2010)
287      IF( ln_landfast_L16 ) THEN   ;   zkt = rn_lf_tensile
288      ELSE                         ;   zkt = 0._wp
289      ENDIF
290      !
291      !------------------------------------------------------------------------------!
292      ! 2) Wind / ocean stress, mass terms, coriolis terms
293      !------------------------------------------------------------------------------!
294      ! sea surface height
295      !    embedded sea ice: compute representative ice top surface
296      !    non-embedded sea ice: use ocean surface for slope calculation
297      zsshdyn(:,:) = ice_var_sshdyn( ssh_m, snwice_mass, snwice_mass_b)
298
299      DO jj = 2, jpjm1
300         DO ji = fs_2, fs_jpim1
301
302            ! ice fraction at U-V points
303            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)
304            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)
305
306            ! Ice/snow mass at U-V points
307            zm1 = ( rhos * vt_s(ji  ,jj  ) + rhoi * vt_i(ji  ,jj  ) )
308            zm2 = ( rhos * vt_s(ji+1,jj  ) + rhoi * vt_i(ji+1,jj  ) )
309            zm3 = ( rhos * vt_s(ji  ,jj+1) + rhoi * vt_i(ji  ,jj+1) )
310            zmassU = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm2 * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1)
311            zmassV = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm3 * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1)
312
313            ! Ocean currents at U-V points
314            v_oceU(ji,jj)   = 0.25_wp * ( v_oce(ji,jj) + v_oce(ji,jj-1) + v_oce(ji+1,jj) + v_oce(ji+1,jj-1) ) * umask(ji,jj,1)
315            u_oceV(ji,jj)   = 0.25_wp * ( u_oce(ji,jj) + u_oce(ji-1,jj) + u_oce(ji,jj+1) + u_oce(ji-1,jj+1) ) * vmask(ji,jj,1)
316
317            ! Coriolis at T points (m*f)
318            zmf(ji,jj)      = zm1 * ff_t(ji,jj)
319
320            ! dt/m at T points (for alpha and beta coefficients)
321            zdt_m(ji,jj)    = zdtevp / MAX( zm1, zmmin )
322           
323            ! m/dt
324            zmU_t(ji,jj)    = zmassU * z1_dtevp
325            zmV_t(ji,jj)    = zmassV * z1_dtevp
326           
327            ! Drag ice-atm.
328            ztaux_ai(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj)
329            ztauy_ai(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj)
330
331            ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points
332            zspgU(ji,jj)    = - zmassU * grav * ( zsshdyn(ji+1,jj) - zsshdyn(ji,jj) ) * r1_e1u(ji,jj)
333            zspgV(ji,jj)    = - zmassV * grav * ( zsshdyn(ji,jj+1) - zsshdyn(ji,jj) ) * r1_e2v(ji,jj)
334
335            ! masks
336            zmsk00x(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) )  ! 0 if no ice
337            zmsk00y(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) )  ! 0 if no ice
338
339            ! switches
340            IF( zmassU <= zmmin .AND. zaU(ji,jj) <= zamin ) THEN   ;   zmsk01x(ji,jj) = 0._wp
341            ELSE                                                   ;   zmsk01x(ji,jj) = 1._wp   ;   ENDIF
342            IF( zmassV <= zmmin .AND. zaV(ji,jj) <= zamin ) THEN   ;   zmsk01y(ji,jj) = 0._wp
343            ELSE                                                   ;   zmsk01y(ji,jj) = 1._wp   ;   ENDIF
344
345         END DO
346      END DO
347      CALL lbc_lnk_multi( 'icedyn_rhg_eap', zmf, 'T', 1., zdt_m, 'T', 1. )
348      !
349      !                                  !== Landfast ice parameterization ==!
350      !
351      IF( ln_landfast_L16 ) THEN         !-- Lemieux 2016
352         DO jj = 2, jpjm1
353            DO ji = fs_2, fs_jpim1
354               ! ice thickness at U-V points
355               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)
356               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)
357               ! ice-bottom stress at U points
358               zvCr = zaU(ji,jj) * rn_lf_depfra * hu_n(ji,jj)
359               ztaux_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvU - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaU(ji,jj) ) )
360               ! ice-bottom stress at V points
361               zvCr = zaV(ji,jj) * rn_lf_depfra * hv_n(ji,jj)
362               ztauy_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvV - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaV(ji,jj) ) )
363               ! ice_bottom stress at T points
364               zvCr = at_i(ji,jj) * rn_lf_depfra * ht_n(ji,jj)
365               tau_icebfr(ji,jj) = - rn_lf_bfr * MAX( 0._wp, vt_i(ji,jj) - zvCr ) * EXP( -rn_crhg * ( 1._wp - at_i(ji,jj) ) )
366            END DO
367         END DO
368         CALL lbc_lnk( 'icedyn_rhg_eap', tau_icebfr(:,:), 'T', 1. )
369         !
370      ELSE                               !-- no landfast
371         DO jj = 2, jpjm1
372            DO ji = fs_2, fs_jpim1
373               ztaux_base(ji,jj) = 0._wp
374               ztauy_base(ji,jj) = 0._wp
375            END DO
376         END DO
377      ENDIF
378
379      !------------------------------------------------------------------------------!
380      ! 3) Solution of the momentum equation, iterative procedure
381      !------------------------------------------------------------------------------!
382      !
383      !                                               ! ==================== !
384      DO jter = 1 , nn_nevp                           !    loop over jter    !
385         !                                            ! ==================== !       
386         l_full_nf_update = jter == nn_nevp   ! false: disable full North fold update (performances) for iter = 1 to nn_nevp-1
387         !
388         ! convergence test
389         IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2  ) THEN
390            DO jj = 1, jpj
391               DO ji = 1, jpi
392                  zu_ice(ji,jj) = u_ice(ji,jj) * umask(ji,jj,1) ! velocity at previous time step
393                  zv_ice(ji,jj) = v_ice(ji,jj) * vmask(ji,jj,1)
394               END DO
395            END DO
396         ENDIF
397
398         ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- !
399         DO jj = 1, jpjm1
400            DO ji = 1, jpim1
401
402               ! shear at F points
403               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)   &
404                  &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
405                  &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
406
407            END DO
408         END DO
409         CALL lbc_lnk( 'icedyn_rhg_eap', zds, 'F', 1._wp )
410
411         DO jj = 2, jpjm1
412            DO ji = 2, jpim1 ! no vector loop
413
414               ! shear**2 at T points (doc eq. A16)
415               zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
416                  &   + 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)  &
417                  &   ) * 0.25_wp * r1_e1e2t(ji,jj)
418             
419               ! divergence at T points
420               zdiv  = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
421                  &    + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
422                  &    ) * r1_e1e2t(ji,jj)
423               zdiv2 = zdiv * zdiv
424               
425               ! tension at T points
426               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)   &
427                  &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
428                  &   ) * r1_e1e2t(ji,jj)
429               zdt2 = zdt * zdt
430
431               ! delta at T points
432               zdelta(ji,jj) = SQRT( zdiv2 + ( zdt2 + zds2 ) * z1_ecc2 ) 
433
434            END DO
435         END DO
436         CALL lbc_lnk( 'icedyn_rhg_eap', zdelta, 'T', 1._wp )
437               
438         ! P/delta at T points
439         DO jj = 1, jpj
440            DO ji = 1, jpi
441               zp_delt(ji,jj) = strength(ji,jj) / ( zdelta(ji,jj) + rn_creepl )
442            END DO
443         END DO
444
445         DO jj = 2, jpj    ! loop ends at jpi,jpj so that no lbc_lnk are needed for zs1 and zs2
446            DO ji = 2, jpi ! no vector loop
447
448                ! shear at T points
449               zdsT = ( zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
450                  &   + zds(ji,jj-1) * e1e2f(ji,jj-1) + zds(ji-1,jj-1) * e1e2f(ji-1,jj-1)  &
451                  &   ) * 0.25_wp * r1_e1e2t(ji,jj)
452               !WRITE(numout,*) 'shear output', ji, jj, zdsT
453               
454              ! divergence at T points (duplication to avoid communications)
455               zdiv  = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
456                  &    + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
457                  &    ) * r1_e1e2t(ji,jj)
458               
459               ! tension at T points (duplication to avoid communications)
460               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)   &
461                  &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
462                  &   ) * r1_e1e2t(ji,jj)
463
464               ! --- anisotropic stress calculation --- !
465               CALL update_stress_rdg (jter, nn_nevp, zdiv, zdt, zdsT, paniso_11(ji,jj), paniso_12(ji,jj), &
466                                       zstressptmp, zstressmtmp, zstress12tmp(ji,jj), strength(ji,jj), zalphar, zalphas)
467
468               ! structure tensor update
469!!$               IF (mod(jter,10) == 0) THEN
470                  CALL calc_ffrac(zstressptmp, zstressmtmp, zstress12tmp(ji,jj), paniso_11(ji,jj), paniso_12(ji,jj), zmresult11,  zmresult12)
471
472!!$                  paniso_11(ji,jj) = (paniso_11(ji,jj)*z1_dtevp_A  + zp5kth + zmresult11) * z1dtevpkth ! implicit
473!!$                  paniso_12(ji,jj) = (paniso_12(ji,jj)*z1_dtevp_A  + zmresult12) * z1dtevpkth ! implicit
474                  paniso_11(ji,jj) = (paniso_11(ji,jj)  + 0.5*2.e-5*zdtevp + zmresult11*zdtevp) / (1. + 2.e-5*zdtevp) ! implicit
475                  paniso_12(ji,jj) = (paniso_12(ji,jj)                     + zmresult12*zdtevp) / (1. + 2.e-5*zdtevp) ! implicit
476!!$               ENDIF
477
478               IF (jter == nn_nevp) THEN
479                  zyield11(ji,jj) = 0.5_wp * (zstressptmp + zstressmtmp)
480                  zyield22(ji,jj) = 0.5_wp * (zstressptmp - zstressmtmp)
481                  zyield12(ji,jj) = zstress12tmp(ji,jj)
482                  prdg_conv(ji,jj) = -min( zalphar, 0._wp)   
483               ENDIF
484
485               ! alpha for aEVP
486               !   gamma = 0.5*P/(delta+creepl) * (c*pi)**2/Area * dt/m
487               !   alpha = beta = sqrt(4*gamma)
488               IF( ln_aEVP ) THEN
489                  zalph1   = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
490                  z1_alph1 = 1._wp / ( zalph1 + 1._wp )
491                  zalph2   = zalph1
492                  z1_alph2 = z1_alph1
493                  ! explicit:
494                  ! z1_alph1 = 1._wp / zalph1
495                  ! z1_alph2 = 1._wp / zalph1
496                  ! zalph1 = zalph1 - 1._wp
497                  ! zalph2 = zalph1
498               ENDIF
499
500               ! stress at T points (zkt/=0 if landfast)
501               !zs1(ji,jj) = ( zs1(ji,jj)*zalph1 + zp_delt(ji,jj) * ( zdiv*(1._wp + zkt) - zdelta(ji,jj)*(1._wp - zkt) ) ) * z1_alph1
502               !zs2(ji,jj) = ( zs2(ji,jj)*zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 * (1._wp + zkt) ) ) * z1_alph2
503!!$               zs1(ji,jj) = ( zs1(ji,jj) + zstressptmp * zalph1  ) * z1_alph1
504!!$               zs2(ji,jj) = ( zs2(ji,jj) + zstressmtmp * zalph1  ) * z1_alph1
505               zs1(ji,jj) = ( zs1(ji,jj) * zalph1 + zstressptmp ) * z1_alph1
506               zs2(ji,jj) = ( zs2(ji,jj) * zalph1 + zstressmtmp ) * z1_alph1
507               !zs1(ji,jj) = ( stressptmp * zs1(ji,jj) + zalph1  ) * z1_alph1
508               !zs2(ji,jj) = ( stressmtmp * zs2(ji,jj) + zalph1  ) * z1_alph1
509               !WRITE(numout,*) 'stress output', ji, jj, zs1(ji,jj)
510            END DO
511         END DO
512         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zstress12tmp, 'T', 1. , paniso_11, 'T', 1. , paniso_12, 'T', 1.)
513
514        ! Save beta at T-points for further computations
515         IF( ln_aEVP ) THEN
516            DO jj = 1, jpj
517               DO ji = 1, jpi
518                  zbeta(ji,jj) = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
519               END DO
520            END DO
521         ENDIF
522
523         DO jj = 1, jpjm1
524            DO ji = 1, jpim1
525               ! stress12tmp at F points
526               zstress12tmpF = ( zstress12tmp(ji,jj+1) * e1e2t(ji,jj+1) + zstress12tmp(ji+1,jj+1) * e1e2t(ji+1,jj+1)  &
527                  &            + zstress12tmp(ji,jj  ) * e1e2t(ji,jj  ) + zstress12tmp(ji+1,jj  ) * e1e2t(ji+1,jj  )  &
528                  &            ) * 0.25_wp * r1_e1e2f(ji,jj)
529
530               ! alpha for aEVP
531               IF( ln_aEVP ) THEN
532                  zalph2   = MAX( zbeta(ji,jj), zbeta(ji+1,jj), zbeta(ji,jj+1), zbeta(ji+1,jj+1) )
533                  z1_alph2 = 1._wp / ( zalph2 + 1._wp )
534                  ! explicit:
535                  ! z1_alph2 = 1._wp / zalph2
536                  ! zalph2 = zalph2 - 1._wp
537               ENDIF
538               
539               ! stress at F points (zkt/=0 if landfast)
540               !zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 * (1._wp + zkt) ) * 0.5_wp ) * z1_alph2
541!!$               zs12(ji,jj) = ( zs12(ji,jj) + zstress12tmpF * zalph1  ) * z1_alph1
542               zs12(ji,jj) = ( zs12(ji,jj) * zalph1 + zstress12tmpF ) * z1_alph1
543               !zs12(ji,jj) = ( stress12tmpF * zs12(ji,jj) + zalph1  ) * z1_alph1
544
545            END DO
546         END DO
547         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zs1, 'T', 1., zs2, 'T', 1., zs12, 'F', 1. )
548
549         ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- !
550         DO jj = 2, jpjm1
551            DO ji = fs_2, fs_jpim1               
552               !                   !--- U points
553               zfU(ji,jj) = 0.5_wp * ( ( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj)                                             &
554                  &                  + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj)    &
555                  &                    ) * r1_e2u(ji,jj)                                                                      &
556                  &                  + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1)  &
557                  &                    ) * 2._wp * r1_e1u(ji,jj)                                                              &
558                  &                  ) * r1_e1e2u(ji,jj)
559               !
560               !                !--- V points
561               zfV(ji,jj) = 0.5_wp * ( ( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj)                                             &
562                  &                  - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj)    &
563                  &                    ) * r1_e1v(ji,jj)                                                                      &
564                  &                  + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj)  &
565                  &                    ) * 2._wp * r1_e2v(ji,jj)                                                              &
566                  &                  ) * r1_e1e2v(ji,jj)
567               !
568               !                !--- ice currents at U-V point
569               v_iceU(ji,jj) = 0.25_wp * ( v_ice(ji,jj) + v_ice(ji,jj-1) + v_ice(ji+1,jj) + v_ice(ji+1,jj-1) ) * umask(ji,jj,1)
570               u_iceV(ji,jj) = 0.25_wp * ( u_ice(ji,jj) + u_ice(ji-1,jj) + u_ice(ji,jj+1) + u_ice(ji-1,jj+1) ) * vmask(ji,jj,1)
571               !
572            END DO
573         END DO
574         !
575         ! --- Computation of ice velocity --- !
576         !  Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta vary as in Kimmritz 2016 & 2017
577         !  Bouillon et al. 2009 (eq 34-35) => stable
578         IF( MOD(jter,2) == 0 ) THEN ! even iterations
579            !
580            DO jj = 2, jpjm1
581               DO ji = fs_2, fs_jpim1
582                  !                 !--- tau_io/(v_oce - v_ice)
583                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
584                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
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  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
590                  zTauB = ztauy_base(ji,jj) / zvel
591                  !                 !--- OceanBottom-to-Ice stress
592                  ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
593                  !
594                  !                 !--- Coriolis at V-points (energy conserving formulation)
595                  zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
596                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
597                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
598                  !
599                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
600                  zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
601                  !
602                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
603                  !                                         1 = sliding friction : TauB < RHS
604                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
605                  !
606                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
607                     zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) )
608                     v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) )         & ! previous velocity
609                        &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
610                        &                                    ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
611                        &            + ( 1._wp - rswitch ) * (  v_ice_b(ji,jj)                                                   & 
612                        &                                     + v_ice  (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
613                        &                                    ) / ( zbetav + 1._wp )                                              &
614                        &             ) * zmsk01y(ji,jj) + v_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01y(ji,jj) )                   & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
615                        &           )   * zmsk00y(ji,jj)
616                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
617                     v_ice(ji,jj) = ( (         rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                                       & ! previous velocity
618                        &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
619                        &                                    ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
620                        &            + ( 1._wp - rswitch ) *   v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
621                        &              ) * zmsk01y(ji,jj) + v_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01y(ji,jj) )                  & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
622                        &            )   * zmsk00y(ji,jj)
623                  ENDIF
624               END DO
625            END DO
626            CALL lbc_lnk( 'icedyn_rhg_eap', v_ice, 'V', -1. )
627            !
628#if defined key_agrif
629!!            CALL agrif_interp_ice( 'V', jter, nn_nevp )
630            CALL agrif_interp_ice( 'V' )
631#endif
632            IF( ln_bdy )   CALL bdy_ice_dyn( 'V' )
633            !
634            DO jj = 2, jpjm1
635               DO ji = fs_2, fs_jpim1         
636                  !                 !--- tau_io/(u_oce - u_ice)
637                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
638                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
639                  !                 !--- Ocean-to-Ice stress
640                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
641                  !
642                  !                 !--- tau_bottom/u_ice
643                  zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
644                  zTauB = ztaux_base(ji,jj) / zvel
645                  !                 !--- OceanBottom-to-Ice stress
646                  ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
647                  !
648                  !                 !--- Coriolis at U-points (energy conserving formulation)
649                  zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
650                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
651                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
652                  !
653                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
654                  zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
655                  !
656                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
657                  !                                         1 = sliding friction : TauB < RHS
658                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
659                  !
660                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
661                     zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) )
662                     u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) )         & ! previous velocity
663                        &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
664                        &                                    ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
665                        &            + ( 1._wp - rswitch ) * (  u_ice_b(ji,jj)                                                   &
666                        &                                     + u_ice  (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
667                        &                                    ) / ( zbetau + 1._wp )                                              &
668                        &             ) * zmsk01x(ji,jj) + u_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01x(ji,jj) )                   & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
669                        &           )   * zmsk00x(ji,jj)
670                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
671                     u_ice(ji,jj) = ( (         rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                                       & ! previous velocity
672                        &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
673                        &                                    ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
674                        &            + ( 1._wp - rswitch ) *   u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
675                        &              ) * zmsk01x(ji,jj) + u_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01x(ji,jj) )                  & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
676                        &            )   * zmsk00x(ji,jj)
677                  ENDIF
678               END DO
679            END DO
680            CALL lbc_lnk( 'icedyn_rhg_eap', u_ice, 'U', -1. )
681            !
682#if defined key_agrif
683!!            CALL agrif_interp_ice( 'U', jter, nn_nevp )
684            CALL agrif_interp_ice( 'U' )
685#endif
686            IF( ln_bdy )   CALL bdy_ice_dyn( 'U' )
687            !
688         ELSE ! odd iterations
689            !
690            DO jj = 2, jpjm1
691               DO ji = fs_2, fs_jpim1
692                  !                 !--- tau_io/(u_oce - u_ice)
693                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
694                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
695                  !                 !--- Ocean-to-Ice stress
696                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
697                  !
698                  !                 !--- tau_bottom/u_ice
699                  zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
700                  zTauB = ztaux_base(ji,jj) / zvel
701                  !                 !--- OceanBottom-to-Ice stress
702                  ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
703                  !
704                  !                 !--- Coriolis at U-points (energy conserving formulation)
705                  zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
706                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
707                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
708                  !
709                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
710                  zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
711                  !
712                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
713                  !                                         1 = sliding friction : TauB < RHS
714                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
715                  !
716                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
717                     zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) )
718                     u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) )         & ! previous velocity
719                        &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
720                        &                                    ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
721                        &            + ( 1._wp - rswitch ) * (  u_ice_b(ji,jj)                                                   &
722                        &                                     + u_ice  (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
723                        &                                    ) / ( zbetau + 1._wp )                                              &
724                        &             ) * zmsk01x(ji,jj) + u_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01x(ji,jj) )                   & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
725                        &           )   * zmsk00x(ji,jj)
726                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
727                     u_ice(ji,jj) = ( (         rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                                       & ! previous velocity
728                        &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
729                        &                                    ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
730                        &            + ( 1._wp - rswitch ) *   u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
731                        &              ) * zmsk01x(ji,jj) + u_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01x(ji,jj) )                  & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
732                        &            )   * zmsk00x(ji,jj)
733                  ENDIF
734               END DO
735            END DO
736            CALL lbc_lnk( 'icedyn_rhg_eap', u_ice, 'U', -1. )
737            !
738#if defined key_agrif
739!!            CALL agrif_interp_ice( 'U', jter, nn_nevp )
740            CALL agrif_interp_ice( 'U' )
741#endif
742            IF( ln_bdy )   CALL bdy_ice_dyn( 'U' )
743            !
744            DO jj = 2, jpjm1
745               DO ji = fs_2, fs_jpim1
746                  !                 !--- tau_io/(v_oce - v_ice)
747                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
748                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
749                  !                 !--- Ocean-to-Ice stress
750                  ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
751                  !
752                  !                 !--- tau_bottom/v_ice
753                  zvel  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
754                  zTauB = ztauy_base(ji,jj) / zvel
755                  !                 !--- OceanBottom-to-Ice stress
756                  ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
757                  !
758                  !                 !--- Coriolis at v-points (energy conserving formulation)
759                  zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
760                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
761                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
762                  !
763                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
764                  zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
765                  !
766                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
767                  !                                         1 = sliding friction : TauB < RHS
768                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
769                  !
770                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
771                     zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) )
772                     v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) )         & ! previous velocity
773                        &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
774                        &                                    ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
775                        &            + ( 1._wp - rswitch ) * (  v_ice_b(ji,jj)                                                   &
776                        &                                     + v_ice  (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
777                        &                                    ) / ( zbetav + 1._wp )                                              & 
778                        &             ) * zmsk01y(ji,jj) + v_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01y(ji,jj) )                   & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
779                        &           )   * zmsk00y(ji,jj)
780                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
781                     v_ice(ji,jj) = ( (         rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                                       & ! previous velocity
782                        &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
783                        &                                    ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
784                        &            + ( 1._wp - rswitch ) *   v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
785                        &              ) * zmsk01y(ji,jj) + v_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01y(ji,jj) )                  & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
786                        &            )   * zmsk00y(ji,jj)
787                  ENDIF
788               END DO
789            END DO
790            CALL lbc_lnk( 'icedyn_rhg_eap', v_ice, 'V', -1. )
791            !
792#if defined key_agrif
793!!            CALL agrif_interp_ice( 'V', jter, nn_nevp )
794            CALL agrif_interp_ice( 'V' )
795#endif
796            IF( ln_bdy )   CALL bdy_ice_dyn( 'V' )
797            !
798         ENDIF
799
800         ! convergence test
801         IF( nn_rhg_chkcvg == 2 )   CALL rhg_cvg( kt, jter, nn_nevp, u_ice, v_ice, zu_ice, zv_ice )
802         !
803         !                                                ! ==================== !
804      END DO                                              !  end loop over jter  !
805      !                                                   ! ==================== !
806      IF( ln_aEVP )   CALL iom_put( 'beta_evp' , zbeta )
807      !
808      CALL lbc_lnk( 'icedyn_rhg_eap', prdg_conv, 'T', 1. )  ! only need this in ridging module after rheology completed
809      !
810      !------------------------------------------------------------------------------!
811      ! 4) Recompute delta, shear and div (inputs for mechanical redistribution)
812      !------------------------------------------------------------------------------!
813      DO jj = 1, jpjm1
814         DO ji = 1, jpim1
815
816            ! shear at F points
817            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)   &
818               &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
819               &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
820
821         END DO
822      END DO           
823     
824      DO jj = 2, jpjm1
825         DO ji = 2, jpim1 ! no vector loop
826           
827            ! tension**2 at T points
828            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)   &
829               &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
830               &   ) * r1_e1e2t(ji,jj)
831            zdt2 = zdt * zdt
832           
833            zten_i(ji,jj) = zdt
834
835            ! shear**2 at T points (doc eq. A16)
836            zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
837               &   + 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)  &
838               &   ) * 0.25_wp * r1_e1e2t(ji,jj)
839           
840            ! shear at T points
841            pshear_i(ji,jj) = SQRT( zdt2 + zds2 )
842
843            ! divergence at T points
844            pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
845               &             + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
846               &             ) * r1_e1e2t(ji,jj)
847           
848            ! delta at T points
849            zfac            = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) ! delta 
850            rswitch         = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zfac ) ) ! 0 if delta=0
851            pdelta_i(ji,jj) = zfac + rn_creepl * rswitch ! delta+creepl
852
853         END DO
854      END DO
855      CALL lbc_lnk_multi( 'icedyn_rhg_eap', pshear_i, 'T', 1., pdivu_i, 'T', 1., pdelta_i, 'T', 1., zten_i, 'T', 1., &
856         &                                  zs1     , 'T', 1., zs2    , 'T', 1., zs12    , 'F', 1. )
857     
858      ! --- Store the stress tensor for the next time step --- !
859      pstress1_i (:,:) = zs1 (:,:)
860      pstress2_i (:,:) = zs2 (:,:)
861      pstress12_i(:,:) = zs12(:,:)
862      !
863
864      !------------------------------------------------------------------------------!
865      ! 5) diagnostics
866      !------------------------------------------------------------------------------!
867      ! --- ice-ocean, ice-atm. & ice-oceanbottom(landfast) stresses --- !
868      IF(  iom_use('utau_oi') .OR. iom_use('vtau_oi') .OR. iom_use('utau_ai') .OR. iom_use('vtau_ai') .OR. &
869         & iom_use('utau_bi') .OR. iom_use('vtau_bi') ) THEN
870         !
871         CALL lbc_lnk_multi( 'icedyn_rhg_eap', ztaux_oi, 'U', -1., ztauy_oi, 'V', -1., ztaux_ai, 'U', -1., ztauy_ai, 'V', -1., &
872            &                                  ztaux_bi, 'U', -1., ztauy_bi, 'V', -1. )
873         !
874         CALL iom_put( 'utau_oi' , ztaux_oi * zmsk00 )
875         CALL iom_put( 'vtau_oi' , ztauy_oi * zmsk00 )
876         CALL iom_put( 'utau_ai' , ztaux_ai * zmsk00 )
877         CALL iom_put( 'vtau_ai' , ztauy_ai * zmsk00 )
878         CALL iom_put( 'utau_bi' , ztaux_bi * zmsk00 )
879         CALL iom_put( 'vtau_bi' , ztauy_bi * zmsk00 )
880      ENDIF
881       
882      ! --- divergence, shear and strength --- !
883      IF( iom_use('icediv') )   CALL iom_put( 'icediv' , pdivu_i  * zmsk00 )   ! divergence
884      IF( iom_use('iceshe') )   CALL iom_put( 'iceshe' , pshear_i * zmsk00 )   ! shear
885      IF( iom_use('icedlt') )   CALL iom_put( 'icedlt' , pdelta_i * zmsk00 )   ! delta
886      IF( iom_use('icestr') )   CALL iom_put( 'icestr' , strength * zmsk00 )   ! strength
887
888      ! --- Stress tensor invariants (SIMIP diags) --- !
889      IF( iom_use('normstr') .OR. iom_use('sheastr') ) THEN
890         !
891         ALLOCATE( zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )
892         !         
893         DO jj = 1, jpj
894            DO ji = 1, jpi
895           
896               ! Ice stresses
897               ! sigma1, sigma2, sigma12 are some useful recombination of the stresses (Hunke and Dukowicz MWR 2002, Bouillon et al., OM2013)
898               ! These are NOT stress tensor components, neither stress invariants, neither stress principal components
899               ! I know, this can be confusing...
900               zfac             =   strength(ji,jj) / ( pdelta_i(ji,jj) + rn_creepl ) 
901               zsig1            =   zfac * ( pdivu_i(ji,jj) - pdelta_i(ji,jj) )
902               zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
903               zsig12           =   zfac * z1_ecc2 * pshear_i(ji,jj)
904               
905               ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008)
906               zsig_I (ji,jj)   =   zsig1 * 0.5_wp                                           ! 1st stress invariant, aka average normal stress, aka negative pressure
907               zsig_II(ji,jj)   =   SQRT ( MAX( 0._wp, zsig2 * zsig2 * 0.25_wp + zsig12 ) )  ! 2nd  ''       '', aka maximum shear stress
908               
909            END DO
910         END DO         
911         !
912         ! Stress tensor invariants (normal and shear stress N/m) - SIMIP diags - definitions following Coon (1974) and Feltham (2008)
913         IF( iom_use('normstr') )   CALL iom_put( 'normstr', zsig_I (:,:) * zmsk00(:,:) ) ! Normal stress
914         IF( iom_use('sheastr') )   CALL iom_put( 'sheastr', zsig_II(:,:) * zmsk00(:,:) ) ! Maximum shear stress
915         
916         DEALLOCATE ( zsig_I, zsig_II )
917         
918      ENDIF
919
920      ! --- Normalized stress tensor principal components --- !
921      ! This are used to plot the normalized yield curve, see Lemieux & Dupont, 2020
922      ! Recommendation 1 : we use ice strength, not replacement pressure
923      ! Recommendation 2 : need to use deformations at PREVIOUS iterate for viscosities
924      IF( iom_use('sig1_pnorm') .OR. iom_use('sig2_pnorm') ) THEN
925         !
926         ALLOCATE( zsig1_p(jpi,jpj) , zsig2_p(jpi,jpj) , zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )         
927         !         
928         DO jj = 1, jpj
929            DO ji = 1, jpi
930           
931               ! Ice stresses computed with **viscosities** (delta, p/delta) at **previous** iterates
932               !                        and **deformations** at current iterates
933               !                        following Lemieux & Dupont (2020)
934               zfac             =   zp_delt(ji,jj)
935               zsig1            =   zfac * ( pdivu_i(ji,jj) - ( zdelta(ji,jj) + rn_creepl ) )
936               zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
937               zsig12           =   zfac * z1_ecc2 * pshear_i(ji,jj)
938               
939               ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008), T-point
940               zsig_I(ji,jj)    =   zsig1 * 0.5_wp                                           ! 1st stress invariant, aka average normal stress, aka negative pressure
941               zsig_II(ji,jj)   =   SQRT ( MAX( 0._wp, zsig2 * zsig2 * 0.25_wp + zsig12 ) )  ! 2nd  ''       '', aka maximum shear stress
942     
943               ! Normalized  principal stresses (used to display the ellipse)
944               z1_strength      =   1._wp / MAX( 1._wp, strength(ji,jj) )
945               zsig1_p(ji,jj)   =   ( zsig_I(ji,jj) + zsig_II(ji,jj) ) * z1_strength
946               zsig2_p(ji,jj)   =   ( zsig_I(ji,jj) - zsig_II(ji,jj) ) * z1_strength
947            END DO
948         END DO               
949         !
950         CALL iom_put( 'sig1_pnorm' , zsig1_p ) 
951         CALL iom_put( 'sig2_pnorm' , zsig2_p ) 
952     
953         DEALLOCATE( zsig1_p , zsig2_p , zsig_I, zsig_II )
954         
955      ENDIF
956
957      ! --- yieldcurve --- !
958      IF( iom_use('yield11') .OR. iom_use('yield12') .OR. iom_use('yield22')) THEN
959
960         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zyield11, 'T', 1., zyield22, 'T', 1., zyield12, 'T', 1. )
961
962         CALL iom_put( 'yield11', zyield11 * zmsk00 )
963         CALL iom_put( 'yield22', zyield22 * zmsk00 )
964         CALL iom_put( 'yield12', zyield12 * zmsk00 )
965      ENDIF
966
967      ! --- anisotropy tensor --- !
968      IF( iom_use('aniso') ) THEN       
969         CALL lbc_lnk( 'icedyn_rhg_eap', paniso_11, 'T', 1. ) 
970         CALL iom_put( 'aniso' , paniso_11 * zmsk00 )
971      ENDIF
972     
973      ! --- SIMIP --- !
974      IF(  iom_use('dssh_dx') .OR. iom_use('dssh_dy') .OR. &
975         & iom_use('corstrx') .OR. iom_use('corstry') .OR. iom_use('intstrx') .OR. iom_use('intstry') ) THEN
976         !
977         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zspgU, 'U', -1., zspgV, 'V', -1., &
978            &                                  zCorU, 'U', -1., zCorV, 'V', -1., zfU, 'U', -1., zfV, 'V', -1. )
979
980         CALL iom_put( 'dssh_dx' , zspgU * zmsk00 )   ! Sea-surface tilt term in force balance (x)
981         CALL iom_put( 'dssh_dy' , zspgV * zmsk00 )   ! Sea-surface tilt term in force balance (y)
982         CALL iom_put( 'corstrx' , zCorU * zmsk00 )   ! Coriolis force term in force balance (x)
983         CALL iom_put( 'corstry' , zCorV * zmsk00 )   ! Coriolis force term in force balance (y)
984         CALL iom_put( 'intstrx' , zfU   * zmsk00 )   ! Internal force term in force balance (x)
985         CALL iom_put( 'intstry' , zfV   * zmsk00 )   ! Internal force term in force balance (y)
986      ENDIF
987
988      IF(  iom_use('xmtrpice') .OR. iom_use('ymtrpice') .OR. &
989         & iom_use('xmtrpsnw') .OR. iom_use('ymtrpsnw') .OR. iom_use('xatrp') .OR. iom_use('yatrp') ) THEN
990         !
991         ALLOCATE( zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) , &
992            &      zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp(jpi,jpj) , zdiag_yatrp(jpi,jpj) )
993         !
994         DO jj = 2, jpjm1
995            DO ji = 2, jpim1
996               ! 2D ice mass, snow mass, area transport arrays (X, Y)
997               zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zmsk00(ji,jj)
998               zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zmsk00(ji,jj)
999
1000               zdiag_xmtrp_ice(ji,jj) = rhoi * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component
1001               zdiag_ymtrp_ice(ji,jj) = rhoi * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) !        ''           Y-   ''
1002
1003               zdiag_xmtrp_snw(ji,jj) = rhos * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
1004               zdiag_ymtrp_snw(ji,jj) = rhos * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) !          ''          Y-   ''
1005
1006               zdiag_xatrp(ji,jj)     = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) )        ! area transport,      X-component
1007               zdiag_yatrp(ji,jj)     = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) )        !        ''            Y-   ''
1008
1009            END DO
1010         END DO
1011
1012         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zdiag_xmtrp_ice, 'U', -1., zdiag_ymtrp_ice, 'V', -1., &
1013            &                                  zdiag_xmtrp_snw, 'U', -1., zdiag_ymtrp_snw, 'V', -1., &
1014            &                                  zdiag_xatrp    , 'U', -1., zdiag_yatrp    , 'V', -1. )
1015
1016         CALL iom_put( 'xmtrpice' , zdiag_xmtrp_ice )   ! X-component of sea-ice mass transport (kg/s)
1017         CALL iom_put( 'ymtrpice' , zdiag_ymtrp_ice )   ! Y-component of sea-ice mass transport
1018         CALL iom_put( 'xmtrpsnw' , zdiag_xmtrp_snw )   ! X-component of snow mass transport (kg/s)
1019         CALL iom_put( 'ymtrpsnw' , zdiag_ymtrp_snw )   ! Y-component of snow mass transport
1020         CALL iom_put( 'xatrp'    , zdiag_xatrp     )   ! X-component of ice area transport
1021         CALL iom_put( 'yatrp'    , zdiag_yatrp     )   ! Y-component of ice area transport
1022
1023         DEALLOCATE( zdiag_xmtrp_ice , zdiag_ymtrp_ice , &
1024            &        zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp , zdiag_yatrp )
1025
1026      ENDIF
1027      !
1028      ! --- convergence tests --- !
1029      IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2 ) THEN
1030         IF( iom_use('uice_cvg') ) THEN
1031            IF( ln_aEVP ) THEN   ! output: beta * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
1032               CALL iom_put( 'uice_cvg', MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * zbeta(:,:) * umask(:,:,1) , &
1033                  &                           ABS( v_ice(:,:) - zv_ice(:,:) ) * zbeta(:,:) * vmask(:,:,1) ) * zmsk15(:,:) )
1034            ELSE                 ! output: nn_nevp * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
1035               CALL iom_put( 'uice_cvg', REAL( nn_nevp ) * MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * umask(:,:,1) , &
1036                  &                                             ABS( v_ice(:,:) - zv_ice(:,:) ) * vmask(:,:,1) ) * zmsk15(:,:) )
1037            ENDIF
1038         ENDIF
1039      ENDIF     
1040      !
1041      DEALLOCATE( zmsk00, zmsk15 )
1042      !
1043   END SUBROUTINE ice_dyn_rhg_eap
1044
1045   
1046   SUBROUTINE rhg_cvg( kt, kiter, kitermax, pu, pv, pub, pvb )
1047      !!----------------------------------------------------------------------
1048      !!                    ***  ROUTINE rhg_cvg  ***
1049      !!                     
1050      !! ** Purpose :   check convergence of oce rheology
1051      !!
1052      !! ** Method  :   create a file ice_cvg.nc containing the convergence of ice velocity
1053      !!                during the sub timestepping of rheology so as:
1054      !!                  uice_cvg = MAX( u(t+1) - u(t) , v(t+1) - v(t) )
1055      !!                This routine is called every sub-iteration, so it is cpu expensive
1056      !!
1057      !! ** Note    :   for the first sub-iteration, uice_cvg is set to 0 (too large otherwise)   
1058      !!----------------------------------------------------------------------
1059      INTEGER ,                 INTENT(in) ::   kt, kiter, kitermax       ! ocean time-step index
1060      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pu, pv, pub, pvb          ! now and before velocities
1061      !!
1062      INTEGER           ::   it, idtime, istatus
1063      INTEGER           ::   ji, jj          ! dummy loop indices
1064      REAL(wp)          ::   zresm           ! local real
1065      CHARACTER(len=20) ::   clname
1066      REAL(wp), DIMENSION(jpi,jpj) ::   zres           ! check convergence
1067      !!----------------------------------------------------------------------
1068
1069      ! create file
1070      IF( kt == nit000 .AND. kiter == 1 ) THEN
1071         !
1072         IF( lwp ) THEN
1073            WRITE(numout,*)
1074            WRITE(numout,*) 'rhg_cvg : ice rheology convergence control'
1075            WRITE(numout,*) '~~~~~~~'
1076         ENDIF
1077         !
1078         IF( lwm ) THEN
1079            clname = 'ice_cvg.nc'
1080            IF( .NOT. Agrif_Root() )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
1081            istatus = NF90_CREATE( TRIM(clname), NF90_CLOBBER, ncvgid )
1082            istatus = NF90_DEF_DIM( ncvgid, 'time'  , NF90_UNLIMITED, idtime )
1083            istatus = NF90_DEF_VAR( ncvgid, 'uice_cvg', NF90_DOUBLE   , (/ idtime /), nvarid )
1084            istatus = NF90_ENDDEF(ncvgid)
1085         ENDIF
1086         !
1087      ENDIF
1088
1089      ! time
1090      it = ( kt - 1 ) * kitermax + kiter
1091     
1092      ! convergence
1093      IF( kiter == 1 ) THEN ! remove the first iteration for calculations of convergence (always very large)
1094         zresm = 0._wp
1095      ELSE
1096         DO jj = 1, jpj
1097            DO ji = 1, jpi
1098               zres(ji,jj) = MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), &
1099                  &               ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * zmsk15(ji,jj)
1100            END DO
1101         END DO
1102         zresm = MAXVAL( zres )
1103         CALL mpp_max( 'icedyn_rhg_evp', zresm )   ! max over the global domain
1104      ENDIF
1105
1106      IF( lwm ) THEN
1107         ! write variables
1108         istatus = NF90_PUT_VAR( ncvgid, nvarid, (/zresm/), (/it/), (/1/) )
1109         ! close file
1110         IF( kt == nitend )   istatus = NF90_CLOSE(ncvgid)
1111      ENDIF
1112     
1113   END SUBROUTINE rhg_cvg
1114
1115
1116   SUBROUTINE update_stress_rdg( ksub, kndte, pdivu, ptension, pshear, pa11, pa12, &
1117      &                          pstressp,  pstressm, pstress12, pstrength, palphar, palphas )
1118      !!---------------------------------------------------------------------
1119      !!                   ***  SUBROUTINE update_stress_rdg  ***
1120      !!                     
1121      !! ** Purpose :   Updates the stress depending on values of strain rate and structure
1122      !!                tensor and for the last subcycle step it computes closing and sliding rate
1123      !!---------------------------------------------------------------------
1124      INTEGER,  INTENT(in   ) ::   ksub, kndte
1125      REAL(wp), INTENT(in   ) ::   pstrength
1126      REAL(wp), INTENT(in   ) ::   pdivu, ptension, pshear
1127      REAL(wp), INTENT(in   ) ::   pa11, pa12
1128      REAL(wp), INTENT(  out) ::   pstressp, pstressm, pstress12
1129      REAL(wp), INTENT(  out) ::   palphar, palphas
1130
1131      INTEGER  ::   kx ,ky, ka
1132
1133      REAL(wp) ::   zstemp11r, zstemp12r, zstemp22r   
1134      REAL(wp) ::   zstemp11s, zstemp12s, zstemp22s 
1135      REAL(wp) ::   za22, zQd11Qd11, zQd11Qd12, zQd12Qd12 
1136      REAL(wp) ::   zQ11Q11, zQ11Q12, zQ12Q12 
1137      REAL(wp) ::   zdtemp11, zdtemp12, zdtemp22 
1138      REAL(wp) ::   zrotstemp11r, zrotstemp12r, zrotstemp22r   
1139      REAL(wp) ::   zrotstemp11s, zrotstemp12s, zrotstemp22s
1140      REAL(wp) ::   zsig11, zsig12, zsig22 
1141      REAL(wp) ::   zsgprm11, zsgprm12, zsgprm22 
1142      REAL(wp) ::   zinvstressconviso 
1143      REAL(wp) ::   zAngle_denom_gamma,  zAngle_denom_alpha 
1144      REAL(wp) ::   zTany_1, zTany_2 
1145      REAL(wp) ::   zx, zy, zdx, zdy, zda, zkxw, kyw, kaw
1146      REAL(wp) ::   zinvdx, zinvdy, zinvda   
1147      REAL(wp) ::   zdtemp1, zdtemp2, zatempprime, zinvsin
1148
1149      REAL(wp), PARAMETER ::   kfriction = 0.45_wp
1150      !!---------------------------------------------------------------------
1151      ! Factor to maintain the same stress as in EVP (see Section 3)
1152      ! Can be set to 1 otherwise
1153!      zinvstressconviso = 1._wp/(1._wp+kfriction*kfriction)
1154      zinvstressconviso = 1._wp
1155 
1156      zinvsin = 1._wp/sin(2._wp*pphi) * zinvstressconviso 
1157      !now uses phi as set in higher code
1158       
1159      ! compute eigenvalues, eigenvectors and angles for structure tensor, strain
1160      ! rates
1161
1162      ! 1) structure tensor
1163      za22 = 1._wp-pa11
1164      zQ11Q11 = 1._wp
1165      zQ12Q12 = rsmall
1166      zQ11Q12 = rsmall
1167 
1168      ! gamma: angle between general coordiantes and principal axis of A
1169      ! here Tan2gamma = 2 a12 / (a11 - a22)
1170      IF((ABS(pa11 - za22) > rsmall).OR.(ABS(pa12) > rsmall)) THEN     
1171         zAngle_denom_gamma = 1._wp/sqrt( ( pa11 - za22 )*( pa11 - za22) + &
1172                              4._wp*pa12*pa12 )
1173   
1174         zQ11Q11 = 0.5_wp + ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Cos^2
1175         zQ12Q12 = 0.5_wp - ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Sin^2
1176         zQ11Q12 = pa12*zAngle_denom_gamma                     !CosSin
1177      ENDIF
1178 
1179      ! rotation Q*atemp*Q^T
1180      zatempprime = zQ11Q11*pa11 + 2.0_wp*zQ11Q12*pa12 + zQ12Q12*za22 
1181     
1182      ! make first principal value the largest
1183      zatempprime = max(zatempprime, 1.0_wp - zatempprime)
1184 
1185      ! 2) strain rate
1186      zdtemp11 = 0.5_wp*(pdivu + ptension)
1187      zdtemp12 = pshear*0.5_wp
1188      zdtemp22 = 0.5_wp*(pdivu - ptension)
1189
1190      ! here Tan2alpha = 2 dtemp12 / (dtemp11 - dtemp22)
1191
1192      zQd11Qd11 = 1.0_wp
1193      zQd12Qd12 = rsmall
1194      zQd11Qd12 = rsmall
1195
1196      IF((ABS( zdtemp11 - zdtemp22) > rsmall).OR. (ABS(zdtemp12) > rsmall)) THEN
1197
1198         zAngle_denom_alpha = 1.0_wp/sqrt( ( zdtemp11 - zdtemp22 )* &
1199                             ( zdtemp11 - zdtemp22 ) + 4.0_wp*zdtemp12*zdtemp12)
1200
1201         zQd11Qd11 = 0.5_wp + ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Cos^2
1202         zQd12Qd12 = 0.5_wp - ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Sin^2
1203         zQd11Qd12 = zdtemp12*zAngle_denom_alpha !CosSin
1204      ENDIF
1205
1206      zdtemp1 = zQd11Qd11*zdtemp11 + 2.0_wp*zQd11Qd12*zdtemp12 + zQd12Qd12*zdtemp22
1207      zdtemp2 = zQd12Qd12*zdtemp11 - 2.0_wp*zQd11Qd12*zdtemp12 + zQd11Qd11*zdtemp22
1208      ! In cos and sin values
1209      zx = 0._wp
1210      IF ((ABS(zdtemp1) > rsmall).OR.(ABS(zdtemp2) > rsmall)) THEN
1211         zx = atan2(zdtemp2,zdtemp1)
1212      ENDIF     
1213               
1214      ! to ensure the angle lies between pi/4 and 9 pi/4
1215      IF (zx < rpi*0.25_wp) zx = zx + rpi*2.0_wp
1216     
1217      ! y: angle between major principal axis of strain rate and structure
1218      ! tensor
1219      ! y = gamma - alpha
1220      ! Expressesed componently with
1221      ! Tany = (Singamma*Cosgamma - Sinalpha*Cosgamma)/(Cos^2gamma - Sin^alpha)
1222     
1223      zTany_1 = zQ11Q12 - zQd11Qd12
1224      zTany_2 = zQ11Q11 - zQd12Qd12
1225     
1226      zy = 0._wp
1227     
1228      IF ((ABS(zTany_1) > rsmall).OR.(ABS(zTany_2) > rsmall)) THEN
1229         zy = atan2(zTany_1,zTany_2)
1230      ENDIF
1231     
1232      ! to make sure y is between 0 and pi
1233      IF (zy > rpi) zy = zy - rpi
1234      IF (zy < 0)  zy = zy + rpi
1235
1236      ! 3) update anisotropic stress tensor
1237      zdx   = rpi/real(nx_yield-1,kind=wp)
1238      zdy   = rpi/real(ny_yield-1,kind=wp)
1239      zda   = 0.5_wp/real(na_yield-1,kind=wp)
1240      zinvdx = 1._wp/zdx
1241      zinvdy = 1._wp/zdy
1242      zinvda = 1._wp/zda
1243
1244      ! % need 8 coords and 8 weights
1245      ! % range in kx
1246      kx  = int((zx-rpi*0.25_wp-rpi)*zinvdx) + 1
1247      !!clem kx  = MAX( 1, MIN( nx_yield-1, INT((zx-rpi*0.25_wp-rpi)*zinvdx) + 1  ) )
1248      zkxw = kx - (zx-rpi*0.25_wp-rpi)*zinvdx 
1249     
1250      ky  = int(zy*zinvdy) + 1
1251      !!clem ky  = MAX( 1, MIN( ny_yield-1, INT(zy*zinvdy) + 1 ) )
1252      kyw = ky - zy*zinvdy 
1253     
1254      ka  = int((zatempprime-0.5_wp)*zinvda) + 1
1255      !!clem ka  = MAX( 1, MIN( na_yield-1, INT((zatempprime-0.5_wp)*zinvda) + 1 ) )
1256      kaw = ka - (zatempprime-0.5_wp)*zinvda
1257     
1258      ! % Determine sigma_r(A1,Zeta,y) and sigma_s (see Section A1 of Tsamados 2013)
1259!!$      zstemp11r =  zkxw  * kyw         * kaw         * s11r(kx  ,ky  ,ka  ) &
1260!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11r(kx+1,ky  ,ka  ) &
1261!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11r(kx  ,ky+1,ka  ) &
1262!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11r(kx  ,ky  ,ka+1) &
1263!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11r(kx+1,ky+1,ka  ) &
1264!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11r(kx+1,ky  ,ka+1) &
1265!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11r(kx  ,ky+1,ka+1) &
1266!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11r(kx+1,ky+1,ka+1)
1267!!$      zstemp12r =  zkxw  * kyw         * kaw         * s12r(kx  ,ky  ,ka  ) &
1268!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12r(kx+1,ky  ,ka  ) &
1269!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12r(kx  ,ky+1,ka  ) &
1270!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12r(kx  ,ky  ,ka+1) &
1271!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12r(kx+1,ky+1,ka  ) &
1272!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12r(kx+1,ky  ,ka+1) &
1273!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12r(kx  ,ky+1,ka+1) &
1274!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12r(kx+1,ky+1,ka+1)
1275!!$      zstemp22r =  zkxw  * kyw         * kaw         * s22r(kx  ,ky  ,ka  ) &
1276!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22r(kx+1,ky  ,ka  ) &
1277!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22r(kx  ,ky+1,ka  ) &
1278!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22r(kx  ,ky  ,ka+1) &
1279!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22r(kx+1,ky+1,ka  ) &
1280!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22r(kx+1,ky  ,ka+1) &
1281!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22r(kx  ,ky+1,ka+1) &
1282!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22r(kx+1,ky+1,ka+1)
1283!!$     
1284!!$      zstemp11s =  zkxw  * kyw         * kaw         * s11s(kx  ,ky  ,ka  ) &
1285!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11s(kx+1,ky  ,ka  ) &
1286!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11s(kx  ,ky+1,ka  ) &
1287!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11s(kx  ,ky  ,ka+1) &
1288!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11s(kx+1,ky+1,ka  ) &
1289!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11s(kx+1,ky  ,ka+1) &
1290!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11s(kx  ,ky+1,ka+1) &
1291!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11s(kx+1,ky+1,ka+1)
1292!!$      zstemp12s =  zkxw  * kyw         * kaw         * s12s(kx  ,ky  ,ka  ) &
1293!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12s(kx+1,ky  ,ka  ) &
1294!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12s(kx  ,ky+1,ka  ) &
1295!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12s(kx  ,ky  ,ka+1) &
1296!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12s(kx+1,ky+1,ka  ) &
1297!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12s(kx+1,ky  ,ka+1) &
1298!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12s(kx  ,ky+1,ka+1) &
1299!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12s(kx+1,ky+1,ka+1)
1300!!$      zstemp22s =  zkxw  * kyw         * kaw         * s22s(kx  ,ky  ,ka  ) &
1301!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22s(kx+1,ky  ,ka  ) &
1302!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22s(kx  ,ky+1,ka  ) &
1303!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22s(kx  ,ky  ,ka+1) &
1304!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22s(kx+1,ky+1,ka  ) &
1305!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22s(kx+1,ky  ,ka+1) &
1306!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22s(kx  ,ky+1,ka+1) &
1307!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22s(kx+1,ky+1,ka+1)
1308
1309      zstemp11r = s11r(kx,ky,ka)
1310      zstemp12r = s12r(kx,ky,ka)
1311      zstemp22r = s22r(kx,ky,ka)
1312      zstemp11s = s11s(kx,ky,ka)
1313      zstemp12s = s12s(kx,ky,ka)
1314      zstemp22s = s22s(kx,ky,ka)
1315     
1316     
1317      ! Calculate mean ice stress over a collection of floes (Equation 3 in
1318      ! Tsamados 2013)
1319
1320      zsig11  = pstrength*(zstemp11r + kfriction*zstemp11s) * zinvsin
1321      zsig12  = pstrength*(zstemp12r + kfriction*zstemp12s) * zinvsin
1322      zsig22  = pstrength*(zstemp22r + kfriction*zstemp22s) * zinvsin
1323
1324      !WRITE(numout,*) 'principal axis stress inside loop', sig11, sig12, sig22
1325
1326      ! Back - rotation of the stress from principal axes into general coordinates
1327
1328      ! Update stress
1329      zsgprm11 = zQ11Q11*zsig11 + zQ12Q12*zsig22 -      2._wp*zQ11Q12 *zsig12
1330      zsgprm12 = zQ11Q12*zsig11 - zQ11Q12*zsig22 + (zQ11Q11 - zQ12Q12)*zsig12
1331      zsgprm22 = zQ12Q12*zsig11 + zQ11Q11*zsig22 +      2._wp*zQ11Q12 *zsig12
1332
1333      pstressp  = zsgprm11 + zsgprm22
1334      pstress12 = zsgprm12
1335      pstressm  = zsgprm11 - zsgprm22
1336
1337      ! Compute ridging and sliding functions in general coordinates
1338      ! (Equation 11 in Tsamados 2013)
1339      IF (ksub == kndte) THEN
1340         zrotstemp11r = zQ11Q11*zstemp11r - 2._wp*zQ11Q12* zstemp12r &
1341                      + zQ12Q12*zstemp22r
1342         zrotstemp12r = zQ11Q11*zstemp12r +       zQ11Q12*(zstemp11r-zstemp22r) &
1343                      - zQ12Q12*zstemp12r
1344         zrotstemp22r = zQ12Q12*zstemp11r + 2._wp*zQ11Q12* zstemp12r & 
1345                      + zQ11Q11*zstemp22r
1346
1347         zrotstemp11s = zQ11Q11*zstemp11s - 2._wp*zQ11Q12* zstemp12s &
1348                      + zQ12Q12*zstemp22s
1349         zrotstemp12s = zQ11Q11*zstemp12s +       zQ11Q12*(zstemp11s-zstemp22s) &
1350                      - zQ12Q12*zstemp12s
1351         zrotstemp22s = zQ12Q12*zstemp11s + 2._wp*zQ11Q12* zstemp12s & 
1352                      + zQ11Q11*zstemp22s
1353
1354         palphar =  zrotstemp11r*zdtemp11 + 2._wp*zrotstemp12r*zdtemp12 &
1355                  + zrotstemp22r*zdtemp22
1356         palphas =  zrotstemp11s*zdtemp11 + 2._wp*zrotstemp12s*zdtemp12 &
1357                  + zrotstemp22s*zdtemp22
1358      ENDIF
1359   END SUBROUTINE update_stress_rdg 
1360
1361!=======================================================================
1362
1363
1364   SUBROUTINE calc_ffrac( pstressp, pstressm, pstress12, pa11, pa12, &
1365      &                   pmresult11, pmresult12 )
1366      !!---------------------------------------------------------------------
1367      !!                   ***  ROUTINE calc_ffrac  ***
1368      !!                     
1369      !! ** Purpose :   Computes term in evolution equation for structure tensor
1370      !!                which determines the ice floe re-orientation due to fracture
1371      !! ** Method :    Eq. 7: Ffrac = -kf(A-S) or = 0 depending on sigma_1 and sigma_2
1372      !!---------------------------------------------------------------------
1373      REAL (wp), INTENT(in)  ::   pstressp, pstressm, pstress12, pa11, pa12
1374      REAL (wp), INTENT(out) ::   pmresult11, pmresult12
1375
1376      ! local variables
1377
1378      REAL (wp) ::   zsigma11, zsigma12, zsigma22  ! stress tensor elements
1379      REAL (wp) ::   zAngle_denom        ! angle with principal component axis
1380      REAL (wp) ::   zsigma_1, zsigma_2   ! principal components of stress
1381      REAL (wp) ::   zQ11, zQ12, zQ11Q11, zQ11Q12, zQ12Q12
1382
1383!!$      REAL (wp), PARAMETER ::   kfrac = 0.0001_wp   ! rate of fracture formation
1384      REAL (wp), PARAMETER ::   kfrac = 1.e-3_wp   ! rate of fracture formation
1385      REAL (wp), PARAMETER ::   threshold = 0.3_wp  ! critical confinement ratio
1386      !!---------------------------------------------------------------
1387      !
1388      zsigma11 = 0.5_wp*(pstressp+pstressm) 
1389      zsigma12 = pstress12 
1390      zsigma22 = 0.5_wp*(pstressp-pstressm) 
1391
1392      ! Here's the change - no longer calculate gamma,
1393      ! use trig formulation, angle signs are kept correct, don't worry
1394   
1395      ! rotate tensor to get into sigma principal axis
1396   
1397      ! here Tan2gamma = 2 sig12 / (sig11 - sig12)
1398      ! add rsmall to the denominator to stop 1/0 errors, makes very little
1399      ! error to the calculated angles < rsmall
1400
1401      zQ11Q11 = 0.1_wp
1402      zQ12Q12 = rsmall
1403      zQ11Q12 = rsmall
1404
1405      IF((ABS( zsigma11 - zsigma22) > rsmall).OR.(ABS(zsigma12) > rsmall)) THEN
1406
1407         zAngle_denom = 1.0_wp/sqrt( ( zsigma11 - zsigma22 )*( zsigma11 - &
1408                       zsigma22 ) + 4.0_wp*zsigma12*zsigma12)
1409
1410         zQ11Q11 = 0.5_wp + ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Cos^2
1411         zQ12Q12 = 0.5_wp - ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Sin^2
1412         zQ11Q12 = zsigma12*zAngle_denom                      !CosSin
1413      ENDIF
1414
1415      zsigma_1 = zQ11Q11*zsigma11 + 2.0_wp*zQ11Q12*zsigma12 + zQ12Q12*zsigma22 ! S(1,1)
1416      zsigma_2 = zQ12Q12*zsigma11 - 2.0_wp*zQ11Q12*zsigma12 + zQ11Q11*zsigma22 ! S(2,2)
1417
1418      ! Pure divergence
1419      IF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 >= 0.0_wp))  THEN
1420         pmresult11 = 0.0_wp
1421         pmresult12 = 0.0_wp
1422
1423      ! Unconfined compression: cracking of blocks not along the axial splitting
1424      ! direction
1425      ! which leads to the loss of their shape, so we again model it through diffusion
1426      ELSEIF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 < 0.0_wp))  THEN
1427         pmresult11 = - kfrac * (pa11 - zQ12Q12)
1428         pmresult12 = - kfrac * (pa12 + zQ11Q12)
1429
1430      ! Shear faulting
1431      ELSEIF (zsigma_2 == 0.0_wp) THEN
1432         pmresult11 = 0.0_wp
1433         pmresult12 = 0.0_wp
1434      ELSEIF ((zsigma_1 <= 0.0_wp).AND.(zsigma_1/zsigma_2 <= threshold)) THEN
1435         pmresult11 = - kfrac * (pa11 - zQ12Q12)
1436         pmresult12 = - kfrac * (pa12 + zQ11Q12)
1437
1438      ! Horizontal spalling
1439      ELSE
1440         pmresult11 = 0.0_wp
1441         pmresult12 = 0.0_wp
1442      ENDIF
1443
1444   END SUBROUTINE calc_ffrac
1445
1446
1447   SUBROUTINE rhg_eap_rst( cdrw, kt )
1448      !!---------------------------------------------------------------------
1449      !!                   ***  ROUTINE rhg_eap_rst  ***
1450      !!                     
1451      !! ** Purpose :   Read or write RHG file in restart file
1452      !!     
1453      !! ** Method  :   use of IOM library
1454      !!----------------------------------------------------------------------
1455      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
1456      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
1457      !
1458      INTEGER  ::   iter                      ! local integer
1459      INTEGER  ::   id1, id2, id3, id4, id5   ! local integers
1460      INTEGER  ::   ix, iy, ip, iz, n, ia     ! local integers
1461   
1462      INTEGER, PARAMETER            ::    nz = 100
1463     
1464      REAL(wp) ::   ainit, xinit, yinit, pinit, zinit
1465      REAL(wp) ::   da, dx, dy, dp, dz, a1
1466
1467      !!clem
1468      REAL(wp) ::   zw1, zw2, zfac, ztemp
1469      REAL(wp) ::   idx, idy, idz
1470
1471      REAL(wp), PARAMETER           ::   eps6 = 1.0e-6_wp
1472      !!----------------------------------------------------------------------
1473      !
1474      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
1475         !                                   ! ---------------
1476         IF( ln_rstart ) THEN                   !* Read the restart file
1477            !
1478            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
1479            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
1480            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
1481            id4 = iom_varid( numrir, 'aniso_11'  , ldstop = .FALSE. )
1482            id5 = iom_varid( numrir, 'aniso_12'  , ldstop = .FALSE. )
1483            !
1484            IF( MIN( id1, id2, id3, id4, id5 ) > 0 ) THEN      ! fields exist
1485               CALL iom_get( numrir, jpdom_autoglo, 'stress1_i' , stress1_i  )
1486               CALL iom_get( numrir, jpdom_autoglo, 'stress2_i' , stress2_i  )
1487               CALL iom_get( numrir, jpdom_autoglo, 'stress12_i', stress12_i )
1488               CALL iom_get( numrir, jpdom_autoglo, 'aniso_11'  , aniso_11   )
1489               CALL iom_get( numrir, jpdom_autoglo, 'aniso_12'  , aniso_12   )
1490            ELSE                                     ! start rheology from rest
1491               IF(lwp) WRITE(numout,*)
1492               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without rheology, set stresses to 0'
1493               stress1_i (:,:) = 0._wp
1494               stress2_i (:,:) = 0._wp
1495               stress12_i(:,:) = 0._wp
1496               aniso_11  (:,:) = 0.5_wp
1497               aniso_12  (:,:) = 0._wp
1498            ENDIF
1499         ELSE                                   !* Start from rest
1500            IF(lwp) WRITE(numout,*)
1501            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set stresses to 0'
1502            stress1_i (:,:) = 0._wp
1503            stress2_i (:,:) = 0._wp
1504            stress12_i(:,:) = 0._wp
1505            aniso_11  (:,:) = 0.5_wp
1506            aniso_12  (:,:) = 0._wp
1507         ENDIF
1508         !
1509
1510         da = 0.5_wp/real(na_yield-1,kind=wp)
1511         ainit = 0.5_wp - da
1512         dx = rpi/real(nx_yield-1,kind=wp)
1513         xinit = rpi + 0.25_wp*rpi - dx
1514         dz = rpi/real(nz,kind=wp)
1515         zinit = -rpi*0.5_wp
1516         dy = rpi/real(ny_yield-1,kind=wp)
1517         yinit = -dy
1518
1519         s11r(:,:,:) = 0._wp
1520         s12r(:,:,:) = 0._wp
1521         s22r(:,:,:) = 0._wp
1522         s11s(:,:,:) = 0._wp
1523         s12s(:,:,:) = 0._wp
1524         s22s(:,:,:) = 0._wp
1525
1526!!$         DO ia=1,na_yield
1527!!$            DO ix=1,nx_yield
1528!!$               DO iy=1,ny_yield
1529!!$                  s11r(ix,iy,ia) = 0._wp
1530!!$                  s12r(ix,iy,ia) = 0._wp
1531!!$                  s22r(ix,iy,ia) = 0._wp
1532!!$                  s11s(ix,iy,ia) = 0._wp
1533!!$                  s12s(ix,iy,ia) = 0._wp
1534!!$                  s22s(ix,iy,ia) = 0._wp
1535!!$                  IF (ia <= na_yield-1) THEN
1536!!$                     DO iz=1,nz
1537!!$                        s11r(ix,iy,ia) = s11r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1538!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1539!!$                           s11kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1540!!$                        s12r(ix,iy,ia) = s12r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1541!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1542!!$                           s12kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1543!!$                        s22r(ix,iy,ia) = s22r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1544!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1545!!$                           s22kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1546!!$                        s11s(ix,iy,ia) = s11s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1547!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1548!!$                           s11ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1549!!$                        s12s(ix,iy,ia) = s12s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1550!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1551!!$                           s12ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1552!!$                        s22s(ix,iy,ia) = s22s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1553!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1554!!$                           s22ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1555!!$                     ENDDO
1556!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1557!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1558!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1559!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1560!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1561!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1562!!$                  ELSE
1563!!$                     s11r(ix,iy,ia) = 0.5_wp*s11kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1564!!$                     s12r(ix,iy,ia) = 0.5_wp*s12kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1565!!$                     s22r(ix,iy,ia) = 0.5_wp*s22kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1566!!$                     s11s(ix,iy,ia) = 0.5_wp*s11ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1567!!$                     s12s(ix,iy,ia) = 0.5_wp*s12ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1568!!$                     s22s(ix,iy,ia) = 0.5_wp*s22ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1569!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1570!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1571!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1572!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1573!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1574!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1575!!$                  ENDIF
1576!!$               ENDDO
1577!!$            ENDDO
1578!!$         ENDDO
1579
1580         !! faster but still very slow => to be improved         
1581         zfac = dz/sin(2._wp*pphi)
1582         DO ia = 1, na_yield-1
1583            zw1 = w1(ainit+ia*da)
1584            zw2 = w2(ainit+ia*da)
1585            DO iz = 1, nz
1586               idz = zinit+iz*dz
1587               ztemp = zw1 * EXP(-zw2*(zinit+iz*dz)*(zinit+iz*dz))
1588               DO iy = 1, ny_yield
1589                  idy = yinit+iy*dy
1590                  DO ix = 1, nx_yield
1591                     idx = xinit+ix*dx
1592                     s11r(ix,iy,ia) = s11r(ix,iy,ia) + ztemp * s11kr(idx,idy,idz)*zfac
1593                     s12r(ix,iy,ia) = s12r(ix,iy,ia) + ztemp * s12kr(idx,idy,idz)*zfac
1594                     s22r(ix,iy,ia) = s22r(ix,iy,ia) + ztemp * s22kr(idx,idy,idz)*zfac 
1595                     s11s(ix,iy,ia) = s11s(ix,iy,ia) + ztemp * s11ks(idx,idy,idz)*zfac
1596                     s12s(ix,iy,ia) = s12s(ix,iy,ia) + ztemp * s12ks(idx,idy,idz)*zfac
1597                     s22s(ix,iy,ia) = s22s(ix,iy,ia) + ztemp * s22ks(idx,idy,idz)*zfac
1598                  END DO
1599               END DO
1600            END DO
1601         END DO
1602
1603         zfac = 1._wp/sin(2._wp*pphi)
1604         ia = na_yield
1605         DO iy = 1, ny_yield
1606            idy = yinit+iy*dy
1607            DO ix = 1, nx_yield
1608               idx = xinit+ix*dx                 
1609               s11r(ix,iy,ia) = 0.5_wp*s11kr(idx,idy,0._wp)*zfac
1610               s12r(ix,iy,ia) = 0.5_wp*s12kr(idx,idy,0._wp)*zfac
1611               s22r(ix,iy,ia) = 0.5_wp*s22kr(idx,idy,0._wp)*zfac
1612               s11s(ix,iy,ia) = 0.5_wp*s11ks(idx,idy,0._wp)*zfac
1613               s12s(ix,iy,ia) = 0.5_wp*s12ks(idx,idy,0._wp)*zfac
1614               s22s(ix,iy,ia) = 0.5_wp*s22ks(idx,idy,0._wp)*zfac
1615            ENDDO
1616         ENDDO
1617         WHERE (ABS(s11r(:,:,:)) < eps6) s11r(:,:,:) = 0._wp
1618         WHERE (ABS(s12r(:,:,:)) < eps6) s12r(:,:,:) = 0._wp
1619         WHERE (ABS(s22r(:,:,:)) < eps6) s22r(:,:,:) = 0._wp
1620         WHERE (ABS(s11s(:,:,:)) < eps6) s11s(:,:,:) = 0._wp
1621         WHERE (ABS(s12s(:,:,:)) < eps6) s12s(:,:,:) = 0._wp
1622         WHERE (ABS(s22s(:,:,:)) < eps6) s22s(:,:,:) = 0._wp
1623
1624
1625      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1626         !                                   ! -------------------
1627         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
1628         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
1629         !
1630         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i  )
1631         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i  )
1632         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
1633         CALL iom_rstput( iter, nitrst, numriw, 'aniso_11'  , aniso_11   )
1634         CALL iom_rstput( iter, nitrst, numriw, 'aniso_12'  , aniso_12   )
1635         !
1636      ENDIF
1637      !
1638   END SUBROUTINE rhg_eap_rst
1639
1640
1641   FUNCTION w1(pa)
1642      !!-------------------------------------------------------------------
1643      !! Function : w1 (see Gaussian function psi in Tsamados et al 2013)
1644      !!-------------------------------------------------------------------
1645      REAL(wp), INTENT(in   ) ::   pa
1646      REAL(wp) ::   w1
1647      !!-------------------------------------------------------------------
1648   
1649      w1 = -   223.87569446_wp &
1650       &   +  2361.21986630_wp*pa &
1651       &   - 10606.56079975_wp*pa*pa &
1652       &   + 26315.50025642_wp*pa*pa*pa &
1653       &   - 38948.30444297_wp*pa*pa*pa*pa &
1654       &   + 34397.72407466_wp*pa*pa*pa*pa*pa &
1655       &   - 16789.98003081_wp*pa*pa*pa*pa*pa*pa &
1656       &   +  3495.82839237_wp*pa*pa*pa*pa*pa*pa*pa
1657   
1658   END FUNCTION w1
1659
1660
1661   FUNCTION w2(pa)
1662      !!-------------------------------------------------------------------
1663      !! Function : w2 (see Gaussian function psi in Tsamados et al 2013)
1664      !!-------------------------------------------------------------------
1665      REAL(wp), INTENT(in   ) ::   pa
1666      REAL(wp) ::   w2
1667      !!-------------------------------------------------------------------
1668   
1669      w2 = -    6670.68911883_wp &
1670       &   +   70222.33061536_wp*pa &
1671       &   -  314871.71525448_wp*pa*pa &
1672       &   +  779570.02793492_wp*pa*pa*pa &
1673       &   - 1151098.82436864_wp*pa*pa*pa*pa &
1674       &   + 1013896.59464498_wp*pa*pa*pa*pa*pa &
1675       &   -  493379.44906738_wp*pa*pa*pa*pa*pa*pa &
1676       &   +  102356.55151800_wp*pa*pa*pa*pa*pa*pa*pa
1677   
1678   END FUNCTION w2
1679
1680   FUNCTION s11kr(px,py,pz) 
1681      !!-------------------------------------------------------------------
1682      !! Function : s11kr
1683      !!-------------------------------------------------------------------
1684      REAL(wp), INTENT(in   ) ::   px,py,pz
1685   
1686      REAL(wp) ::   s11kr, zpih
1687   
1688      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1689      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1690      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1691      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1692      REAL(wp) ::   zd11, zd12, zd22
1693      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1694      REAL(wp) ::   zHen1t2, zHen2t1
1695      !!-------------------------------------------------------------------
1696   
1697      zpih = 0.5_wp*rpi
1698   
1699      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1700      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1701      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1702      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1703      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1704      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1705      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1706      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1707      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1708      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1709      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1710      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1711      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1712      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1713      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1714      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1715   ! In expression of tensor d, with this formulatin d(x)=-d(x+pi)
1716   ! Solution, when diagonalizing always check sgn(a11-a22) if > then keep x else
1717   ! x=x-pi/2
1718      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1719      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1720      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1721      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1722      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1723      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1724   
1725      IF (-zIIn1t2>=rsmall) THEN
1726      zHen1t2 = 1._wp
1727      ELSE
1728      zHen1t2 = 0._wp
1729      ENDIF
1730   
1731      IF (-zIIn2t1>=rsmall) THEN
1732      zHen2t1 = 1._wp
1733      ELSE
1734      zHen2t1 = 0._wp
1735      ENDIF
1736   
1737      s11kr = (- zHen1t2 * zn1t2i11 - zHen2t1 * zn2t1i11)
1738
1739   END FUNCTION s11kr
1740
1741   FUNCTION s12kr(px,py,pz)
1742      !!-------------------------------------------------------------------
1743      !! Function : s12kr
1744      !!-------------------------------------------------------------------
1745      REAL(wp), INTENT(in   ) ::   px,py,pz
1746
1747      REAL(wp) ::   s12kr, zs12r0, zs21r0, zpih
1748
1749      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1750      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1751      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1752      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1753      REAL(wp) ::   zd11, zd12, zd22
1754      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1755      REAL(wp) ::   zHen1t2, zHen2t1
1756      !!-------------------------------------------------------------------
1757      zpih = 0.5_wp*rpi
1758   
1759      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1760      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1761      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1762      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1763      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1764      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1765      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1766      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1767      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1768      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1769      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1770      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1771      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1772      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1773      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1774      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1775      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1776      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1777      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1778      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1779      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1780      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1781   
1782      IF (-zIIn1t2>=rsmall) THEN
1783      zHen1t2 = 1._wp
1784      ELSE
1785      zHen1t2 = 0._wp
1786      ENDIF
1787   
1788      IF (-zIIn2t1>=rsmall) THEN
1789      zHen2t1 = 1._wp
1790      ELSE
1791      zHen2t1 = 0._wp
1792      ENDIF
1793   
1794      zs12r0 = (- zHen1t2 * zn1t2i12 - zHen2t1 * zn2t1i12)
1795      zs21r0 = (- zHen1t2 * zn1t2i21 - zHen2t1 * zn2t1i21)
1796      s12kr=0.5_wp*(zs12r0+zs21r0)
1797   
1798   END FUNCTION s12kr
1799
1800   FUNCTION s22kr(px,py,pz)
1801      !!-------------------------------------------------------------------
1802      !! Function : s22kr
1803      !!-------------------------------------------------------------------
1804      REAL(wp), INTENT(in   ) ::   px,py,pz
1805
1806      REAL(wp) ::   s22kr, zpih
1807
1808      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1809      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1810      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1811      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1812      REAL(wp) ::   zd11, zd12, zd22
1813      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1814      REAL(wp) ::   zHen1t2, zHen2t1
1815      !!-------------------------------------------------------------------
1816      zpih = 0.5_wp*rpi
1817
1818      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1819      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1820      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1821      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1822      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1823      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1824      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1825      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1826      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1827      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1828      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1829      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1830      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1831      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1832      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1833      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1834      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1835      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1836      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1837      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1838      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1839      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1840
1841      IF (-zIIn1t2>=rsmall) THEN
1842      zHen1t2 = 1._wp
1843      ELSE
1844      zHen1t2 = 0._wp
1845      ENDIF
1846
1847      IF (-zIIn2t1>=rsmall) THEN
1848      zHen2t1 = 1._wp
1849      ELSE
1850      zHen2t1 = 0._wp
1851      ENDIF
1852
1853      s22kr = (- zHen1t2 * zn1t2i22 - zHen2t1 * zn2t1i22)
1854
1855   END FUNCTION s22kr
1856
1857   FUNCTION s11ks(px,py,pz)
1858      !!-------------------------------------------------------------------
1859      !! Function : s11ks
1860      !!-------------------------------------------------------------------
1861      REAL(wp), INTENT(in   ) ::   px,py,pz
1862
1863      REAL(wp) ::   s11ks, zpih
1864
1865      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1866      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1867      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1868      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1869      REAL(wp) ::   zd11, zd12, zd22
1870      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1871      REAL(wp) ::   zHen1t2, zHen2t1
1872      !!-------------------------------------------------------------------
1873      zpih = 0.5_wp*rpi
1874
1875      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1876      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1877      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1878      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1879      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1880      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1881      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1882      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1883      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1884      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1885      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1886      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1887      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1888      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1889      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1890      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1891      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1892      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1893      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1894      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1895      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1896      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1897
1898      IF (-zIIn1t2>=rsmall) THEN
1899      zHen1t2 = 1._wp
1900      ELSE
1901      zHen1t2 = 0._wp
1902      ENDIF
1903
1904      IF (-zIIn2t1>=rsmall) THEN
1905      zHen2t1 = 1._wp
1906      ELSE
1907      zHen2t1 = 0._wp
1908      ENDIF
1909
1910      s11ks = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i11 + zHen2t1 * zt2t1i11)
1911
1912   END FUNCTION s11ks
1913
1914   FUNCTION s12ks(px,py,pz)
1915      !!-------------------------------------------------------------------
1916      !! Function : s12ks
1917      !!-------------------------------------------------------------------
1918      REAL(wp), INTENT(in   ) ::   px,py,pz
1919
1920      REAL(wp) ::   s12ks, zs12s0, zs21s0, zpih
1921
1922      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1923      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1924      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1925      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1926      REAL(wp) ::   zd11, zd12, zd22
1927      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1928      REAL(wp) ::   zHen1t2, zHen2t1
1929      !!-------------------------------------------------------------------
1930      zpih = 0.5_wp*rpi
1931
1932      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1933      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1934      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1935      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1936      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1937      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1938      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1939      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1940      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1941      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1942      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1943      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1944      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1945      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1946      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1947      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1948      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1949      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1950      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1951      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1952      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1953      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1954
1955      IF (-zIIn1t2>=rsmall) THEN
1956      zHen1t2 = 1._wp
1957      ELSE
1958      zHen1t2 = 0._wp
1959      ENDIF
1960
1961      IF (-zIIn2t1>=rsmall) THEN
1962      zHen2t1 = 1._wp
1963      ELSE
1964      zHen2t1 = 0._wp
1965      ENDIF
1966
1967      zs12s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i12 + zHen2t1 * zt2t1i12)
1968      zs21s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i21 + zHen2t1 * zt2t1i21)
1969      s12ks=0.5_wp*(zs12s0+zs21s0)
1970
1971   END FUNCTION s12ks
1972
1973   FUNCTION s22ks(px,py,pz) 
1974      !!-------------------------------------------------------------------
1975      !! Function : s22ks
1976      !!-------------------------------------------------------------------
1977      REAL(wp), INTENT(in   ) ::   px,py,pz
1978
1979      REAL(wp) ::   s22ks, zpih
1980
1981      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1982      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1983      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1984      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1985      REAL(wp) ::   zd11, zd12, zd22
1986      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1987      REAL(wp) ::   zHen1t2, zHen2t1
1988      !!-------------------------------------------------------------------
1989      zpih = 0.5_wp*rpi
1990
1991      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1992      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1993      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1994      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1995      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1996      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1997      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1998      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1999      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
2000      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
2001      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
2002      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
2003      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
2004      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
2005      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
2006      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
2007      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
2008      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
2009      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
2010      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
2011      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
2012      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
2013
2014      IF (-zIIn1t2>=rsmall) THEN
2015      zHen1t2 = 1._wp
2016      ELSE
2017      zHen1t2 = 0._wp
2018      ENDIF
2019
2020      IF (-zIIn2t1>=rsmall) THEN
2021      zHen2t1 = 1._wp
2022      ELSE
2023      zHen2t1 = 0._wp
2024      ENDIF
2025
2026      s22ks = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i22 + zHen2t1 * zt2t1i22)
2027
2028   END FUNCTION s22ks
2029
2030#else
2031   !!----------------------------------------------------------------------
2032   !!   Default option         Empty module           NO SI3 sea-ice model
2033   !!----------------------------------------------------------------------
2034#endif
2035
2036   !!==============================================================================
2037END MODULE icedyn_rhg_eap
Note: See TracBrowser for help on using the repository browser.