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

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
icedyn_rhg_eap.F90 in NEMO/branches/2019/dev_r11842_SI3-10_EAP/tests/ICE_RHEO/MY_SRC – NEMO

source: NEMO/branches/2019/dev_r11842_SI3-10_EAP/tests/ICE_RHEO/MY_SRC/icedyn_rhg_eap.F90 @ 13746

Last change on this file since 13746 was 13746, checked in by stefryn, 3 years ago

update test case

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