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/UKMO/NEMO_4.0.4_EAP_rheology_fix/src/ICE – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_EAP_rheology_fix/src/ICE/icedyn_rhg_eap.F90

Last change on this file was 15803, checked in by annkeen, 2 years ago

Change default ice_cvg.nc output to be mean rather than maximum

File size: 107.5 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, Y. Aksenov, C. Rousset)  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, zshear                  ! tension, shear
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) = zmsk00(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) = zmsk00(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_eap( 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           
841            ! maximum shear rate at T points (includes tension, output only)
842            pshear_i(ji,jj) = SQRT( zdt2 + zds2 )
843
844            ! shear at T-points
845            zshear(ji,jj)   = SQRT( zds2 )
846
847            ! divergence at T points
848            pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
849               &             + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
850               &             ) * r1_e1e2t(ji,jj)
851           
852            ! delta at T points
853            zfac            = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) ! delta 
854            rswitch         = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zfac ) ) ! 0 if delta=0
855            pdelta_i(ji,jj) = zfac + rn_creepl * rswitch ! delta+creepl
856
857         END DO
858      END DO
859      CALL lbc_lnk_multi( 'icedyn_rhg_eap', pshear_i, 'T', 1., pdivu_i, 'T', 1., pdelta_i, 'T', 1., zten_i, 'T', 1., &
860         &                                  zs1     , 'T', 1., zs2    , 'T', 1., zs12    , 'F', 1., zshear, 'T', 1. )
861     
862      ! --- Store the stress tensor for the next time step --- !
863      pstress1_i (:,:) = zs1 (:,:)
864      pstress2_i (:,:) = zs2 (:,:)
865      pstress12_i(:,:) = zs12(:,:)
866      !
867
868      !------------------------------------------------------------------------------!
869      ! 5) diagnostics
870      !------------------------------------------------------------------------------!
871      ! --- ice-ocean, ice-atm. & ice-oceanbottom(landfast) stresses --- !
872      IF(  iom_use('utau_oi') .OR. iom_use('vtau_oi') .OR. iom_use('utau_ai') .OR. iom_use('vtau_ai') .OR. &
873         & iom_use('utau_bi') .OR. iom_use('vtau_bi') ) THEN
874         !
875         CALL lbc_lnk_multi( 'icedyn_rhg_eap', ztaux_oi, 'U', -1., ztauy_oi, 'V', -1., ztaux_ai, 'U', -1., ztauy_ai, 'V', -1., &
876            &                                  ztaux_bi, 'U', -1., ztauy_bi, 'V', -1. )
877         !
878         CALL iom_put( 'utau_oi' , ztaux_oi * zmsk00 )
879         CALL iom_put( 'vtau_oi' , ztauy_oi * zmsk00 )
880         CALL iom_put( 'utau_ai' , ztaux_ai * zmsk00 )
881         CALL iom_put( 'vtau_ai' , ztauy_ai * zmsk00 )
882         CALL iom_put( 'utau_bi' , ztaux_bi * zmsk00 )
883         CALL iom_put( 'vtau_bi' , ztauy_bi * zmsk00 )
884      ENDIF
885       
886      ! --- divergence, shear and strength --- !
887      IF( iom_use('icediv') )   CALL iom_put( 'icediv' , pdivu_i  * zmsk00 )   ! divergence
888      IF( iom_use('iceshe') )   CALL iom_put( 'iceshe' , pshear_i * zmsk00 )   ! shear
889      IF( iom_use('icedlt') )   CALL iom_put( 'icedlt' , pdelta_i * zmsk00 )   ! delta
890      IF( iom_use('icestr') )   CALL iom_put( 'icestr' , strength * zmsk00 )   ! strength
891
892      ! --- Stress tensor invariants (SIMIP diags) --- !
893      IF( iom_use('normstr') .OR. iom_use('sheastr') ) THEN
894         !
895         ALLOCATE( zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )
896         !         
897         DO jj = 1, jpj
898            DO ji = 1, jpi
899           
900               ! Ice stresses
901               ! sigma1, sigma2, sigma12 are some useful recombination of the stresses (Hunke and Dukowicz MWR 2002, Bouillon et al., OM2013)
902               ! These are NOT stress tensor components, neither stress invariants, neither stress principal components
903               ! I know, this can be confusing...
904               zfac             =   strength(ji,jj) / ( pdelta_i(ji,jj) + rn_creepl ) 
905               zsig1            =   zfac * ( pdivu_i(ji,jj) - pdelta_i(ji,jj) )
906               zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
907               zsig12           =   zfac * z1_ecc2 * zshear(ji,jj) * 0.5_wp
908               
909               ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008)
910               zsig_I (ji,jj)   =   zsig1 * 0.5_wp                                           ! 1st stress invariant, aka average normal stress, aka negative pressure
911               zsig_II(ji,jj)   =   SQRT( zsig2 * zsig2 * 0.25_wp + zsig12 * zsig12 )        ! 2nd  ''       ''    , aka maximum shear stress
912               
913            END DO
914         END DO         
915         !
916         ! Stress tensor invariants (normal and shear stress N/m) - SIMIP diags - definitions following Coon (1974) and Feltham (2008)
917         IF( iom_use('normstr') )   CALL iom_put( 'normstr', zsig_I (:,:) * zmsk00(:,:) ) ! Normal stress
918         IF( iom_use('sheastr') )   CALL iom_put( 'sheastr', zsig_II(:,:) * zmsk00(:,:) ) ! Maximum shear stress
919         
920         DEALLOCATE ( zsig_I, zsig_II )
921         
922      ENDIF
923
924      ! --- Normalized stress tensor principal components --- !
925      ! This are used to plot the normalized yield curve, see Lemieux & Dupont, 2020
926      ! Recommendation 1 : we use ice strength, not replacement pressure
927      ! Recommendation 2 : need to use deformations at PREVIOUS iterate for viscosities
928      IF( iom_use('sig1_pnorm') .OR. iom_use('sig2_pnorm') ) THEN
929         !
930         ALLOCATE( zsig1_p(jpi,jpj) , zsig2_p(jpi,jpj) , zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )         
931         !         
932         DO jj = 1, jpj
933            DO ji = 1, jpi
934           
935               ! Ice stresses computed with **viscosities** (delta, p/delta) at **previous** iterates
936               !                        and **deformations** at current iterates
937               !                        following Lemieux & Dupont (2020)
938               zfac             =   strength(ji,jj) / ( pdelta_i(ji,jj) + rn_creepl )
939               zsig1            =   zfac * ( pdivu_i(ji,jj) - zdelta(ji,jj) )
940               zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
941               zsig12           =   zfac * z1_ecc2 * zshear(ji,jj) * 0.5_wp
942               
943               ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008), T-point
944               zsig_I(ji,jj)    =   zsig1 * 0.5_wp                                           ! 1st stress invariant, aka average normal stress, aka negative pressure
945               zsig_II(ji,jj)   =   SQRT( zsig2 * zsig2 * 0.25_wp + zsig12 * zsig12 )        ! 2nd  ''       ''    , aka maximum shear stress
946     
947               ! Normalized  principal stresses (used to display the ellipse)
948               z1_strength      =   1._wp / MAX( 1._wp, strength(ji,jj) )
949               zsig1_p(ji,jj)   =   ( zsig_I(ji,jj) + zsig_II(ji,jj) ) * z1_strength
950               zsig2_p(ji,jj)   =   ( zsig_I(ji,jj) - zsig_II(ji,jj) ) * z1_strength
951            END DO
952         END DO               
953         !
954         IF( iom_use('sig1_pnorm') )  CALL iom_put( 'sig1_pnorm' , zsig1_p(:,:) * zmsk00(:,:) ) 
955         IF( iom_use('sig2_pnorm') )  CALL iom_put( 'sig2_pnorm' , zsig2_p(:,:) * zmsk00(:,:) ) 
956     
957         DEALLOCATE( zsig1_p , zsig2_p , zsig_I, zsig_II )
958         
959      ENDIF
960
961      ! --- yieldcurve --- !
962      IF( iom_use('yield11') .OR. iom_use('yield12') .OR. iom_use('yield22')) THEN
963
964         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zyield11, 'T', 1., zyield22, 'T', 1., zyield12, 'T', 1. )
965
966         CALL iom_put( 'yield11', zyield11 * zmsk00 )
967         CALL iom_put( 'yield22', zyield22 * zmsk00 )
968         CALL iom_put( 'yield12', zyield12 * zmsk00 )
969      ENDIF
970
971      ! --- anisotropy tensor --- !
972      IF( iom_use('aniso') ) THEN       
973         CALL lbc_lnk( 'icedyn_rhg_eap', paniso_11, 'T', 1. ) 
974         CALL iom_put( 'aniso' , paniso_11 * zmsk00 )
975      ENDIF
976     
977      ! --- SIMIP --- !
978      IF(  iom_use('dssh_dx') .OR. iom_use('dssh_dy') .OR. &
979         & iom_use('corstrx') .OR. iom_use('corstry') .OR. iom_use('intstrx') .OR. iom_use('intstry') ) THEN
980         !
981         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zspgU, 'U', -1., zspgV, 'V', -1., &
982            &                                  zCorU, 'U', -1., zCorV, 'V', -1., zfU, 'U', -1., zfV, 'V', -1. )
983
984         CALL iom_put( 'dssh_dx' , zspgU * zmsk00 )   ! Sea-surface tilt term in force balance (x)
985         CALL iom_put( 'dssh_dy' , zspgV * zmsk00 )   ! Sea-surface tilt term in force balance (y)
986         CALL iom_put( 'corstrx' , zCorU * zmsk00 )   ! Coriolis force term in force balance (x)
987         CALL iom_put( 'corstry' , zCorV * zmsk00 )   ! Coriolis force term in force balance (y)
988         CALL iom_put( 'intstrx' , zfU   * zmsk00 )   ! Internal force term in force balance (x)
989         CALL iom_put( 'intstry' , zfV   * zmsk00 )   ! Internal force term in force balance (y)
990      ENDIF
991
992      IF(  iom_use('xmtrpice') .OR. iom_use('ymtrpice') .OR. &
993         & iom_use('xmtrpsnw') .OR. iom_use('ymtrpsnw') .OR. iom_use('xatrp') .OR. iom_use('yatrp') ) THEN
994         !
995         ALLOCATE( zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) , &
996            &      zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp(jpi,jpj) , zdiag_yatrp(jpi,jpj) )
997         !
998         DO jj = 2, jpjm1
999            DO ji = 2, jpim1
1000               ! 2D ice mass, snow mass, area transport arrays (X, Y)
1001               zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zmsk00(ji,jj)
1002               zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zmsk00(ji,jj)
1003
1004               zdiag_xmtrp_ice(ji,jj) = rhoi * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component
1005               zdiag_ymtrp_ice(ji,jj) = rhoi * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) !        ''           Y-   ''
1006
1007               zdiag_xmtrp_snw(ji,jj) = rhos * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
1008               zdiag_ymtrp_snw(ji,jj) = rhos * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) !          ''          Y-   ''
1009
1010               zdiag_xatrp(ji,jj)     = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) )        ! area transport,      X-component
1011               zdiag_yatrp(ji,jj)     = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) )        !        ''            Y-   ''
1012
1013            END DO
1014         END DO
1015
1016         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zdiag_xmtrp_ice, 'U', -1., zdiag_ymtrp_ice, 'V', -1., &
1017            &                                  zdiag_xmtrp_snw, 'U', -1., zdiag_ymtrp_snw, 'V', -1., &
1018            &                                  zdiag_xatrp    , 'U', -1., zdiag_yatrp    , 'V', -1. )
1019
1020         CALL iom_put( 'xmtrpice' , zdiag_xmtrp_ice )   ! X-component of sea-ice mass transport (kg/s)
1021         CALL iom_put( 'ymtrpice' , zdiag_ymtrp_ice )   ! Y-component of sea-ice mass transport
1022         CALL iom_put( 'xmtrpsnw' , zdiag_xmtrp_snw )   ! X-component of snow mass transport (kg/s)
1023         CALL iom_put( 'ymtrpsnw' , zdiag_ymtrp_snw )   ! Y-component of snow mass transport
1024         CALL iom_put( 'xatrp'    , zdiag_xatrp     )   ! X-component of ice area transport
1025         CALL iom_put( 'yatrp'    , zdiag_yatrp     )   ! Y-component of ice area transport
1026
1027         DEALLOCATE( zdiag_xmtrp_ice , zdiag_ymtrp_ice , &
1028            &        zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp , zdiag_yatrp )
1029
1030      ENDIF
1031      !
1032      ! --- convergence tests --- !
1033      IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2 ) THEN
1034         IF( iom_use('uice_cvg') ) THEN
1035            IF( ln_aEVP ) THEN   ! output: beta * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
1036               CALL iom_put( 'uice_cvg', MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * zbeta(:,:) * umask(:,:,1) , &
1037                  &                           ABS( v_ice(:,:) - zv_ice(:,:) ) * zbeta(:,:) * vmask(:,:,1) ) * zmsk15(:,:) )
1038            ELSE                 ! output: nn_nevp * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
1039               CALL iom_put( 'uice_cvg', REAL( nn_nevp ) * MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * umask(:,:,1) , &
1040                  &                                             ABS( v_ice(:,:) - zv_ice(:,:) ) * vmask(:,:,1) ) * zmsk15(:,:) )
1041            ENDIF
1042         ENDIF
1043      ENDIF     
1044      !
1045      DEALLOCATE( zmsk00, zmsk15 )
1046      !
1047   END SUBROUTINE ice_dyn_rhg_eap
1048
1049   
1050   SUBROUTINE rhg_cvg_eap( kt, kiter, kitermax, pu, pv, pub, pvb )
1051      !!----------------------------------------------------------------------
1052      !!                    ***  ROUTINE rhg_cvg_eap  ***
1053      !!                     
1054      !! ** Purpose :   check convergence of oce rheology
1055      !!
1056      !! ** Method  :   create a file ice_cvg.nc containing the convergence of ice velocity
1057      !!                during the sub timestepping of rheology so as:
1058      !!                  uice_cvg = MAX( u(t+1) - u(t) , v(t+1) - v(t) )
1059      !!                This routine is called every sub-iteration, so it is cpu expensive
1060      !!
1061      !! ** Note    :   for the first sub-iteration, uice_cvg is set to 0 (too large otherwise)   
1062      !!----------------------------------------------------------------------
1063      INTEGER ,                 INTENT(in) ::   kt, kiter, kitermax       ! ocean time-step index
1064      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pu, pv, pub, pvb          ! now and before velocities
1065      !!
1066      INTEGER           ::   it, idtime, istatus
1067      INTEGER           ::   ji, jj          ! dummy loop indices
1068      REAL(wp)          ::   zresm           ! local real
1069      CHARACTER(len=20) ::   clname
1070      LOGICAL           ::   ll_maxcvg
1071      REAL(wp), DIMENSION(jpi,jpj) ::   zres1           ! check convergence
1072      REAL(wp), DIMENSION(jpi,jpj) ::   zres2           ! check convergence
1073      REAL(wp), DIMENSION(2)         ::   ztmp
1074      !!----------------------------------------------------------------------
1075      ll_maxcvg = .FALSE.
1076
1077      ! create file
1078      IF( kt == nit000 .AND. kiter == 1 ) THEN
1079         !
1080         IF( lwp ) THEN
1081            WRITE(numout,*)
1082            WRITE(numout,*) 'rhg_cvg_eap : ice rheology convergence control'
1083            WRITE(numout,*) '~~~~~~~'
1084         ENDIF
1085         !
1086         IF( lwm ) THEN
1087            clname = 'ice_cvg.nc'
1088            IF( .NOT. Agrif_Root() )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
1089            istatus = NF90_CREATE( TRIM(clname), NF90_CLOBBER, ncvgid )
1090            istatus = NF90_DEF_DIM( ncvgid, 'time'  , NF90_UNLIMITED, idtime )
1091            istatus = NF90_DEF_VAR( ncvgid, 'uice_cvg', NF90_DOUBLE   , (/ idtime /), nvarid )
1092            istatus = NF90_ENDDEF(ncvgid)
1093         ENDIF
1094         !
1095      ENDIF
1096
1097      ! time
1098      it = ( kt - nit000 ) * kitermax + kiter
1099
1100     ! convergence
1101      IF( kiter == 1 ) THEN ! remove the first iteration for calculations of convergence (always very large)
1102         zresm = 0._wp
1103      ELSE
1104         zresm = 0._wp
1105         IF( ll_maxcvg ) THEN   ! error max over the domain
1106           DO jj = 1, jpj
1107              DO ji = 1, jpi
1108                 zres1(ji,jj) = MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), &
1109                  &               ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * zmsk15(ji,jj)
1110                 zres2(ji,jj) = 0._wp
1111              END DO
1112           END DO
1113           zresm = MAXVAL( zres1 )
1114           CALL mpp_max( 'icedyn_rhg_evp', zresm )   ! max over the global domain
1115         ELSE                    ! error averaged over the domain
1116           DO jj = 1, jpj
1117              DO ji = 1, jpi
1118                 zres1(ji,jj) = MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), &
1119                  &                 ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * zmsk15(ji,jj)
1120                 zres2(ji,jj) = zmsk15(ji,jj)
1121              END DO
1122           END DO
1123           ztmp(1) = glob_sum( 'icedyn_rhg_evp', zres1 )
1124           ztmp(2) = glob_sum( 'icedyn_rhg_evp', zres2 )
1125           IF( ztmp(2) /= 0._wp )   zresm = ztmp(1) / ztmp(2)
1126         ENDIF
1127      ENDIF
1128
1129
1130      IF( lwm ) THEN
1131         ! write variables
1132         istatus = NF90_PUT_VAR( ncvgid, nvarid, (/zresm/), (/it/), (/1/) )
1133         ! close file
1134         IF( kt == nitend - nn_fsbc + 1 .AND. kiter == kitermax )   istatus = NF90_CLOSE(ncvgid)
1135      ENDIF
1136     
1137   END SUBROUTINE rhg_cvg_eap
1138
1139   SUBROUTINE update_stress_rdg( ksub, kndte, pdivu, ptension, pshear, pa11, pa12, &
1140      &                          pstressp,  pstressm, pstress12, pstrength, palphar, palphas )
1141      !!---------------------------------------------------------------------
1142      !!                   ***  SUBROUTINE update_stress_rdg  ***
1143      !!                     
1144      !! ** Purpose :   Updates the stress depending on values of strain rate and structure
1145      !!                tensor and for the last subcycle step it computes closing and sliding rate
1146      !!---------------------------------------------------------------------
1147      INTEGER,  INTENT(in   ) ::   ksub, kndte
1148      REAL(wp), INTENT(in   ) ::   pstrength
1149      REAL(wp), INTENT(in   ) ::   pdivu, ptension, pshear
1150      REAL(wp), INTENT(in   ) ::   pa11, pa12
1151      REAL(wp), INTENT(  out) ::   pstressp, pstressm, pstress12
1152      REAL(wp), INTENT(  out) ::   palphar, palphas
1153
1154      INTEGER  ::   kx ,ky, ka
1155
1156      REAL(wp) ::   zstemp11r, zstemp12r, zstemp22r   
1157      REAL(wp) ::   zstemp11s, zstemp12s, zstemp22s 
1158      REAL(wp) ::   za22, zQd11Qd11, zQd11Qd12, zQd12Qd12 
1159      REAL(wp) ::   zQ11Q11, zQ11Q12, zQ12Q12 
1160      REAL(wp) ::   zdtemp11, zdtemp12, zdtemp22 
1161      REAL(wp) ::   zrotstemp11r, zrotstemp12r, zrotstemp22r   
1162      REAL(wp) ::   zrotstemp11s, zrotstemp12s, zrotstemp22s
1163      REAL(wp) ::   zsig11, zsig12, zsig22 
1164      REAL(wp) ::   zsgprm11, zsgprm12, zsgprm22 
1165      REAL(wp) ::   zinvstressconviso 
1166      REAL(wp) ::   zAngle_denom_gamma,  zAngle_denom_alpha 
1167      REAL(wp) ::   zTany_1, zTany_2 
1168      REAL(wp) ::   zx, zy, zdx, zdy, zda, zkxw, kyw, kaw
1169      REAL(wp) ::   zinvdx, zinvdy, zinvda   
1170      REAL(wp) ::   zdtemp1, zdtemp2, zatempprime, zinvsin
1171
1172      REAL(wp), PARAMETER ::   kfriction = 0.45_wp
1173      !!---------------------------------------------------------------------
1174      ! Factor to maintain the same stress as in EVP (see Section 3)
1175      ! Can be set to 1 otherwise
1176!      zinvstressconviso = 1._wp/(1._wp+kfriction*kfriction)
1177      zinvstressconviso = 1._wp
1178 
1179      zinvsin = 1._wp/sin(2._wp*pphi) * zinvstressconviso 
1180      !now uses phi as set in higher code
1181       
1182      ! compute eigenvalues, eigenvectors and angles for structure tensor, strain
1183      ! rates
1184
1185      ! 1) structure tensor
1186      za22 = 1._wp-pa11
1187      zQ11Q11 = 1._wp
1188      zQ12Q12 = rsmall
1189      zQ11Q12 = rsmall
1190 
1191      ! gamma: angle between general coordiantes and principal axis of A
1192      ! here Tan2gamma = 2 a12 / (a11 - a22)
1193      IF((ABS(pa11 - za22) > rsmall).OR.(ABS(pa12) > rsmall)) THEN     
1194         zAngle_denom_gamma = 1._wp/sqrt( ( pa11 - za22 )*( pa11 - za22) + &
1195                              4._wp*pa12*pa12 )
1196   
1197         zQ11Q11 = 0.5_wp + ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Cos^2
1198         zQ12Q12 = 0.5_wp - ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Sin^2
1199         zQ11Q12 = pa12*zAngle_denom_gamma                     !CosSin
1200      ENDIF
1201 
1202      ! rotation Q*atemp*Q^T
1203      zatempprime = zQ11Q11*pa11 + 2.0_wp*zQ11Q12*pa12 + zQ12Q12*za22 
1204     
1205      ! make first principal value the largest
1206      zatempprime = max(zatempprime, 1.0_wp - zatempprime)
1207 
1208      ! 2) strain rate
1209      zdtemp11 = 0.5_wp*(pdivu + ptension)
1210      zdtemp12 = pshear*0.5_wp
1211      zdtemp22 = 0.5_wp*(pdivu - ptension)
1212
1213      ! here Tan2alpha = 2 dtemp12 / (dtemp11 - dtemp22)
1214
1215      zQd11Qd11 = 1.0_wp
1216      zQd12Qd12 = rsmall
1217      zQd11Qd12 = rsmall
1218
1219      IF((ABS( zdtemp11 - zdtemp22) > rsmall).OR. (ABS(zdtemp12) > rsmall)) THEN
1220
1221         zAngle_denom_alpha = 1.0_wp/sqrt( ( zdtemp11 - zdtemp22 )* &
1222                             ( zdtemp11 - zdtemp22 ) + 4.0_wp*zdtemp12*zdtemp12)
1223
1224         zQd11Qd11 = 0.5_wp + ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Cos^2
1225         zQd12Qd12 = 0.5_wp - ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Sin^2
1226         zQd11Qd12 = zdtemp12*zAngle_denom_alpha !CosSin
1227      ENDIF
1228
1229      zdtemp1 = zQd11Qd11*zdtemp11 + 2.0_wp*zQd11Qd12*zdtemp12 + zQd12Qd12*zdtemp22
1230      zdtemp2 = zQd12Qd12*zdtemp11 - 2.0_wp*zQd11Qd12*zdtemp12 + zQd11Qd11*zdtemp22
1231      ! In cos and sin values
1232      zx = 0._wp
1233      IF ((ABS(zdtemp1) > rsmall).OR.(ABS(zdtemp2) > rsmall)) THEN
1234         zx = atan2(zdtemp2,zdtemp1)
1235      ENDIF     
1236               
1237      ! to ensure the angle lies between pi/4 and 9 pi/4
1238      IF (zx < rpi*0.25_wp) zx = zx + rpi*2.0_wp
1239     
1240      ! y: angle between major principal axis of strain rate and structure
1241      ! tensor
1242      ! y = gamma - alpha
1243      ! Expressesed componently with
1244      ! Tany = (Singamma*Cosgamma - Sinalpha*Cosgamma)/(Cos^2gamma - Sin^alpha)
1245     
1246      zTany_1 = zQ11Q12 - zQd11Qd12
1247      zTany_2 = zQ11Q11 - zQd12Qd12
1248     
1249      zy = 0._wp
1250     
1251      IF ((ABS(zTany_1) > rsmall).OR.(ABS(zTany_2) > rsmall)) THEN
1252         zy = atan2(zTany_1,zTany_2)
1253      ENDIF
1254     
1255      ! to make sure y is between 0 and pi
1256      IF (zy > rpi) zy = zy - rpi
1257      IF (zy < 0)  zy = zy + rpi
1258
1259      ! 3) update anisotropic stress tensor
1260      zdx   = rpi/real(nx_yield-1,kind=wp)
1261      zdy   = rpi/real(ny_yield-1,kind=wp)
1262      zda   = 0.5_wp/real(na_yield-1,kind=wp)
1263      zinvdx = 1._wp/zdx
1264      zinvdy = 1._wp/zdy
1265      zinvda = 1._wp/zda
1266
1267      ! % need 8 coords and 8 weights
1268      ! % range in kx
1269      kx  = int((zx-rpi*0.25_wp-rpi)*zinvdx) + 1
1270      !!clem kx  = MAX( 1, MIN( nx_yield-1, INT((zx-rpi*0.25_wp-rpi)*zinvdx) + 1  ) )
1271      zkxw = kx - (zx-rpi*0.25_wp-rpi)*zinvdx 
1272     
1273      ky  = int(zy*zinvdy) + 1
1274      !!clem ky  = MAX( 1, MIN( ny_yield-1, INT(zy*zinvdy) + 1 ) )
1275      kyw = ky - zy*zinvdy 
1276     
1277      ka  = int((zatempprime-0.5_wp)*zinvda) + 1
1278      !!clem ka  = MAX( 1, MIN( na_yield-1, INT((zatempprime-0.5_wp)*zinvda) + 1 ) )
1279      kaw = ka - (zatempprime-0.5_wp)*zinvda
1280     
1281      ! % Determine sigma_r(A1,Zeta,y) and sigma_s (see Section A1 of Tsamados 2013)
1282!!$      zstemp11r =  zkxw  * kyw         * kaw         * s11r(kx  ,ky  ,ka  ) &
1283!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11r(kx+1,ky  ,ka  ) &
1284!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11r(kx  ,ky+1,ka  ) &
1285!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11r(kx  ,ky  ,ka+1) &
1286!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11r(kx+1,ky+1,ka  ) &
1287!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11r(kx+1,ky  ,ka+1) &
1288!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11r(kx  ,ky+1,ka+1) &
1289!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11r(kx+1,ky+1,ka+1)
1290!!$      zstemp12r =  zkxw  * kyw         * kaw         * s12r(kx  ,ky  ,ka  ) &
1291!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12r(kx+1,ky  ,ka  ) &
1292!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12r(kx  ,ky+1,ka  ) &
1293!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12r(kx  ,ky  ,ka+1) &
1294!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12r(kx+1,ky+1,ka  ) &
1295!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12r(kx+1,ky  ,ka+1) &
1296!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12r(kx  ,ky+1,ka+1) &
1297!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12r(kx+1,ky+1,ka+1)
1298!!$      zstemp22r =  zkxw  * kyw         * kaw         * s22r(kx  ,ky  ,ka  ) &
1299!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22r(kx+1,ky  ,ka  ) &
1300!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22r(kx  ,ky+1,ka  ) &
1301!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22r(kx  ,ky  ,ka+1) &
1302!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22r(kx+1,ky+1,ka  ) &
1303!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22r(kx+1,ky  ,ka+1) &
1304!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22r(kx  ,ky+1,ka+1) &
1305!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22r(kx+1,ky+1,ka+1)
1306!!$     
1307!!$      zstemp11s =  zkxw  * kyw         * kaw         * s11s(kx  ,ky  ,ka  ) &
1308!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11s(kx+1,ky  ,ka  ) &
1309!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11s(kx  ,ky+1,ka  ) &
1310!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11s(kx  ,ky  ,ka+1) &
1311!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11s(kx+1,ky+1,ka  ) &
1312!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11s(kx+1,ky  ,ka+1) &
1313!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11s(kx  ,ky+1,ka+1) &
1314!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11s(kx+1,ky+1,ka+1)
1315!!$      zstemp12s =  zkxw  * kyw         * kaw         * s12s(kx  ,ky  ,ka  ) &
1316!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12s(kx+1,ky  ,ka  ) &
1317!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12s(kx  ,ky+1,ka  ) &
1318!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12s(kx  ,ky  ,ka+1) &
1319!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12s(kx+1,ky+1,ka  ) &
1320!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12s(kx+1,ky  ,ka+1) &
1321!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12s(kx  ,ky+1,ka+1) &
1322!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12s(kx+1,ky+1,ka+1)
1323!!$      zstemp22s =  zkxw  * kyw         * kaw         * s22s(kx  ,ky  ,ka  ) &
1324!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22s(kx+1,ky  ,ka  ) &
1325!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22s(kx  ,ky+1,ka  ) &
1326!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22s(kx  ,ky  ,ka+1) &
1327!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22s(kx+1,ky+1,ka  ) &
1328!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22s(kx+1,ky  ,ka+1) &
1329!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22s(kx  ,ky+1,ka+1) &
1330!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22s(kx+1,ky+1,ka+1)
1331
1332      zstemp11r = s11r(kx,ky,ka)
1333      zstemp12r = s12r(kx,ky,ka)
1334      zstemp22r = s22r(kx,ky,ka)
1335      zstemp11s = s11s(kx,ky,ka)
1336      zstemp12s = s12s(kx,ky,ka)
1337      zstemp22s = s22s(kx,ky,ka)
1338     
1339     
1340      ! Calculate mean ice stress over a collection of floes (Equation 3 in
1341      ! Tsamados 2013)
1342
1343      zsig11  = pstrength*(zstemp11r + kfriction*zstemp11s) * zinvsin
1344      zsig12  = pstrength*(zstemp12r + kfriction*zstemp12s) * zinvsin
1345      zsig22  = pstrength*(zstemp22r + kfriction*zstemp22s) * zinvsin
1346
1347      !WRITE(numout,*) 'principal axis stress inside loop', sig11, sig12, sig22
1348
1349      ! Back - rotation of the stress from principal axes into general coordinates
1350
1351      ! Update stress
1352      zsgprm11 = zQ11Q11*zsig11 + zQ12Q12*zsig22 -      2._wp*zQ11Q12 *zsig12
1353      zsgprm12 = zQ11Q12*zsig11 - zQ11Q12*zsig22 + (zQ11Q11 - zQ12Q12)*zsig12
1354      zsgprm22 = zQ12Q12*zsig11 + zQ11Q11*zsig22 +      2._wp*zQ11Q12 *zsig12
1355
1356      pstressp  = zsgprm11 + zsgprm22
1357      pstress12 = zsgprm12
1358      pstressm  = zsgprm11 - zsgprm22
1359
1360      ! Compute ridging and sliding functions in general coordinates
1361      ! (Equation 11 in Tsamados 2013)
1362      IF (ksub == kndte) THEN
1363         zrotstemp11r = zQ11Q11*zstemp11r - 2._wp*zQ11Q12* zstemp12r &
1364                      + zQ12Q12*zstemp22r
1365         zrotstemp12r = zQ11Q11*zstemp12r +       zQ11Q12*(zstemp11r-zstemp22r) &
1366                      - zQ12Q12*zstemp12r
1367         zrotstemp22r = zQ12Q12*zstemp11r + 2._wp*zQ11Q12* zstemp12r & 
1368                      + zQ11Q11*zstemp22r
1369
1370         zrotstemp11s = zQ11Q11*zstemp11s - 2._wp*zQ11Q12* zstemp12s &
1371                      + zQ12Q12*zstemp22s
1372         zrotstemp12s = zQ11Q11*zstemp12s +       zQ11Q12*(zstemp11s-zstemp22s) &
1373                      - zQ12Q12*zstemp12s
1374         zrotstemp22s = zQ12Q12*zstemp11s + 2._wp*zQ11Q12* zstemp12s & 
1375                      + zQ11Q11*zstemp22s
1376
1377         palphar =  zrotstemp11r*zdtemp11 + 2._wp*zrotstemp12r*zdtemp12 &
1378                  + zrotstemp22r*zdtemp22
1379         palphas =  zrotstemp11s*zdtemp11 + 2._wp*zrotstemp12s*zdtemp12 &
1380                  + zrotstemp22s*zdtemp22
1381      ENDIF
1382   END SUBROUTINE update_stress_rdg 
1383
1384!=======================================================================
1385
1386
1387   SUBROUTINE calc_ffrac( pstressp, pstressm, pstress12, pa11, pa12, &
1388      &                   pmresult11, pmresult12 )
1389      !!---------------------------------------------------------------------
1390      !!                   ***  ROUTINE calc_ffrac  ***
1391      !!                     
1392      !! ** Purpose :   Computes term in evolution equation for structure tensor
1393      !!                which determines the ice floe re-orientation due to fracture
1394      !! ** Method :    Eq. 7: Ffrac = -kf(A-S) or = 0 depending on sigma_1 and sigma_2
1395      !!---------------------------------------------------------------------
1396      REAL (wp), INTENT(in)  ::   pstressp, pstressm, pstress12, pa11, pa12
1397      REAL (wp), INTENT(out) ::   pmresult11, pmresult12
1398
1399      ! local variables
1400
1401      REAL (wp) ::   zsigma11, zsigma12, zsigma22  ! stress tensor elements
1402      REAL (wp) ::   zAngle_denom        ! angle with principal component axis
1403      REAL (wp) ::   zsigma_1, zsigma_2   ! principal components of stress
1404      REAL (wp) ::   zQ11, zQ12, zQ11Q11, zQ11Q12, zQ12Q12
1405
1406!!$      REAL (wp), PARAMETER ::   kfrac = 0.0001_wp   ! rate of fracture formation
1407      REAL (wp), PARAMETER ::   kfrac = 1.e-3_wp   ! rate of fracture formation
1408      REAL (wp), PARAMETER ::   threshold = 0.3_wp  ! critical confinement ratio
1409      !!---------------------------------------------------------------
1410      !
1411      zsigma11 = 0.5_wp*(pstressp+pstressm) 
1412      zsigma12 = pstress12 
1413      zsigma22 = 0.5_wp*(pstressp-pstressm) 
1414
1415      ! Here's the change - no longer calculate gamma,
1416      ! use trig formulation, angle signs are kept correct, don't worry
1417   
1418      ! rotate tensor to get into sigma principal axis
1419   
1420      ! here Tan2gamma = 2 sig12 / (sig11 - sig12)
1421      ! add rsmall to the denominator to stop 1/0 errors, makes very little
1422      ! error to the calculated angles < rsmall
1423
1424      zQ11Q11 = 0.1_wp
1425      zQ12Q12 = rsmall
1426      zQ11Q12 = rsmall
1427
1428      IF((ABS( zsigma11 - zsigma22) > rsmall).OR.(ABS(zsigma12) > rsmall)) THEN
1429
1430         zAngle_denom = 1.0_wp/sqrt( ( zsigma11 - zsigma22 )*( zsigma11 - &
1431                       zsigma22 ) + 4.0_wp*zsigma12*zsigma12)
1432
1433         zQ11Q11 = 0.5_wp + ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Cos^2
1434         zQ12Q12 = 0.5_wp - ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Sin^2
1435         zQ11Q12 = zsigma12*zAngle_denom                      !CosSin
1436      ENDIF
1437
1438      zsigma_1 = zQ11Q11*zsigma11 + 2.0_wp*zQ11Q12*zsigma12 + zQ12Q12*zsigma22 ! S(1,1)
1439      zsigma_2 = zQ12Q12*zsigma11 - 2.0_wp*zQ11Q12*zsigma12 + zQ11Q11*zsigma22 ! S(2,2)
1440
1441      ! Pure divergence
1442      IF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 >= 0.0_wp))  THEN
1443         pmresult11 = 0.0_wp
1444         pmresult12 = 0.0_wp
1445
1446      ! Unconfined compression: cracking of blocks not along the axial splitting
1447      ! direction
1448      ! which leads to the loss of their shape, so we again model it through diffusion
1449      ELSEIF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 < 0.0_wp))  THEN
1450         pmresult11 = - kfrac * (pa11 - zQ12Q12)
1451         pmresult12 = - kfrac * (pa12 + zQ11Q12)
1452
1453      ! Shear faulting
1454      ELSEIF (zsigma_2 == 0.0_wp) THEN
1455         pmresult11 = 0.0_wp
1456         pmresult12 = 0.0_wp
1457      ELSEIF ((zsigma_1 <= 0.0_wp).AND.(zsigma_1/zsigma_2 <= threshold)) THEN
1458         pmresult11 = - kfrac * (pa11 - zQ12Q12)
1459         pmresult12 = - kfrac * (pa12 + zQ11Q12)
1460
1461      ! Horizontal spalling
1462      ELSE
1463         pmresult11 = 0.0_wp
1464         pmresult12 = 0.0_wp
1465      ENDIF
1466
1467   END SUBROUTINE calc_ffrac
1468
1469
1470   SUBROUTINE rhg_eap_rst( cdrw, kt )
1471      !!---------------------------------------------------------------------
1472      !!                   ***  ROUTINE rhg_eap_rst  ***
1473      !!                     
1474      !! ** Purpose :   Read or write RHG file in restart file
1475      !!     
1476      !! ** Method  :   use of IOM library
1477      !!----------------------------------------------------------------------
1478      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
1479      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
1480      !
1481      INTEGER  ::   iter                      ! local integer
1482      INTEGER  ::   id1, id2, id3, id4, id5   ! local integers
1483      INTEGER  ::   ix, iy, ip, iz, n, ia     ! local integers
1484   
1485      INTEGER, PARAMETER            ::    nz = 100
1486     
1487      REAL(wp) ::   ainit, xinit, yinit, pinit, zinit
1488      REAL(wp) ::   da, dx, dy, dp, dz, a1
1489
1490      !!clem
1491      REAL(wp) ::   zw1, zw2, zfac, ztemp
1492      REAL(wp) ::   idx, idy, idz
1493
1494      REAL(wp), PARAMETER           ::   eps6 = 1.0e-6_wp
1495      !!----------------------------------------------------------------------
1496      !
1497      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
1498         !                                   ! ---------------
1499         IF( ln_rstart ) THEN                   !* Read the restart file
1500            !
1501            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
1502            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
1503            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
1504            id4 = iom_varid( numrir, 'aniso_11'  , ldstop = .FALSE. )
1505            id5 = iom_varid( numrir, 'aniso_12'  , ldstop = .FALSE. )
1506            !
1507            IF( MIN( id1, id2, id3, id4, id5 ) > 0 ) THEN      ! fields exist
1508               CALL iom_get( numrir, jpdom_autoglo, 'stress1_i' , stress1_i  )
1509               CALL iom_get( numrir, jpdom_autoglo, 'stress2_i' , stress2_i  )
1510               CALL iom_get( numrir, jpdom_autoglo, 'stress12_i', stress12_i )
1511               CALL iom_get( numrir, jpdom_autoglo, 'aniso_11'  , aniso_11   )
1512               CALL iom_get( numrir, jpdom_autoglo, 'aniso_12'  , aniso_12   )
1513            ELSE                                     ! start rheology from rest
1514               IF(lwp) WRITE(numout,*)
1515               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without rheology, set stresses to 0'
1516               stress1_i (:,:) = 0._wp
1517               stress2_i (:,:) = 0._wp
1518               stress12_i(:,:) = 0._wp
1519               aniso_11  (:,:) = 0.5_wp
1520               aniso_12  (:,:) = 0._wp
1521            ENDIF
1522         ELSE                                   !* Start from rest
1523            IF(lwp) WRITE(numout,*)
1524            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set stresses to 0'
1525            stress1_i (:,:) = 0._wp
1526            stress2_i (:,:) = 0._wp
1527            stress12_i(:,:) = 0._wp
1528            aniso_11  (:,:) = 0.5_wp
1529            aniso_12  (:,:) = 0._wp
1530         ENDIF
1531         !
1532
1533         da = 0.5_wp/real(na_yield-1,kind=wp)
1534         ainit = 0.5_wp - da
1535         dx = rpi/real(nx_yield-1,kind=wp)
1536         xinit = rpi + 0.25_wp*rpi - dx
1537         dz = rpi/real(nz,kind=wp)
1538         zinit = -rpi*0.5_wp
1539         dy = rpi/real(ny_yield-1,kind=wp)
1540         yinit = -dy
1541
1542         s11r(:,:,:) = 0._wp
1543         s12r(:,:,:) = 0._wp
1544         s22r(:,:,:) = 0._wp
1545         s11s(:,:,:) = 0._wp
1546         s12s(:,:,:) = 0._wp
1547         s22s(:,:,:) = 0._wp
1548
1549!!$         DO ia=1,na_yield
1550!!$            DO ix=1,nx_yield
1551!!$               DO iy=1,ny_yield
1552!!$                  s11r(ix,iy,ia) = 0._wp
1553!!$                  s12r(ix,iy,ia) = 0._wp
1554!!$                  s22r(ix,iy,ia) = 0._wp
1555!!$                  s11s(ix,iy,ia) = 0._wp
1556!!$                  s12s(ix,iy,ia) = 0._wp
1557!!$                  s22s(ix,iy,ia) = 0._wp
1558!!$                  IF (ia <= na_yield-1) THEN
1559!!$                     DO iz=1,nz
1560!!$                        s11r(ix,iy,ia) = s11r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1561!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1562!!$                           s11kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1563!!$                        s12r(ix,iy,ia) = s12r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1564!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1565!!$                           s12kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1566!!$                        s22r(ix,iy,ia) = s22r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1567!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1568!!$                           s22kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1569!!$                        s11s(ix,iy,ia) = s11s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1570!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1571!!$                           s11ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1572!!$                        s12s(ix,iy,ia) = s12s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1573!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1574!!$                           s12ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1575!!$                        s22s(ix,iy,ia) = s22s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1576!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1577!!$                           s22ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
1578!!$                     ENDDO
1579!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1580!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1581!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1582!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1583!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1584!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1585!!$                  ELSE
1586!!$                     s11r(ix,iy,ia) = 0.5_wp*s11kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1587!!$                     s12r(ix,iy,ia) = 0.5_wp*s12kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1588!!$                     s22r(ix,iy,ia) = 0.5_wp*s22kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1589!!$                     s11s(ix,iy,ia) = 0.5_wp*s11ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1590!!$                     s12s(ix,iy,ia) = 0.5_wp*s12ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1591!!$                     s22s(ix,iy,ia) = 0.5_wp*s22ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
1592!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1593!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1594!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1595!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1596!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1597!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1598!!$                  ENDIF
1599!!$               ENDDO
1600!!$            ENDDO
1601!!$         ENDDO
1602
1603         !! faster but still very slow => to be improved         
1604         zfac = dz/sin(2._wp*pphi)
1605         DO ia = 1, na_yield-1
1606            zw1 = w1(ainit+ia*da)
1607            zw2 = w2(ainit+ia*da)
1608            DO iz = 1, nz
1609               idz = zinit+iz*dz
1610               ztemp = zw1 * EXP(-zw2*(zinit+iz*dz)*(zinit+iz*dz))
1611               DO iy = 1, ny_yield
1612                  idy = yinit+iy*dy
1613                  DO ix = 1, nx_yield
1614                     idx = xinit+ix*dx
1615                     s11r(ix,iy,ia) = s11r(ix,iy,ia) + ztemp * s11kr(idx,idy,idz)*zfac
1616                     s12r(ix,iy,ia) = s12r(ix,iy,ia) + ztemp * s12kr(idx,idy,idz)*zfac
1617                     s22r(ix,iy,ia) = s22r(ix,iy,ia) + ztemp * s22kr(idx,idy,idz)*zfac 
1618                     s11s(ix,iy,ia) = s11s(ix,iy,ia) + ztemp * s11ks(idx,idy,idz)*zfac
1619                     s12s(ix,iy,ia) = s12s(ix,iy,ia) + ztemp * s12ks(idx,idy,idz)*zfac
1620                     s22s(ix,iy,ia) = s22s(ix,iy,ia) + ztemp * s22ks(idx,idy,idz)*zfac
1621                  END DO
1622               END DO
1623            END DO
1624         END DO
1625
1626         zfac = 1._wp/sin(2._wp*pphi)
1627         ia = na_yield
1628         DO iy = 1, ny_yield
1629            idy = yinit+iy*dy
1630            DO ix = 1, nx_yield
1631               idx = xinit+ix*dx                 
1632               s11r(ix,iy,ia) = 0.5_wp*s11kr(idx,idy,0._wp)*zfac
1633               s12r(ix,iy,ia) = 0.5_wp*s12kr(idx,idy,0._wp)*zfac
1634               s22r(ix,iy,ia) = 0.5_wp*s22kr(idx,idy,0._wp)*zfac
1635               s11s(ix,iy,ia) = 0.5_wp*s11ks(idx,idy,0._wp)*zfac
1636               s12s(ix,iy,ia) = 0.5_wp*s12ks(idx,idy,0._wp)*zfac
1637               s22s(ix,iy,ia) = 0.5_wp*s22ks(idx,idy,0._wp)*zfac
1638            ENDDO
1639         ENDDO
1640         WHERE (ABS(s11r(:,:,:)) < eps6) s11r(:,:,:) = 0._wp
1641         WHERE (ABS(s12r(:,:,:)) < eps6) s12r(:,:,:) = 0._wp
1642         WHERE (ABS(s22r(:,:,:)) < eps6) s22r(:,:,:) = 0._wp
1643         WHERE (ABS(s11s(:,:,:)) < eps6) s11s(:,:,:) = 0._wp
1644         WHERE (ABS(s12s(:,:,:)) < eps6) s12s(:,:,:) = 0._wp
1645         WHERE (ABS(s22s(:,:,:)) < eps6) s22s(:,:,:) = 0._wp
1646
1647
1648      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1649         !                                   ! -------------------
1650         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
1651         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
1652         !
1653         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i  )
1654         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i  )
1655         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
1656         CALL iom_rstput( iter, nitrst, numriw, 'aniso_11'  , aniso_11   )
1657         CALL iom_rstput( iter, nitrst, numriw, 'aniso_12'  , aniso_12   )
1658         !
1659      ENDIF
1660      !
1661   END SUBROUTINE rhg_eap_rst
1662
1663
1664   FUNCTION w1(pa)
1665      !!-------------------------------------------------------------------
1666      !! Function : w1 (see Gaussian function psi in Tsamados et al 2013)
1667      !!-------------------------------------------------------------------
1668      REAL(wp), INTENT(in   ) ::   pa
1669      REAL(wp) ::   w1
1670      !!-------------------------------------------------------------------
1671   
1672      w1 = -   223.87569446_wp &
1673       &   +  2361.21986630_wp*pa &
1674       &   - 10606.56079975_wp*pa*pa &
1675       &   + 26315.50025642_wp*pa*pa*pa &
1676       &   - 38948.30444297_wp*pa*pa*pa*pa &
1677       &   + 34397.72407466_wp*pa*pa*pa*pa*pa &
1678       &   - 16789.98003081_wp*pa*pa*pa*pa*pa*pa &
1679       &   +  3495.82839237_wp*pa*pa*pa*pa*pa*pa*pa
1680   
1681   END FUNCTION w1
1682
1683
1684   FUNCTION w2(pa)
1685      !!-------------------------------------------------------------------
1686      !! Function : w2 (see Gaussian function psi in Tsamados et al 2013)
1687      !!-------------------------------------------------------------------
1688      REAL(wp), INTENT(in   ) ::   pa
1689      REAL(wp) ::   w2
1690      !!-------------------------------------------------------------------
1691   
1692      w2 = -    6670.68911883_wp &
1693       &   +   70222.33061536_wp*pa &
1694       &   -  314871.71525448_wp*pa*pa &
1695       &   +  779570.02793492_wp*pa*pa*pa &
1696       &   - 1151098.82436864_wp*pa*pa*pa*pa &
1697       &   + 1013896.59464498_wp*pa*pa*pa*pa*pa &
1698       &   -  493379.44906738_wp*pa*pa*pa*pa*pa*pa &
1699       &   +  102356.55151800_wp*pa*pa*pa*pa*pa*pa*pa
1700   
1701   END FUNCTION w2
1702
1703   FUNCTION s11kr(px,py,pz) 
1704      !!-------------------------------------------------------------------
1705      !! Function : s11kr
1706      !!-------------------------------------------------------------------
1707      REAL(wp), INTENT(in   ) ::   px,py,pz
1708   
1709      REAL(wp) ::   s11kr, zpih
1710   
1711      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1712      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1713      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1714      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1715      REAL(wp) ::   zd11, zd12, zd22
1716      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1717      REAL(wp) ::   zHen1t2, zHen2t1
1718      !!-------------------------------------------------------------------
1719   
1720      zpih = 0.5_wp*rpi
1721   
1722      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1723      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1724      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1725      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1726      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1727      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1728      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1729      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1730      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1731      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1732      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1733      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1734      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1735      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1736      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1737      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1738   ! In expression of tensor d, with this formulatin d(x)=-d(x+pi)
1739   ! Solution, when diagonalizing always check sgn(a11-a22) if > then keep x else
1740   ! x=x-pi/2
1741      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1742      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1743      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1744      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1745      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1746      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1747   
1748      IF (-zIIn1t2>=rsmall) THEN
1749      zHen1t2 = 1._wp
1750      ELSE
1751      zHen1t2 = 0._wp
1752      ENDIF
1753   
1754      IF (-zIIn2t1>=rsmall) THEN
1755      zHen2t1 = 1._wp
1756      ELSE
1757      zHen2t1 = 0._wp
1758      ENDIF
1759   
1760      s11kr = (- zHen1t2 * zn1t2i11 - zHen2t1 * zn2t1i11)
1761
1762   END FUNCTION s11kr
1763
1764   FUNCTION s12kr(px,py,pz)
1765      !!-------------------------------------------------------------------
1766      !! Function : s12kr
1767      !!-------------------------------------------------------------------
1768      REAL(wp), INTENT(in   ) ::   px,py,pz
1769
1770      REAL(wp) ::   s12kr, zs12r0, zs21r0, zpih
1771
1772      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1773      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1774      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1775      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1776      REAL(wp) ::   zd11, zd12, zd22
1777      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1778      REAL(wp) ::   zHen1t2, zHen2t1
1779      !!-------------------------------------------------------------------
1780      zpih = 0.5_wp*rpi
1781   
1782      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1783      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1784      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1785      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1786      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1787      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1788      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1789      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1790      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1791      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1792      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1793      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1794      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1795      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1796      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1797      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1798      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1799      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1800      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1801      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1802      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1803      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1804   
1805      IF (-zIIn1t2>=rsmall) THEN
1806      zHen1t2 = 1._wp
1807      ELSE
1808      zHen1t2 = 0._wp
1809      ENDIF
1810   
1811      IF (-zIIn2t1>=rsmall) THEN
1812      zHen2t1 = 1._wp
1813      ELSE
1814      zHen2t1 = 0._wp
1815      ENDIF
1816   
1817      zs12r0 = (- zHen1t2 * zn1t2i12 - zHen2t1 * zn2t1i12)
1818      zs21r0 = (- zHen1t2 * zn1t2i21 - zHen2t1 * zn2t1i21)
1819      s12kr=0.5_wp*(zs12r0+zs21r0)
1820   
1821   END FUNCTION s12kr
1822
1823   FUNCTION s22kr(px,py,pz)
1824      !!-------------------------------------------------------------------
1825      !! Function : s22kr
1826      !!-------------------------------------------------------------------
1827      REAL(wp), INTENT(in   ) ::   px,py,pz
1828
1829      REAL(wp) ::   s22kr, zpih
1830
1831      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1832      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1833      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1834      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1835      REAL(wp) ::   zd11, zd12, zd22
1836      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1837      REAL(wp) ::   zHen1t2, zHen2t1
1838      !!-------------------------------------------------------------------
1839      zpih = 0.5_wp*rpi
1840
1841      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1842      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1843      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1844      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1845      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1846      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1847      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1848      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1849      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1850      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1851      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1852      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1853      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1854      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1855      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1856      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1857      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1858      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1859      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1860      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1861      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1862      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1863
1864      IF (-zIIn1t2>=rsmall) THEN
1865      zHen1t2 = 1._wp
1866      ELSE
1867      zHen1t2 = 0._wp
1868      ENDIF
1869
1870      IF (-zIIn2t1>=rsmall) THEN
1871      zHen2t1 = 1._wp
1872      ELSE
1873      zHen2t1 = 0._wp
1874      ENDIF
1875
1876      s22kr = (- zHen1t2 * zn1t2i22 - zHen2t1 * zn2t1i22)
1877
1878   END FUNCTION s22kr
1879
1880   FUNCTION s11ks(px,py,pz)
1881      !!-------------------------------------------------------------------
1882      !! Function : s11ks
1883      !!-------------------------------------------------------------------
1884      REAL(wp), INTENT(in   ) ::   px,py,pz
1885
1886      REAL(wp) ::   s11ks, zpih
1887
1888      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1889      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1890      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1891      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1892      REAL(wp) ::   zd11, zd12, zd22
1893      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1894      REAL(wp) ::   zHen1t2, zHen2t1
1895      !!-------------------------------------------------------------------
1896      zpih = 0.5_wp*rpi
1897
1898      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1899      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1900      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1901      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1902      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1903      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1904      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1905      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1906      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1907      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1908      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1909      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1910      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1911      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1912      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1913      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1914      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1915      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1916      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1917      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1918      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1919      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1920
1921      IF (-zIIn1t2>=rsmall) THEN
1922      zHen1t2 = 1._wp
1923      ELSE
1924      zHen1t2 = 0._wp
1925      ENDIF
1926
1927      IF (-zIIn2t1>=rsmall) THEN
1928      zHen2t1 = 1._wp
1929      ELSE
1930      zHen2t1 = 0._wp
1931      ENDIF
1932
1933      s11ks = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i11 + zHen2t1 * zt2t1i11)
1934
1935   END FUNCTION s11ks
1936
1937   FUNCTION s12ks(px,py,pz)
1938      !!-------------------------------------------------------------------
1939      !! Function : s12ks
1940      !!-------------------------------------------------------------------
1941      REAL(wp), INTENT(in   ) ::   px,py,pz
1942
1943      REAL(wp) ::   s12ks, zs12s0, zs21s0, zpih
1944
1945      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
1946      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
1947      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
1948      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
1949      REAL(wp) ::   zd11, zd12, zd22
1950      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
1951      REAL(wp) ::   zHen1t2, zHen2t1
1952      !!-------------------------------------------------------------------
1953      zpih = 0.5_wp*rpi
1954
1955      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
1956      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
1957      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
1958      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
1959      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
1960      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
1961      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
1962      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
1963      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
1964      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
1965      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
1966      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
1967      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
1968      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
1969      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
1970      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
1971      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
1972      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
1973      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
1974      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
1975      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
1976      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
1977
1978      IF (-zIIn1t2>=rsmall) THEN
1979      zHen1t2 = 1._wp
1980      ELSE
1981      zHen1t2 = 0._wp
1982      ENDIF
1983
1984      IF (-zIIn2t1>=rsmall) THEN
1985      zHen2t1 = 1._wp
1986      ELSE
1987      zHen2t1 = 0._wp
1988      ENDIF
1989
1990      zs12s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i12 + zHen2t1 * zt2t1i12)
1991      zs21s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i21 + zHen2t1 * zt2t1i21)
1992      s12ks=0.5_wp*(zs12s0+zs21s0)
1993
1994   END FUNCTION s12ks
1995
1996   FUNCTION s22ks(px,py,pz) 
1997      !!-------------------------------------------------------------------
1998      !! Function : s22ks
1999      !!-------------------------------------------------------------------
2000      REAL(wp), INTENT(in   ) ::   px,py,pz
2001
2002      REAL(wp) ::   s22ks, zpih
2003
2004      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
2005      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
2006      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
2007      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
2008      REAL(wp) ::   zd11, zd12, zd22
2009      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
2010      REAL(wp) ::   zHen1t2, zHen2t1
2011      !!-------------------------------------------------------------------
2012      zpih = 0.5_wp*rpi
2013
2014      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
2015      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
2016      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
2017      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
2018      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
2019      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
2020      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
2021      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
2022      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
2023      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
2024      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
2025      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
2026      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
2027      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
2028      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
2029      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
2030      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
2031      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
2032      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
2033      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
2034      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
2035      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22
2036
2037      IF (-zIIn1t2>=rsmall) THEN
2038      zHen1t2 = 1._wp
2039      ELSE
2040      zHen1t2 = 0._wp
2041      ENDIF
2042
2043      IF (-zIIn2t1>=rsmall) THEN
2044      zHen2t1 = 1._wp
2045      ELSE
2046      zHen2t1 = 0._wp
2047      ENDIF
2048
2049      s22ks = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i22 + zHen2t1 * zt2t1i22)
2050
2051   END FUNCTION s22ks
2052
2053#else
2054   !!----------------------------------------------------------------------
2055   !!   Default option         Empty module           NO SI3 sea-ice model
2056   !!----------------------------------------------------------------------
2057#endif
2058
2059   !!==============================================================================
2060END MODULE icedyn_rhg_eap
Note: See TracBrowser for help on using the repository browser.