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 @ 12263

Last change on this file since 12263 was 12263, checked in by stefryn, 4 years ago

add test case for rheology, will change lateral bounday conditions later

File size: 88.0 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   IMPLICIT NONE
46   PRIVATE
47
48   PUBLIC   ice_dyn_rhg_eap   ! called by icedyn_rhg.F90
49   PUBLIC   rhg_eap_rst       ! called by icedyn_rhg.F90
50
51   REAL(wp), PARAMETER           ::   phi = 3.141592653589793_wp/12._wp    ! diamond shaped floe smaller angle (default phi = 30 deg)
52
53   ! look-up table for calculating structure tensor
54   INTEGER, PARAMETER ::   nx_yield = 41
55   INTEGER, PARAMETER ::   ny_yield = 41
56   INTEGER, PARAMETER ::   na_yield = 21
57
58   REAL(wp), DIMENSION(nx_yield, ny_yield, na_yield) ::   s11r, s12r, s22r, s11s, s12s, s22s
59
60   !! * Substitutions
61#  include "vectopt_loop_substitute.h90"
62   !!----------------------------------------------------------------------
63   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
64   !! $Id: icedyn_rhg_eap.F90 11536 2019-09-11 13:54:18Z smasson $
65   !! Software governed by the CeCILL license (see ./LICENSE)
66   !!----------------------------------------------------------------------
67CONTAINS
68
69   SUBROUTINE ice_dyn_rhg_eap( kt, pstress1_i, pstress2_i, pstress12_i, pshear_i, pdivu_i, pdelta_i, paniso_11, paniso_12 )
70      !!-------------------------------------------------------------------
71      !!                 ***  SUBROUTINE ice_dyn_rhg_eap  ***
72      !!                             EAP-C-grid
73      !!
74      !! ** purpose : determines sea ice drift from wind stress, ice-ocean
75      !!  stress and sea-surface slope. Ice-ice interaction is described by
76      !!  a non-linear elasto-anisotropic-plastic (EAP) law including shear
77      !!  strength and a bulk rheology . 
78      !!
79      !!  The points in the C-grid look like this, dear reader
80      !!
81      !!                              (ji,jj)
82      !!                                 |
83      !!                                 |
84      !!                      (ji-1,jj)  |  (ji,jj)
85      !!                             ---------   
86      !!                            |         |
87      !!                            | (ji,jj) |------(ji,jj)
88      !!                            |         |
89      !!                             ---------   
90      !!                     (ji-1,jj-1)     (ji,jj-1)
91      !!
92      !! ** Inputs  : - wind forcing (stress), oceanic currents
93      !!                ice total volume (vt_i) per unit area
94      !!                snow total volume (vt_s) per unit area
95      !!
96      !! ** Action  : - compute u_ice, v_ice : the components of the
97      !!                sea-ice velocity vector
98      !!              - compute delta_i, shear_i, divu_i, which are inputs
99      !!                of the ice thickness distribution
100      !!
101      !! ** Steps   : 0) compute mask at F point
102      !!              1) Compute ice snow mass, ice strength
103      !!              2) Compute wind, oceanic stresses, mass terms and
104      !!                 coriolis terms of the momentum equation
105      !!              3) Solve the momentum equation (iterative procedure)
106      !!              4) Recompute delta, shear and divergence
107      !!                 (which are inputs of the ITD) & store stress
108      !!                 for the next time step
109      !!              5) Diagnostics including charge ellipse
110      !!
111      !! ** Notes   : There is the possibility to use aEVP from the nice work of Kimmritz et al. (2016 & 2017)
112      !!              by setting up ln_aEVP=T (i.e. changing alpha and beta parameters).
113      !!              This is an upgraded version of mEVP from Bouillon et al. 2013
114      !!              (i.e. more stable and better convergence)
115      !!
116      !! References : Hunke and Dukowicz, JPO97
117      !!              Bouillon et al., Ocean Modelling 2009
118      !!              Bouillon et al., Ocean Modelling 2013
119      !!              Kimmritz et al., Ocean Modelling 2016 & 2017
120      !!-------------------------------------------------------------------
121      INTEGER                 , INTENT(in   ) ::   kt                                    ! time step
122      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pstress1_i, pstress2_i, pstress12_i   !
123      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pshear_i  , pdivu_i   , pdelta_i      !
124      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   paniso_11 , paniso_12                 ! structure tensor components
125      !!
126      INTEGER ::   ji, jj, ij0, ij1, ii0, ii1     ! dummy loop indices
127      INTEGER ::   jter         ! local integers
128      !
129      REAL(wp) ::   zrhoco                                              ! rau0 * rn_cio
130      REAL(wp) ::   zdtevp, z1_dtevp                                    ! time step for subcycling
131      REAL(wp) ::   ecc2, z1_ecc2                                       ! square of yield ellipse eccenticity
132      REAL(wp) ::   zalph1, z1_alph1                                    ! alpha coef from Bouillon 2009 or Kimmritz 2017
133      REAL(wp) ::   zm1, zm2, zm3, zmassU, zmassV, zvU, zvV             ! ice/snow mass and volume
134      REAL(wp) ::   zdelta, zp_delf, zds2, zdt, zdt2, zdiv, zdiv2, zdsT ! temporary scalars
135      REAL(wp) ::   zTauO, zTauB, zRHS, zvel                            ! temporary scalars
136      REAL(wp) ::   zkt                                                 ! isotropic tensile strength for landfast ice
137      REAL(wp) ::   zvCr                                                ! critical ice volume above which ice is landfast
138      !
139      REAL(wp) ::   zresm                                               ! Maximal error on ice velocity
140      REAL(wp) ::   zintb, zintn                                        ! dummy argument
141      REAL(wp) ::   zfac_x, zfac_y
142      REAL(wp) ::   zshear, zdum1, zdum2
143      REAL(wp) ::   stressptmp, stressmtmp, stress12tmpF                ! anisotropic stress tensor components
144      REAL(wp) ::   alphar, alphas                                      ! for mechanical redistribution
145      REAL(wp) ::   invw                                                ! for test case
146      !
147      REAL(wp), DIMENSION(jpi,jpj) ::   stress12tmp                     ! anisotropic stress tensor component for regridding
148      REAL(wp), DIMENSION(jpi,jpj) ::   zp_delt                         ! P/delta at T points
149      REAL(wp), DIMENSION(jpi,jpj) ::   zbeta                           ! beta coef from Kimmritz 2017
150      !
151      REAL(wp), DIMENSION(jpi,jpj) ::   zdt_m                           ! (dt / ice-snow_mass) on T points
152      REAL(wp), DIMENSION(jpi,jpj) ::   zaU  , zaV                      ! ice fraction on U/V points
153      REAL(wp), DIMENSION(jpi,jpj) ::   zmU_t, zmV_t                    ! (ice-snow_mass / dt) on U/V points
154      REAL(wp), DIMENSION(jpi,jpj) ::   zmf                             ! coriolis parameter at T points
155      REAL(wp), DIMENSION(jpi,jpj) ::   v_oceU, u_oceV, v_iceU, u_iceV  ! ocean/ice u/v component on V/U points                           
156      !
157      REAL(wp), DIMENSION(jpi,jpj) ::   zds                             ! shear
158      REAL(wp), DIMENSION(jpi,jpj) ::   zs1, zs2, zs12                  ! stress tensor components
159!!$      REAL(wp), DIMENSION(jpi,jpj) ::   zu_ice, zv_ice, zresr           ! check convergence
160      REAL(wp), DIMENSION(jpi,jpj) ::   zsshdyn                         ! array used for the calculation of ice surface slope:
161      !                                                                 !    ocean surface (ssh_m) if ice is not embedded
162      !                                                                 !    ice bottom surface if ice is embedded   
163      REAL(wp), DIMENSION(jpi,jpj) ::   zfU  , zfV                      ! internal stresses
164      REAL(wp), DIMENSION(jpi,jpj) ::   zspgU, zspgV                    ! surface pressure gradient at U/V points
165      REAL(wp), DIMENSION(jpi,jpj) ::   zCorU, zCorV                    ! Coriolis stress array
166      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_ai, ztauy_ai              ! ice-atm. stress at U-V points
167      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_oi, ztauy_oi              ! ice-ocean stress at U-V points
168      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_bi, ztauy_bi              ! ice-OceanBottom stress at U-V points (landfast)
169      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_base, ztauy_base          ! ice-bottom stress at U-V points (landfast)
170      !
171      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk01x, zmsk01y                ! dummy arrays
172      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk00x, zmsk00y                ! mask for ice presence
173      REAL(wp), DIMENSION(jpi,jpj) ::   zfmask, zwf                     ! mask at F points for the ice
174
175      REAL(wp), PARAMETER          ::   zepsi  = 1.0e-20_wp             ! tolerance parameter
176      REAL(wp), PARAMETER          ::   zmmin  = 1._wp                  ! ice mass (kg/m2)  below which ice velocity becomes very small
177      REAL(wp), PARAMETER          ::   zamin  = 0.001_wp               ! ice concentration below which ice velocity becomes very small
178      !! --- diags
179      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk00
180      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zsig1, zsig2, zsig3
181      !! --- SIMIP diags
182      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s)
183      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s)
184      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s)
185      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s)
186      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xatrp     ! X-component of area transport (m2/s)
187      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_yatrp     ! Y-component of area transport (m2/s)     
188      !!-------------------------------------------------------------------
189
190      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_dyn_rhg_eap: EAP sea-ice rheology'
191      !
192!!gm for Clem:  OPTIMIZATION:  I think zfmask can be computed one for all at the initialization....
193      !------------------------------------------------------------------------------!
194      ! 0) mask at F points for the ice
195      !------------------------------------------------------------------------------!
196      ! ocean/land mask
197      DO jj = 1, jpjm1
198         DO ji = 1, jpim1      ! NO vector opt.
199            zfmask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
200         END DO
201      END DO
202      CALL lbc_lnk( 'icedyn_rhg_eap', zfmask, 'F', 1._wp )
203
204      ! Lateral boundary conditions on velocity (modify zfmask)
205      zwf(:,:) = zfmask(:,:)
206      DO jj = 2, jpjm1
207         DO ji = fs_2, fs_jpim1   ! vector opt.
208            IF( zfmask(ji,jj) == 0._wp ) THEN
209               zfmask(ji,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,jj), zwf(ji,jj+1), zwf(ji-1,jj), zwf(ji,jj-1) ) )
210            ENDIF
211         END DO
212      END DO
213      DO jj = 2, jpjm1
214         IF( zfmask(1,jj) == 0._wp ) THEN
215            zfmask(1  ,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(2,jj), zwf(1,jj+1), zwf(1,jj-1) ) )
216         ENDIF
217         IF( zfmask(jpi,jj) == 0._wp ) THEN
218            zfmask(jpi,jj) = rn_ishlat * MIN( 1._wp , MAX( zwf(jpi,jj+1), zwf(jpim1,jj), zwf(jpi,jj-1) ) )
219         ENDIF
220      END DO
221      DO ji = 2, jpim1
222         IF( zfmask(ji,1) == 0._wp ) THEN
223            zfmask(ji,1  ) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,1), zwf(ji,2), zwf(ji-1,1) ) )
224         ENDIF
225         IF( zfmask(ji,jpj) == 0._wp ) THEN
226            zfmask(ji,jpj) = rn_ishlat * MIN( 1._wp , MAX( zwf(ji+1,jpj), zwf(ji-1,jpj), zwf(ji,jpjm1) ) )
227         ENDIF
228      END DO
229      CALL lbc_lnk( 'icedyn_rhg_eap', zfmask, 'F', 1._wp )
230
231      !------------------------------------------------------------------------------!
232      ! 1) define some variables and initialize arrays
233      !------------------------------------------------------------------------------!
234      zrhoco = rau0 * rn_cio 
235!extra code for test case boundary conditions
236      invw=1._wp/(zrhoco*0.5_wp)
237
238      ! ecc2: square of yield ellipse eccenticrity
239      ecc2    = rn_ecc * rn_ecc
240      z1_ecc2 = 1._wp / ecc2
241
242      ! Time step for subcycling
243      zdtevp   = rdt_ice / REAL( nn_nevp )
244      z1_dtevp = 1._wp / zdtevp
245
246      ! alpha parameters (Bouillon 2009)
247      IF( .NOT. ln_aEVP ) THEN
248         zalph1 = ( 2._wp * rn_relast * rdt_ice ) * z1_dtevp
249         z1_alph1 = 1._wp / ( zalph1 + 1._wp )
250      ENDIF
251         
252      ! Initialise stress tensor
253      zs1 (:,:) = pstress1_i (:,:) 
254      zs2 (:,:) = pstress2_i (:,:)
255      zs12(:,:) = pstress12_i(:,:)
256
257      ! Ice strength
258      CALL ice_strength
259
260      ! landfast param from Lemieux(2016): add isotropic tensile strength (following Konig Beatty and Holland, 2010)
261      IF( ln_landfast_L16 ) THEN   ;   zkt = rn_tensile
262      ELSE                         ;   zkt = 0._wp
263      ENDIF
264      !
265      !------------------------------------------------------------------------------!
266      ! 2) Wind / ocean stress, mass terms, coriolis terms
267      !------------------------------------------------------------------------------!
268      ! sea surface height
269      !    embedded sea ice: compute representative ice top surface
270      !    non-embedded sea ice: use ocean surface for slope calculation
271      zsshdyn(:,:) = ice_var_sshdyn( ssh_m, snwice_mass, snwice_mass_b)
272
273      DO jj = 2, jpjm1
274         DO ji = fs_2, fs_jpim1
275
276            ! ice fraction at U-V points
277            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)
278            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)
279
280            ! Ice/snow mass at U-V points
281            zm1 = ( rhos * vt_s(ji  ,jj  ) + rhoi * vt_i(ji  ,jj  ) )
282            zm2 = ( rhos * vt_s(ji+1,jj  ) + rhoi * vt_i(ji+1,jj  ) )
283            zm3 = ( rhos * vt_s(ji  ,jj+1) + rhoi * vt_i(ji  ,jj+1) )
284            zmassU = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm2 * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1)
285            zmassV = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm3 * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1)
286
287            ! Ocean currents at U-V points
288            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)
289            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)
290
291            ! Coriolis at T points (m*f)
292            zmf(ji,jj)      = zm1 * ff_t(ji,jj)
293
294            ! dt/m at T points (for alpha and beta coefficients)
295            zdt_m(ji,jj)    = zdtevp / MAX( zm1, zmmin )
296           
297            ! m/dt
298            zmU_t(ji,jj)    = zmassU * z1_dtevp
299            zmV_t(ji,jj)    = zmassV * z1_dtevp
300           
301            ! Drag ice-atm.
302            ztaux_ai(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj)
303            ztauy_ai(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj)
304
305            ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points
306            zspgU(ji,jj)    = - zmassU * grav * ( zsshdyn(ji+1,jj) - zsshdyn(ji,jj) ) * r1_e1u(ji,jj)
307            zspgV(ji,jj)    = - zmassV * grav * ( zsshdyn(ji,jj+1) - zsshdyn(ji,jj) ) * r1_e2v(ji,jj)
308
309            ! masks
310            zmsk00x(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) )  ! 0 if no ice
311            zmsk00y(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) )  ! 0 if no ice
312
313            ! switches
314            IF( zmassU <= zmmin .AND. zaU(ji,jj) <= zamin ) THEN   ;   zmsk01x(ji,jj) = 0._wp
315            ELSE                                                   ;   zmsk01x(ji,jj) = 1._wp   ;   ENDIF
316            IF( zmassV <= zmmin .AND. zaV(ji,jj) <= zamin ) THEN   ;   zmsk01y(ji,jj) = 0._wp
317            ELSE                                                   ;   zmsk01y(ji,jj) = 1._wp   ;   ENDIF
318
319         END DO
320      END DO
321      CALL lbc_lnk_multi( 'icedyn_rhg_eap', zmf, 'T', 1., zdt_m, 'T', 1. )
322      !
323      !                                  !== Landfast ice parameterization ==!
324      !
325      IF( ln_landfast_L16 ) THEN         !-- Lemieux 2016
326         DO jj = 2, jpjm1
327            DO ji = fs_2, fs_jpim1
328               ! ice thickness at U-V points
329               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)
330               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)
331               ! ice-bottom stress at U points
332               zvCr = zaU(ji,jj) * rn_depfra * hu_n(ji,jj)
333               ztaux_base(ji,jj) = - rn_icebfr * MAX( 0._wp, zvU - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaU(ji,jj) ) )
334               ! ice-bottom stress at V points
335               zvCr = zaV(ji,jj) * rn_depfra * hv_n(ji,jj)
336               ztauy_base(ji,jj) = - rn_icebfr * MAX( 0._wp, zvV - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaV(ji,jj) ) )
337               ! ice_bottom stress at T points
338               zvCr = at_i(ji,jj) * rn_depfra * ht_n(ji,jj)
339               tau_icebfr(ji,jj) = - rn_icebfr * MAX( 0._wp, vt_i(ji,jj) - zvCr ) * EXP( -rn_crhg * ( 1._wp - at_i(ji,jj) ) )
340            END DO
341         END DO
342         CALL lbc_lnk( 'icedyn_rhg_eap', tau_icebfr(:,:), 'T', 1. )
343         !
344      ELSE                               !-- no landfast
345         DO jj = 2, jpjm1
346            DO ji = fs_2, fs_jpim1
347               ztaux_base(ji,jj) = 0._wp
348               ztauy_base(ji,jj) = 0._wp
349            END DO
350         END DO
351      ENDIF
352
353      !------------------------------------------------------------------------------!
354      ! 3) Solution of the momentum equation, iterative procedure
355      !------------------------------------------------------------------------------!
356      !
357      !                                               ! ==================== !
358      DO jter = 1 , nn_nevp                           !    loop over jter    !
359         !                                            ! ==================== !       
360         l_full_nf_update = jter == nn_nevp   ! false: disable full North fold update (performances) for iter = 1 to nn_nevp-1
361         !
362!!$         IF(ln_ctl) THEN   ! Convergence test
363!!$            DO jj = 1, jpjm1
364!!$               zu_ice(:,jj) = u_ice(:,jj) ! velocity at previous time step
365!!$               zv_ice(:,jj) = v_ice(:,jj)
366!!$            END DO
367!!$         ENDIF
368
369         ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- !
370         DO jj = 1, jpjm1         ! loops start at 1 since there is no boundary condition (lbc_lnk) at i=1 and j=1 for F points
371            DO ji = 1, jpim1
372
373               ! shear at F points
374               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)   &
375                  &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
376                  &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
377
378            END DO
379         END DO
380         CALL lbc_lnk( 'icedyn_rhg_eap', zds, 'F', 1. )
381
382         DO jj = 2, jpj    ! loop to jpi,jpj to avoid making a communication for zs1,zs2,zs12
383            DO ji = 2, jpi ! no vector loop
384
385               ! shear at T points
386               zdsT = ( zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
387                  &   + zds(ji,jj-1) * e1e2f(ji,jj-1) + zds(ji-1,jj-1) * e1e2f(ji-1,jj-1)  &
388                  &   ) * 0.25_wp * r1_e1e2t(ji,jj)
389               !WRITE(numout,*) 'shear output', ji, jj, zdsT
390               
391
392               ! shear**2 at T points (doc eq. A16)
393               zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
394                  &   + 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)  &
395                  &   ) * 0.25_wp * r1_e1e2t(ji,jj)
396             
397               ! divergence at T points
398               zdiv  = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
399                  &    + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
400                  &    ) * r1_e1e2t(ji,jj)
401               zdiv2 = zdiv * zdiv
402               
403               ! tension at T points
404               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)   &
405                  &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
406                  &   ) * r1_e1e2t(ji,jj)
407               zdt2 = zdt * zdt
408               
409               ! --- anisotropic stress calculation --- !
410               CALL update_stress_rdg (jter, nn_nevp, phi, zdiv, zdt, &
411                                       zdsT, paniso_11(ji,jj), paniso_12(ji,jj), &
412                                       stressptmp, stressmtmp, &
413                                       stress12tmp(ji,jj), strength(ji,jj), &
414                                       alphar, alphas)
415               IF (jter == nn_nevp) THEN
416                  rdg_conv(ji,jj) = -min( alphar, 0._wp)   
417               ENDIF
418
419               ! delta at T points
420               zdelta = SQRT( zdiv2 + ( zdt2 + zds2 ) * z1_ecc2 ) 
421
422               ! P/delta at T points
423               zp_delt(ji,jj) = strength(ji,jj) / ( zdelta + rn_creepl )
424
425               ! alpha & beta for aEVP
426               !   gamma = 0.5*P/(delta+creepl) * (c*pi)**2/Area * dt/m
427               !   alpha = beta = sqrt(4*gamma)
428               IF( ln_aEVP ) THEN
429                  zalph1   = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
430                  z1_alph1 = 1._wp / ( zalph1 + 1._wp )
431               ENDIF
432               
433               ! stress at T points (zkt/=0 if landfast)
434               !zs1(ji,jj) = ( zs1(ji,jj) * zalph1 + zp_delt(ji,jj) * ( zdiv * (1._wp + zkt) - zdelta * (1._wp - zkt) ) ) * z1_alph1
435               !zs2(ji,jj) = ( zs2(ji,jj) * zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 * (1._wp + zkt) ) ) * z1_alph2
436               zs1(ji,jj) = ( zs1(ji,jj) + stressptmp * zalph1  ) * z1_alph1
437               zs2(ji,jj) = ( zs2(ji,jj) + stressmtmp * zalph1  ) * z1_alph1
438               !zs1(ji,jj) = ( stressptmp * zs1(ji,jj) + zalph1  ) * z1_alph1
439               !zs2(ji,jj) = ( stressmtmp * zs2(ji,jj) + zalph1  ) * z1_alph1
440               !WRITE(numout,*) 'stress output', ji, jj, zs1(ji,jj)
441            END DO
442         END DO
443         !CALL lbc_lnk( 'icedyn_rhg_eap', zp_delt, 'T', 1. )
444         CALL lbc_lnk( 'icedyn_rhg_eap', stress12tmp, 'T', 1. )
445
446         DO jj = 1, jpjm1
447            DO ji = 1, jpim1
448               ! stress12tmp at F points
449               stress12tmpF = ( stress12tmp(ji,jj+1) * e1e2t(ji,jj+1) + stress12tmp(ji+1,jj+1) * e1e2t(ji+1,jj+1)  &
450                  &   + stress12tmp(ji,jj  ) * e1e2t(ji,jj  ) + stress12tmp(ji+1,jj  ) * e1e2t(ji+1,jj  )  &
451                  &   ) * 0.25_wp * r1_e1e2f(ji,jj)
452
453               ! alpha & beta for aEVP
454               IF( ln_aEVP ) THEN
455                  zalph1   = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
456                  z1_alph1 = 1._wp / ( zalph1 + 1._wp )
457               ENDIF
458               
459               ! stress at F points (zkt/=0 if landfast)
460               !zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 * (1._wp + zkt) ) * 0.5_wp ) * z1_alph2
461               zs12(ji,jj) = ( zs12(ji,jj) + stress12tmpF * zalph1  ) * z1_alph1
462               !zs12(ji,jj) = ( stress12tmpF * zs12(ji,jj) + zalph1  ) * z1_alph1
463
464            END DO
465         END DO
466
467         ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- !
468         DO jj = 2, jpjm1
469            DO ji = fs_2, fs_jpim1               
470               !                   !--- U points
471               zfU(ji,jj) = 0.5_wp * ( ( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj)                                             &
472                  &                  + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj)    &
473                  &                    ) * r1_e2u(ji,jj)                                                                      &
474                  &                  + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1)  &
475                  &                    ) * 2._wp * r1_e1u(ji,jj)                                                              &
476                  &                  ) * r1_e1e2u(ji,jj)
477               !
478               !                !--- V points
479               zfV(ji,jj) = 0.5_wp * ( ( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj)                                             &
480                  &                  - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj)    &
481                  &                    ) * r1_e1v(ji,jj)                                                                      &
482                  &                  + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj)  &
483                  &                    ) * 2._wp * r1_e2v(ji,jj)                                                              &
484                  &                  ) * r1_e1e2v(ji,jj)
485               !
486               !                !--- ice currents at U-V point
487               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)
488               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)
489               !
490            END DO
491         END DO
492         !
493         ! --- Computation of ice velocity --- !
494         !  Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta vary as in Kimmritz 2016 & 2017
495         !  Bouillon et al. 2009 (eq 34-35) => stable
496         IF( MOD(jter,2) == 0 ) THEN ! even iterations
497            !
498            DO jj = 2, jpjm1
499               DO ji = fs_2, fs_jpim1
500                  !                 !--- tau_io/(v_oce - v_ice)
501                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
502                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
503                  !                 !--- Ocean-to-Ice stress
504                  ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
505                  !
506                  !                 !--- tau_bottom/v_ice
507                  zvel  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
508                  zTauB = ztauy_base(ji,jj) / zvel
509                  !                 !--- OceanBottom-to-Ice stress
510                  ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
511                  !
512                  !                 !--- Coriolis at V-points (energy conserving formulation)
513                  zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
514                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
515                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
516                  !
517                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
518                  zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
519                  !
520                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
521                  !                                         1 = sliding friction : TauB < RHS
522                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
523                  !
524                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
525                     v_ice(ji,jj) = ( (          rswitch * ( zmV_t(ji,jj) * ( zbeta(ji,jj) * v_ice(ji,jj) + v_ice_b(ji,jj) )       & ! previous velocity
526                        &                                  + zRHS + zTauO * v_ice(ji,jj) )                                         & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
527                        &                                  / MAX( zepsi, zmV_t(ji,jj) * ( zbeta(ji,jj) + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
528                        &               + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )           & ! static friction => slow decrease to v=0
529                        &             ) * 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
530                        &           )   * zmsk00y(ji,jj)
531                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
532                     v_ice(ji,jj) = ( (           rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                                       & ! previous velocity
533                        &                                     + zRHS + zTauO * v_ice(ji,jj) )                                      & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
534                        &                                     / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                         & ! m/dt + tau_io(only ice part) + landfast
535                        &                + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )          & ! static friction => slow decrease to v=0
536                        &              ) * 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
537                        &            )   * zmsk00y(ji,jj)
538                  ENDIF
539!extra code for test case boundary conditions
540                  IF (mjg(jj)<25 .or. mjg(jj)>975 .or. mig(ji)<25 .or. mig(ji)>975) THEN
541                     v_ice(ji,jj) = invw*(ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj))
542                  END IF
543
544               END DO
545            END DO
546            CALL lbc_lnk( 'icedyn_rhg_eap', v_ice, 'V', -1. )
547            !
548#if defined key_agrif
549!!            CALL agrif_interp_ice( 'V', jter, nn_nevp )
550            CALL agrif_interp_ice( 'V' )
551#endif
552            IF( ln_bdy )   CALL bdy_ice_dyn( 'V' )
553            !
554            DO jj = 2, jpjm1
555               DO ji = fs_2, fs_jpim1         
556                  !                 !--- tau_io/(u_oce - u_ice)
557                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
558                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
559                  !                 !--- Ocean-to-Ice stress
560                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
561                  !
562                  !                 !--- tau_bottom/u_ice
563                  zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
564                  zTauB = ztaux_base(ji,jj) / zvel
565                  !                 !--- OceanBottom-to-Ice stress
566                  ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
567                  !
568                  !                 !--- Coriolis at U-points (energy conserving formulation)
569                  zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
570                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
571                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
572                  !
573                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
574                  zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
575                  !
576                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
577                  !                                         1 = sliding friction : TauB < RHS
578                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
579                  !
580                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
581                     u_ice(ji,jj) = ( (          rswitch * ( zmU_t(ji,jj) * ( zbeta(ji,jj) * u_ice(ji,jj) + u_ice_b(ji,jj) )       & ! previous velocity
582                        &                                  + zRHS + zTauO * u_ice(ji,jj) )                                         & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
583                        &                                  / MAX( zepsi, zmU_t(ji,jj) * ( zbeta(ji,jj) + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
584                        &               + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )           & ! static friction => slow decrease to v=0
585                        &             ) * 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
586                        &           )   * zmsk00x(ji,jj)
587                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
588                     u_ice(ji,jj) = ( (           rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                                       & ! previous velocity
589                        &                                     + zRHS + zTauO * u_ice(ji,jj) )                                      & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
590                        &                                     / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                         & ! m/dt + tau_io(only ice part) + landfast
591                        &                + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )          & ! static friction => slow decrease to v=0
592                        &              ) * 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
593                        &            )   * zmsk00x(ji,jj)
594                  ENDIF
595!extra code for test case boundary conditions
596                  IF (mjg(jj)<25 .or. mjg(jj)>975 .or. mig(ji)<25 .or. mig(ji)>975) THEN
597                     u_ice(ji,jj) = invw*(ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj))
598                  END IF
599
600               END DO
601            END DO
602            CALL lbc_lnk( 'icedyn_rhg_eap', u_ice, 'U', -1. )
603            !
604#if defined key_agrif
605!!            CALL agrif_interp_ice( 'U', jter, nn_nevp )
606            CALL agrif_interp_ice( 'U' )
607#endif
608            IF( ln_bdy )   CALL bdy_ice_dyn( 'U' )
609            !
610         ELSE ! odd iterations
611            !
612            DO jj = 2, jpjm1
613               DO ji = fs_2, fs_jpim1
614                  !                 !--- tau_io/(u_oce - u_ice)
615                  zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
616                     &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
617                  !                 !--- Ocean-to-Ice stress
618                  ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
619                  !
620                  !                 !--- tau_bottom/u_ice
621                  zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
622                  zTauB = ztaux_base(ji,jj) / zvel
623                  !                 !--- OceanBottom-to-Ice stress
624                  ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
625                  !
626                  !                 !--- Coriolis at U-points (energy conserving formulation)
627                  zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
628                     &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
629                     &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
630                  !
631                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
632                  zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
633                  !
634                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
635                  !                                         1 = sliding friction : TauB < RHS
636                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
637                  !
638                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
639                     u_ice(ji,jj) = ( (          rswitch * ( zmU_t(ji,jj) * ( zbeta(ji,jj) * u_ice(ji,jj) + u_ice_b(ji,jj) )       & ! previous velocity
640                        &                                  + zRHS + zTauO * u_ice(ji,jj) )                                         & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
641                        &                                  / MAX( zepsi, zmU_t(ji,jj) * ( zbeta(ji,jj) + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
642                        &               + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )           & ! static friction => slow decrease to v=0
643                        &             ) * 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
644                        &           )   * zmsk00x(ji,jj)
645                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
646                     u_ice(ji,jj) = ( (           rswitch   * ( zmU_t(ji,jj)  * u_ice(ji,jj)                                       & ! previous velocity
647                        &                                     + zRHS + zTauO * u_ice(ji,jj) )                                      & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
648                        &                                     / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                         & ! m/dt + tau_io(only ice part) + landfast
649                        &                + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )          & ! static friction => slow decrease to v=0
650                        &              ) * 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
651                        &            )   * zmsk00x(ji,jj)
652                  ENDIF
653!extra code for test case boundary conditions
654                  IF (mjg(jj)<25 .or. mjg(jj)>975 .or. mig(ji)<25 .or. mig(ji)>975) THEN
655                     u_ice(ji,jj) = invw*(ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj))
656                  END IF
657
658               END DO
659            END DO
660            CALL lbc_lnk( 'icedyn_rhg_eap', u_ice, 'U', -1. )
661            !
662#if defined key_agrif
663!!            CALL agrif_interp_ice( 'U', jter, nn_nevp )
664            CALL agrif_interp_ice( 'U' )
665#endif
666            IF( ln_bdy )   CALL bdy_ice_dyn( 'U' )
667            !
668            DO jj = 2, jpjm1
669               DO ji = fs_2, fs_jpim1
670                  !                 !--- tau_io/(v_oce - v_ice)
671                  zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
672                     &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
673                  !                 !--- Ocean-to-Ice stress
674                  ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
675                  !
676                  !                 !--- tau_bottom/v_ice
677                  zvel  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
678                  zTauB = ztauy_base(ji,jj) / zvel
679                  !                 !--- OceanBottom-to-Ice stress
680                  ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
681                  !
682                  !                 !--- Coriolis at v-points (energy conserving formulation)
683                  zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
684                     &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
685                     &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
686                  !
687                  !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
688                  zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
689                  !
690                  !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
691                  !                                         1 = sliding friction : TauB < RHS
692                  rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
693                  !
694                  IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
695                     v_ice(ji,jj) = ( (          rswitch * ( zmV_t(ji,jj) * ( zbeta(ji,jj) * v_ice(ji,jj) + v_ice_b(ji,jj) )       & ! previous velocity
696                        &                                  + zRHS + zTauO * v_ice(ji,jj) )                                         & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
697                        &                                  / MAX( zepsi, zmV_t(ji,jj) * ( zbeta(ji,jj) + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
698                        &               + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )           & ! static friction => slow decrease to v=0
699                        &             ) * 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
700                        &           )   * zmsk00y(ji,jj)
701                  ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
702                     v_ice(ji,jj) = ( (           rswitch   * ( zmV_t(ji,jj)  * v_ice(ji,jj)                                       & ! previous velocity
703                        &                                     + zRHS + zTauO * v_ice(ji,jj) )                                      & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
704                        &                                     / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                         & ! m/dt + tau_io(only ice part) + landfast
705                        &                + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lfrelax )          & ! static friction => slow decrease to v=0
706                        &              ) * 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
707                        &            )   * zmsk00y(ji,jj)
708                  ENDIF
709!extra code for test case boundary conditions
710                  IF (mjg(jj)<25 .or. mjg(jj)>975 .or. mig(ji)<25 .or. mig(ji)>975) THEN
711                     v_ice(ji,jj) = invw*(ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj))
712                  END IF
713
714               END DO
715            END DO
716            CALL lbc_lnk( 'icedyn_rhg_eap', v_ice, 'V', -1. )
717            !
718#if defined key_agrif
719!!            CALL agrif_interp_ice( 'V', jter, nn_nevp )
720            CALL agrif_interp_ice( 'V' )
721#endif
722            IF( ln_bdy )   CALL bdy_ice_dyn( 'V' )
723            !
724         ENDIF
725
726!!$         IF(ln_ctl) THEN   ! Convergence test
727!!$            DO jj = 2 , jpjm1
728!!$               zresr(:,jj) = MAX( ABS( u_ice(:,jj) - zu_ice(:,jj) ), ABS( v_ice(:,jj) - zv_ice(:,jj) ) )
729!!$            END DO
730!!$            zresm = MAXVAL( zresr( 1:jpi, 2:jpjm1 ) )
731!!$            CALL mpp_max( 'icedyn_rhg_eap', zresm )   ! max over the global domain
732!!$         ENDIF
733         !
734         !                                                ! ==================== !
735      END DO                                              !  end loop over jter  !
736      !                                                   ! ==================== !
737      !
738      !------------------------------------------------------------------------------!
739      ! 4) Recompute delta, shear and div (inputs for mechanical redistribution)
740      !------------------------------------------------------------------------------!
741      DO jj = 1, jpjm1
742         DO ji = 1, jpim1
743
744            ! shear at F points
745            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)   &
746               &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
747               &         ) * r1_e1e2f(ji,jj) * zfmask(ji,jj)
748
749         END DO
750      END DO           
751     
752      DO jj = 2, jpjm1
753         DO ji = 2, jpim1 ! no vector loop
754           
755            ! tension**2 at T points
756            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)   &
757               &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
758               &   ) * r1_e1e2t(ji,jj)
759            zdt2 = zdt * zdt
760           
761            ! shear**2 at T points (doc eq. A16)
762            zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
763               &   + 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)  &
764               &   ) * 0.25_wp * r1_e1e2t(ji,jj)
765           
766            ! shear at T points
767            pshear_i(ji,jj) = SQRT( zdt2 + zds2 )
768
769            ! divergence at T points
770            pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
771               &             + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
772               &             ) * r1_e1e2t(ji,jj)
773           
774            ! delta at T points
775            zdelta         = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) 
776            rswitch        = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zdelta ) ) ! 0 if delta=0
777            pdelta_i(ji,jj) = zdelta + rn_creepl * rswitch
778
779         END DO
780      END DO
781      CALL lbc_lnk_multi( 'icedyn_rhg_eap', pshear_i, 'T', 1., pdivu_i, 'T', 1., pdelta_i, 'T', 1. )
782      !CALL lbc_lnk_multi( 'icedyn_rhg_eap', pshear_i, 'T', 1., pdivu_i, 'T', 1. )
783     
784      ! --- Store the stress tensor for the next time step --- !
785      CALL lbc_lnk_multi( 'icedyn_rhg_eap', zs1, 'T', 1., zs2, 'T', 1., zs12, 'F', 1. )
786      pstress1_i (:,:) = zs1 (:,:)
787      pstress2_i (:,:) = zs2 (:,:)
788      pstress12_i(:,:) = zs12(:,:)
789      !
790
791      !------------------------------------------------------------------------------!
792      ! 5) diagnostics
793      !------------------------------------------------------------------------------!
794      DO jj = 1, jpj
795         DO ji = 1, jpi
796            zmsk00(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) ! 1 if ice, 0 if no ice
797         END DO
798      END DO
799
800      ! --- ice-ocean, ice-atm. & ice-oceanbottom(landfast) stresses --- !
801      IF(  iom_use('utau_oi') .OR. iom_use('vtau_oi') .OR. iom_use('utau_ai') .OR. iom_use('vtau_ai') .OR. &
802         & iom_use('utau_bi') .OR. iom_use('vtau_bi') ) THEN
803         !
804         CALL lbc_lnk_multi( 'icedyn_rhg_eap', ztaux_oi, 'U', -1., ztauy_oi, 'V', -1., ztaux_ai, 'U', -1., ztauy_ai, 'V', -1., &
805            &                                  ztaux_bi, 'U', -1., ztauy_bi, 'V', -1. )
806         !
807         CALL iom_put( 'utau_oi' , ztaux_oi * zmsk00 )
808         CALL iom_put( 'vtau_oi' , ztauy_oi * zmsk00 )
809         CALL iom_put( 'utau_ai' , ztaux_ai * zmsk00 )
810         CALL iom_put( 'vtau_ai' , ztauy_ai * zmsk00 )
811         CALL iom_put( 'utau_bi' , ztaux_bi * zmsk00 )
812         CALL iom_put( 'vtau_bi' , ztauy_bi * zmsk00 )
813      ENDIF
814       
815      ! --- divergence, shear and strength --- !
816      IF( iom_use('icediv') )   CALL iom_put( 'icediv' , pdivu_i  * zmsk00 )   ! divergence
817      IF( iom_use('iceshe') )   CALL iom_put( 'iceshe' , pshear_i * zmsk00 )   ! shear
818      IF( iom_use('icestr') )   CALL iom_put( 'icestr' , strength * zmsk00 )   ! strength
819
820      ! --- stress tensor --- !
821      IF( iom_use('isig1') .OR. iom_use('isig2') .OR. iom_use('isig3') .OR. iom_use('normstr') .OR. iom_use('sheastr') ) THEN
822         !
823         ALLOCATE( zsig1(jpi,jpj) , zsig2(jpi,jpj) , zsig3(jpi,jpj) )
824         !         
825         DO jj = 2, jpjm1
826            DO ji = 2, jpim1
827               zdum1 = ( zmsk00(ji-1,jj) * pstress12_i(ji-1,jj) + zmsk00(ji  ,jj-1) * pstress12_i(ji  ,jj-1) +  &  ! stress12_i at T-point
828                  &      zmsk00(ji  ,jj) * pstress12_i(ji  ,jj) + zmsk00(ji-1,jj-1) * pstress12_i(ji-1,jj-1) )  &
829                  &    / MAX( 1._wp, zmsk00(ji-1,jj) + zmsk00(ji,jj-1) + zmsk00(ji,jj) + zmsk00(ji-1,jj-1) )
830
831               zshear = SQRT( pstress2_i(ji,jj) * pstress2_i(ji,jj) + 4._wp * zdum1 * zdum1 ) ! shear stress 
832
833               zdum2 = zmsk00(ji,jj) / MAX( 1._wp, strength(ji,jj) )
834
835!!               zsig1(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) + zshear ) ! principal stress (y-direction, see Hunke & Dukowicz 2002)
836!!               zsig2(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) - zshear ) ! principal stress (x-direction, see Hunke & Dukowicz 2002)
837!!               zsig3(ji,jj) = zdum2**2 * ( ( pstress1_i(ji,jj) + strength(ji,jj) )**2 + ( rn_ecc * zshear )**2 ) ! quadratic relation linking compressive stress to shear stress
838!!                                                                                                               ! (scheme converges if this value is ~1, see Bouillon et al 2009 (eq. 11))
839               zsig1(ji,jj) = 0.5_wp * zdum2 * ( pstress1_i(ji,jj) )          ! compressive stress, see Bouillon et al. 2015
840               zsig2(ji,jj) = 0.5_wp * zdum2 * ( zshear )                     ! shear stress
841               zsig3(ji,jj) = zdum2**2 * ( ( pstress1_i(ji,jj) + strength(ji,jj) )**2 + ( rn_ecc * zshear )**2 )
842            END DO
843         END DO
844         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zsig1, 'T', 1., zsig2, 'T', 1., zsig3, 'T', 1. )
845         !
846         CALL iom_put( 'isig1' , zsig1 )
847         CALL iom_put( 'isig2' , zsig2 )
848         CALL iom_put( 'isig3' , zsig3 )
849         !
850         ! Stress tensor invariants (normal and shear stress N/m)
851         IF( iom_use('normstr') )   CALL iom_put( 'normstr' ,       ( zs1(:,:) + zs2(:,:) )                       * zmsk00(:,:) ) ! Normal stress
852         IF( iom_use('sheastr') )   CALL iom_put( 'sheastr' , SQRT( ( zs1(:,:) - zs2(:,:) )**2 + 4*zs12(:,:)**2 ) * zmsk00(:,:) ) ! Shear stress
853
854         DEALLOCATE( zsig1 , zsig2 , zsig3 )
855      ENDIF
856
857      ! --- anisotropy tensor --- !
858      IF( iom_use('aniso') ) THEN
859         CALL iom_put( 'aniso' , aniso_11 * zmsk00 )
860      ENDIF
861     
862      ! --- SIMIP --- !
863      IF(  iom_use('dssh_dx') .OR. iom_use('dssh_dy') .OR. &
864         & iom_use('corstrx') .OR. iom_use('corstry') .OR. iom_use('intstrx') .OR. iom_use('intstry') ) THEN
865         !
866         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zspgU, 'U', -1., zspgV, 'V', -1., &
867            &                                  zCorU, 'U', -1., zCorV, 'V', -1., zfU, 'U', -1., zfV, 'V', -1. )
868
869         CALL iom_put( 'dssh_dx' , zspgU * zmsk00 )   ! Sea-surface tilt term in force balance (x)
870         CALL iom_put( 'dssh_dy' , zspgV * zmsk00 )   ! Sea-surface tilt term in force balance (y)
871         CALL iom_put( 'corstrx' , zCorU * zmsk00 )   ! Coriolis force term in force balance (x)
872         CALL iom_put( 'corstry' , zCorV * zmsk00 )   ! Coriolis force term in force balance (y)
873         CALL iom_put( 'intstrx' , zfU   * zmsk00 )   ! Internal force term in force balance (x)
874         CALL iom_put( 'intstry' , zfV   * zmsk00 )   ! Internal force term in force balance (y)
875      ENDIF
876
877      IF(  iom_use('xmtrpice') .OR. iom_use('ymtrpice') .OR. &
878         & iom_use('xmtrpsnw') .OR. iom_use('ymtrpsnw') .OR. iom_use('xatrp') .OR. iom_use('yatrp') ) THEN
879         !
880         ALLOCATE( zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) , &
881            &      zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp(jpi,jpj) , zdiag_yatrp(jpi,jpj) )
882         !
883         DO jj = 2, jpjm1
884            DO ji = 2, jpim1
885               ! 2D ice mass, snow mass, area transport arrays (X, Y)
886               zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zmsk00(ji,jj)
887               zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zmsk00(ji,jj)
888
889               zdiag_xmtrp_ice(ji,jj) = rhoi * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component
890               zdiag_ymtrp_ice(ji,jj) = rhoi * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) !        ''           Y-   ''
891
892               zdiag_xmtrp_snw(ji,jj) = rhos * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
893               zdiag_ymtrp_snw(ji,jj) = rhos * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) !          ''          Y-   ''
894
895               zdiag_xatrp(ji,jj)     = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) )        ! area transport,      X-component
896               zdiag_yatrp(ji,jj)     = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) )        !        ''            Y-   ''
897
898            END DO
899         END DO
900
901         CALL lbc_lnk_multi( 'icedyn_rhg_eap', zdiag_xmtrp_ice, 'U', -1., zdiag_ymtrp_ice, 'V', -1., &
902            &                                  zdiag_xmtrp_snw, 'U', -1., zdiag_ymtrp_snw, 'V', -1., &
903            &                                  zdiag_xatrp    , 'U', -1., zdiag_yatrp    , 'V', -1. )
904
905         CALL iom_put( 'xmtrpice' , zdiag_xmtrp_ice )   ! X-component of sea-ice mass transport (kg/s)
906         CALL iom_put( 'ymtrpice' , zdiag_ymtrp_ice )   ! Y-component of sea-ice mass transport
907         CALL iom_put( 'xmtrpsnw' , zdiag_xmtrp_snw )   ! X-component of snow mass transport (kg/s)
908         CALL iom_put( 'ymtrpsnw' , zdiag_ymtrp_snw )   ! Y-component of snow mass transport
909         CALL iom_put( 'xatrp'    , zdiag_xatrp     )   ! X-component of ice area transport
910         CALL iom_put( 'yatrp'    , zdiag_yatrp     )   ! Y-component of ice area transport
911
912         DEALLOCATE( zdiag_xmtrp_ice , zdiag_ymtrp_ice , &
913            &        zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp , zdiag_yatrp )
914
915      ENDIF
916      !
917   END SUBROUTINE ice_dyn_rhg_eap
918
919   SUBROUTINE update_stress_rdg (ksub, ndte, phi, divu, tension, &
920                                 shear, a11, a12, &
921                                 stressp,  stressm, &
922                                 stress12, strength, &
923                                 alphar, alphas)
924      !!---------------------------------------------------------------------
925      !!                   ***  SUBROUTINE rhg_eap_rst  ***
926      !! Updates the stress depending on values of strain rate and structure
927      !! tensor and for ksub=ndte it computes closing and sliding rate
928      !!---------------------------------------------------------------------
929      INTEGER,  INTENT(in   ) ::   ksub, ndte
930      REAL(wp), INTENT(in   ) ::   phi, strength
931      REAL(wp), INTENT(in   ) ::   divu, tension, shear
932      REAL(wp), INTENT(in   ) ::   a11, a12
933      REAL(wp), INTENT(  out) ::   stressp, stressm, stress12
934      REAL(wp), INTENT(  out) ::   alphar, alphas
935
936      INTEGER  ::   kx ,ky, ka
937
938      REAL(wp) ::   stemp11r, stemp12r, stemp22r   
939      REAL(wp) ::   stemp11s, stemp12s, stemp22s 
940      REAL(wp) ::   a22, Qd11Qd11, Qd11Qd12, Qd12Qd12 
941      REAL(wp) ::   Q11Q11, Q11Q12, Q12Q12 
942      REAL(wp) ::   dtemp11, dtemp12, dtemp22 
943      REAL(wp) ::   rotstemp11r, rotstemp12r, rotstemp22r   
944      REAL(wp) ::   rotstemp11s, rotstemp12s, rotstemp22s
945      REAL(wp) ::   sig11, sig12, sig22 
946      REAL(wp) ::   sgprm11, sgprm12, sgprm22 
947      REAL(wp) ::   invstressconviso 
948      REAL(wp) ::   Angle_denom_gamma,  Angle_denom_alpha 
949      REAL(wp) ::   Tany_1, Tany_2 
950      REAL(wp) ::   x, y, dx, dy, da, kxw, kyw, kaw
951      REAL(wp) ::   invdx, invdy, invda   
952      REAL(wp) ::   dtemp1, dtemp2, atempprime, a, invsin
953
954      REAL(wp), PARAMETER ::   kfriction = 0.45_wp
955      !!---------------------------------------------------------------------
956      ! Factor to maintain the same stress as in EVP (see Section 3)
957      ! Can be set to 1 otherwise
958      invstressconviso = 1._wp/(1._wp+kfriction*kfriction)
959 
960      invsin = 1._wp/sin(2._wp*phi) * invstressconviso 
961      !now uses phi as set in higher code
962       
963      ! compute eigenvalues, eigenvectors and angles for structure tensor, strain
964      ! rates
965
966      ! 1) structure tensor
967      a22 = 1._wp-a11
968      Q11Q11 = 1._wp
969      Q12Q12 = rsmall
970      Q11Q12 = rsmall
971 
972      ! gamma: angle between general coordiantes and principal axis of A
973      ! here Tan2gamma = 2 a12 / (a11 - a22)
974      if((ABS(a11 - a22) > rsmall).or.(ABS(a12) > rsmall)) then     
975         Angle_denom_gamma = 1._wp/sqrt( ( a11 - a22 )*( a11 - a22) + &
976                             4._wp*a12*a12 )
977   
978         Q11Q11 = 0.5_wp + ( a11 - a22 )*0.5_wp*Angle_denom_gamma   !Cos^2
979         Q12Q12 = 0.5_wp - ( a11 - a22 )*0.5_wp*Angle_denom_gamma   !Sin^2
980         Q11Q12 = a12*Angle_denom_gamma                     !CosSin
981      endif
982 
983      ! rotation Q*atemp*Q^T
984      atempprime = Q11Q11*a11 + 2.0_wp*Q11Q12*a12 + Q12Q12*a22 
985     
986      ! make first principal value the largest
987      atempprime = max(atempprime, 1.0_wp - atempprime)
988 
989      ! 2) strain rate
990      dtemp11 = 0.5_wp*(divu + tension)
991      dtemp12 = shear*0.5_wp
992      dtemp22 = 0.5_wp*(divu - tension)
993
994      !WRITE(numout,*) 'div, tension, shear inside loop', ksub, divu, tension, shear
995
996      ! here Tan2alpha = 2 dtemp12 / (dtemp11 - dtemp22)
997
998      Qd11Qd11 = 1.0_wp
999      Qd12Qd12 = rsmall
1000      Qd11Qd12 = rsmall
1001
1002      if((ABS( dtemp11 - dtemp22) > rsmall).or. (ABS(dtemp12) > rsmall)) then
1003
1004         Angle_denom_alpha = 1.0_wp/sqrt( ( dtemp11 - dtemp22 )* &
1005                             ( dtemp11 - dtemp22 ) + 4.0_wp*dtemp12*dtemp12)
1006
1007         Qd11Qd11 = 0.5_wp + ( dtemp11 - dtemp22 )*0.5_wp*Angle_denom_alpha !Cos^2
1008         Qd12Qd12 = 0.5_wp - ( dtemp11 - dtemp22 )*0.5_wp*Angle_denom_alpha !Sin^2
1009         Qd11Qd12 = dtemp12*Angle_denom_alpha !CosSin
1010      endif
1011
1012      dtemp1 = Qd11Qd11*dtemp11 + 2.0_wp*Qd11Qd12*dtemp12 + Qd12Qd12*dtemp22
1013      dtemp2 = Qd12Qd12*dtemp11 - 2.0_wp*Qd11Qd12*dtemp12 + Qd11Qd11*dtemp22
1014      ! In cos and sin values
1015      x = 0._wp
1016      if ((ABS(dtemp1) > rsmall).or.(ABS(dtemp2) > rsmall)) then
1017         x = atan2(dtemp2,dtemp1)
1018      endif     
1019               
1020      ! to ensure the angle lies between pi/4 and 9 pi/4
1021      if (x < rpi*0.25_wp) x = x + rpi*2.0_wp
1022     
1023      ! y: angle between major principal axis of strain rate and structure
1024      ! tensor
1025      ! y = gamma - alpha
1026      ! Expressesed componently with
1027      ! Tany = (Singamma*Cosgamma - Sinalpha*Cosgamma)/(Cos^2gamma - Sin^alpha)
1028     
1029      Tany_1 = Q11Q12 - Qd11Qd12
1030      Tany_2 = Q11Q11 - Qd12Qd12
1031     
1032      y = 0._wp
1033     
1034      if ((ABS(Tany_1) > rsmall).or.(ABS(Tany_2) > rsmall)) then
1035         y = atan2(Tany_1,Tany_2)
1036      endif
1037     
1038      ! to make sure y is between 0 and pi
1039      if (y > rpi) y = y - rpi
1040      if (y < 0)  y = y + rpi
1041
1042      ! 3) update anisotropic stress tensor
1043      dx   = rpi/real(nx_yield-1,kind=wp)
1044      dy   = rpi/real(ny_yield-1,kind=wp)
1045      da   = 0.5_wp/real(na_yield-1,kind=wp)
1046      invdx = 1._wp/dx
1047      invdy = 1._wp/dy
1048      invda = 1._wp/da
1049
1050      ! % need 8 coords and 8 weights
1051      ! % range in kx
1052      kx  = int((x-rpi*0.25_wp-rpi)*invdx) + 1
1053      kxw = kx - (x-rpi*0.25_wp-rpi)*invdx 
1054     
1055      ky  = int(y*invdy) + 1
1056      kyw = ky - y*invdy 
1057     
1058      ka  = int((atempprime-0.5_wp)*invda) + 1
1059      kaw = ka - (atempprime-0.5_wp)*invda
1060     
1061      ! % Determine sigma_r(A1,Zeta,y) and sigma_s (see Section A1 of Tsamados 2013)
1062      stemp11r =  kxw   * kyw         * kaw         * s11r(kx  ,ky  ,ka  ) &
1063        & + (1._wp-kxw) * kyw         * kaw         * s11r(kx+1,ky  ,ka  ) &
1064        & + kxw         * (1._wp-kyw) * kaw         * s11r(kx  ,ky+1,ka  ) &
1065        & + kxw         * kyw         * (1._wp-kaw) * s11r(kx  ,ky  ,ka+1) &
1066        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s11r(kx+1,ky+1,ka  ) &
1067        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s11r(kx+1,ky  ,ka+1) &
1068        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s11r(kx  ,ky+1,ka+1) &
1069        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s11r(kx+1,ky+1,ka+1)
1070      stemp12r =  kxw   * kyw         * kaw         * s12r(kx  ,ky  ,ka  ) &
1071        & + (1._wp-kxw) * kyw         * kaw         * s12r(kx+1,ky  ,ka  ) &
1072        & + kxw         * (1._wp-kyw) * kaw         * s12r(kx  ,ky+1,ka  ) &
1073        & + kxw         * kyw         * (1._wp-kaw) * s12r(kx  ,ky  ,ka+1) &
1074        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s12r(kx+1,ky+1,ka  ) &
1075        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s12r(kx+1,ky  ,ka+1) &
1076        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s12r(kx  ,ky+1,ka+1) &
1077        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s12r(kx+1,ky+1,ka+1)
1078      stemp22r = kxw    * kyw         * kaw         * s22r(kx  ,ky  ,ka  ) &
1079        & + (1._wp-kxw) * kyw         * kaw         * s22r(kx+1,ky  ,ka  ) &
1080        & + kxw         * (1._wp-kyw) * kaw         * s22r(kx  ,ky+1,ka  ) &
1081        & + kxw         * kyw         * (1._wp-kaw) * s22r(kx  ,ky  ,ka+1) &
1082        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s22r(kx+1,ky+1,ka  ) &
1083        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s22r(kx+1,ky  ,ka+1) &
1084        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s22r(kx  ,ky+1,ka+1) &
1085        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s22r(kx+1,ky+1,ka+1)
1086     
1087      stemp11s =  kxw   * kyw         * kaw         * s11s(kx  ,ky  ,ka  ) &
1088        & + (1._wp-kxw) * kyw         * kaw         * s11s(kx+1,ky  ,ka  ) &
1089        & + kxw         * (1._wp-kyw) * kaw         * s11s(kx  ,ky+1,ka  ) &
1090        & + kxw         * kyw         * (1._wp-kaw) * s11s(kx  ,ky  ,ka+1) &
1091        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s11s(kx+1,ky+1,ka  ) &
1092        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s11s(kx+1,ky  ,ka+1) &
1093        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s11s(kx  ,ky+1,ka+1) &
1094        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s11s(kx+1,ky+1,ka+1)
1095      stemp12s =  kxw   * kyw         * kaw         * s12s(kx  ,ky  ,ka  ) &
1096        & + (1._wp-kxw) * kyw         * kaw         * s12s(kx+1,ky  ,ka  ) &
1097        & + kxw         * (1._wp-kyw) * kaw         * s12s(kx  ,ky+1,ka  ) &
1098        & + kxw         * kyw         * (1._wp-kaw) * s12s(kx  ,ky  ,ka+1) &
1099        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s12s(kx+1,ky+1,ka  ) &
1100        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s12s(kx+1,ky  ,ka+1) &
1101        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s12s(kx  ,ky+1,ka+1) &
1102        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s12s(kx+1,ky+1,ka+1)
1103      stemp22s =  kxw   * kyw         * kaw         * s22s(kx  ,ky  ,ka  ) &
1104        & + (1._wp-kxw) * kyw         * kaw         * s22s(kx+1,ky  ,ka  ) &
1105        & + kxw         * (1._wp-kyw) * kaw         * s22s(kx  ,ky+1,ka  ) &
1106        & + kxw         * kyw         * (1._wp-kaw) * s22s(kx  ,ky  ,ka+1) &
1107        & + (1._wp-kxw) * (1._wp-kyw) * kaw         * s22s(kx+1,ky+1,ka  ) &
1108        & + (1._wp-kxw) * kyw         * (1._wp-kaw) * s22s(kx+1,ky  ,ka+1) &
1109        & + kxw         * (1._wp-kyw) * (1._wp-kaw) * s22s(kx  ,ky+1,ka+1) &
1110        & + (1._wp-kxw) * (1._wp-kyw) * (1._wp-kaw) * s22s(kx+1,ky+1,ka+1)
1111                   
1112      ! Calculate mean ice stress over a collection of floes (Equation 3 in
1113      ! Tsamados 2013)
1114
1115      sig11  = strength*(stemp11r + kfriction*stemp11s) * invsin
1116      sig12  = strength*(stemp12r + kfriction*stemp12s) * invsin
1117      sig22  = strength*(stemp22r + kfriction*stemp22s) * invsin
1118
1119      !WRITE(numout,*) 'strength inside loop', strength
1120      !WRITE(numout,*) 'principal axis stress inside loop', sig11, sig12, sig22
1121
1122      ! Back - rotation of the stress from principal axes into general coordinates
1123
1124      ! Update stress
1125      sgprm11 = Q11Q11*sig11 + Q12Q12*sig22 -        2._wp*Q11Q12 *sig12
1126      sgprm12 = Q11Q12*sig11 - Q11Q12*sig22 + (Q11Q11 - Q12Q12)*sig12
1127      sgprm22 = Q12Q12*sig11 + Q11Q11*sig22 +        2._wp*Q11Q12 *sig12
1128
1129      stressp  = sgprm11 + sgprm22
1130      stress12 = sgprm12
1131      stressm  = sgprm11 - sgprm22
1132
1133      !WRITE(numout,*) 'stress output inside loop', ksub, stressp
1134
1135      ! Compute ridging and sliding functions in general coordinates
1136      ! (Equation 11 in Tsamados 2013)
1137      if (ksub == ndte) then
1138         rotstemp11r = Q11Q11*stemp11r - 2._wp*Q11Q12* stemp12r &
1139                     + Q12Q12*stemp22r
1140         rotstemp12r = Q11Q11*stemp12r +    Q11Q12*(stemp11r-stemp22r) &
1141                     - Q12Q12*stemp12r
1142         rotstemp22r = Q12Q12*stemp11r + 2._wp*Q11Q12* stemp12r & 
1143                     + Q11Q11*stemp22r
1144
1145         rotstemp11s = Q11Q11*stemp11s - 2._wp*Q11Q12* stemp12s &
1146                     + Q12Q12*stemp22s
1147         rotstemp12s = Q11Q11*stemp12s +    Q11Q12*(stemp11s-stemp22s) &
1148                     - Q12Q12*stemp12s
1149         rotstemp22s = Q12Q12*stemp11s + 2._wp*Q11Q12* stemp12s & 
1150                     + Q11Q11*stemp22s
1151
1152         alphar =  rotstemp11r*dtemp11 + 2._wp*rotstemp12r*dtemp12 &
1153                 + rotstemp22r*dtemp22
1154         alphas =  rotstemp11s*dtemp11 + 2._wp*rotstemp12s*dtemp12 &
1155                 + rotstemp22s*dtemp22
1156      endif
1157   END SUBROUTINE update_stress_rdg 
1158
1159   SUBROUTINE rhg_eap_rst( cdrw, kt )
1160      !!---------------------------------------------------------------------
1161      !!                   ***  ROUTINE rhg_eap_rst  ***
1162      !!                     
1163      !! ** Purpose :   Read or write RHG file in restart file
1164      !!     
1165      !! ** Method  :   use of IOM library
1166      !!----------------------------------------------------------------------
1167      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
1168      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
1169      !
1170      INTEGER  ::   iter                      ! local integer
1171      INTEGER  ::   id1, id2, id3, id4, id5   ! local integers
1172      INTEGER  ::   ix, iy, ip, iz, n, ia     ! local integers
1173   
1174      INTEGER, PARAMETER            ::    nz = 100
1175     
1176      REAL(wp) ::   ainit, xinit, yinit, pinit, zinit
1177      REAL(wp) ::   da, dx, dy, dp, dz, a1
1178
1179      REAL(wp), PARAMETER           ::   eps6 = 1.0e-6_wp
1180      !!----------------------------------------------------------------------
1181      !
1182      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
1183         !                                   ! ---------------
1184         IF( ln_rstart ) THEN                   !* Read the restart file
1185            !
1186            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
1187            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
1188            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
1189            id4 = iom_varid( numrir, 'aniso_11' , ldstop = .FALSE. )
1190            id5 = iom_varid( numrir, 'aniso_12' , ldstop = .FALSE. )
1191            !
1192            IF( MIN( id1, id2, id3, id4, id5 ) > 0 ) THEN      ! fields exist
1193               CALL iom_get( numrir, jpdom_autoglo, 'stress1_i' , stress1_i  )
1194               CALL iom_get( numrir, jpdom_autoglo, 'stress2_i' , stress2_i  )
1195               CALL iom_get( numrir, jpdom_autoglo, 'stress12_i', stress12_i )
1196               CALL iom_get( numrir, jpdom_autoglo, 'aniso_11' , aniso_11 )
1197               CALL iom_get( numrir, jpdom_autoglo, 'aniso_12' , aniso_12 )
1198            ELSE                                     ! start rheology from rest
1199               IF(lwp) WRITE(numout,*)
1200               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without rheology, set stresses to 0'
1201               stress1_i (:,:) = 0._wp
1202               stress2_i (:,:) = 0._wp
1203               stress12_i(:,:) = 0._wp
1204               aniso_11 (:,:) = 0.5_wp
1205               aniso_12 (:,:) = 0._wp
1206            ENDIF
1207         ELSE                                   !* Start from rest
1208            IF(lwp) WRITE(numout,*)
1209            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set stresses to 0'
1210            stress1_i (:,:) = 0._wp
1211            stress2_i (:,:) = 0._wp
1212            stress12_i(:,:) = 0._wp
1213            aniso_11 (:,:) = 0.5_wp
1214            aniso_12 (:,:) = 0._wp
1215         ENDIF
1216         !
1217
1218      da = 0.5_wp/real(na_yield-1,kind=wp)
1219      ainit = 0.5_wp - da
1220      dx = rpi/real(nx_yield-1,kind=wp)
1221      xinit = rpi + 0.25_wp*rpi - dx
1222      dz = rpi/real(nz,kind=wp)
1223      zinit = -rpi*0.5_wp
1224      dy = rpi/real(ny_yield-1,kind=wp)
1225      yinit = -dy
1226
1227      DO ia=1,na_yield
1228       DO ix=1,nx_yield
1229        DO iy=1,ny_yield
1230         s11r(ix,iy,ia) = 0._wp
1231         s12r(ix,iy,ia) = 0._wp
1232         s22r(ix,iy,ia) = 0._wp
1233         s11s(ix,iy,ia) = 0._wp
1234         s12s(ix,iy,ia) = 0._wp
1235         s22s(ix,iy,ia) = 0._wp
1236         IF (ia <= na_yield-1) then
1237          DO iz=1,nz
1238          s11r(ix,iy,ia) = s11r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1239           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1240           s11kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi)
1241          s12r(ix,iy,ia) = s12r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1242           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1243           s12kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi)
1244          s22r(ix,iy,ia) = s22r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1245           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1246           s22kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi) 
1247          s11s(ix,iy,ia) = s11s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1248           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1249           s11ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi)
1250          s12s(ix,iy,ia) = s12s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1251           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1252           s12ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi)
1253          s22s(ix,iy,ia) = s22s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
1254           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
1255           s22ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz,phi)*dz/sin(2._wp*phi)
1256          ENDDO
1257          IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1258          IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1259          IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1260          IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1261          IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1262          IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1263         ELSE
1264          s11r(ix,iy,ia) = 0.5_wp*s11kr(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1265          s12r(ix,iy,ia) = 0.5_wp*s12kr(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1266          s22r(ix,iy,ia) = 0.5_wp*s22kr(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1267          s11s(ix,iy,ia) = 0.5_wp*s11ks(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1268          s12s(ix,iy,ia) = 0.5_wp*s12ks(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1269          s22s(ix,iy,ia) = 0.5_wp*s22ks(xinit+ix*dx,yinit+iy*dy,0._wp,phi)/sin(2._wp*phi)
1270          IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
1271          IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
1272          IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
1273          IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
1274          IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
1275          IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
1276         ENDIF
1277        ENDDO
1278       ENDDO
1279      ENDDO
1280
1281      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1282         !                                   ! -------------------
1283         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
1284         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
1285         !
1286         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i  )
1287         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i  )
1288         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
1289         CALL iom_rstput( iter, nitrst, numriw, 'aniso_11' , aniso_11 )
1290         CALL iom_rstput( iter, nitrst, numriw, 'aniso_12' , aniso_12 )
1291         !
1292      ENDIF
1293      !
1294   END SUBROUTINE rhg_eap_rst
1295
1296
1297   FUNCTION w1(a)
1298      !!-------------------------------------------------------------------
1299      !! Function : w1 (see Gaussian function psi in Tsamados et al 2013)
1300      !!-------------------------------------------------------------------
1301      REAL(wp), INTENT(in   ) ::   a
1302      REAL(wp) ::   w1
1303      !!-------------------------------------------------------------------
1304   
1305      w1 = - 223.87569446_wp &
1306           + 2361.2198663_wp*a &
1307           - 10606.56079975_wp*a*a &
1308           + 26315.50025642_wp*a*a*a &
1309           - 38948.30444297_wp*a*a*a*a &
1310           + 34397.72407466_wp*a*a*a*a*a &
1311           - 16789.98003081_wp*a*a*a*a*a*a &
1312           + 3495.82839237_wp*a*a*a*a*a*a*a
1313   
1314   END FUNCTION w1
1315
1316
1317   FUNCTION w2(a)
1318      !!-------------------------------------------------------------------
1319      !! Function : w2 (see Gaussian function psi in Tsamados et al 2013)
1320      !!-------------------------------------------------------------------
1321      REAL(wp), INTENT(in   ) ::   a
1322      REAL(wp) ::   w2
1323      !!-------------------------------------------------------------------
1324   
1325      w2 = - 6670.68911883_wp &
1326           + 70222.33061536_wp*a &
1327           - 314871.71525448_wp*a*a &
1328           + 779570.02793492_wp*a*a*a &
1329           - 1151098.82436864_wp*a*a*a*a &
1330           + 1013896.59464498_wp*a*a*a*a*a &
1331           - 493379.44906738_wp*a*a*a*a*a*a &
1332           + 102356.551518_wp*a*a*a*a*a*a*a
1333   
1334   END FUNCTION w2
1335
1336   FUNCTION s11kr(x,y,z,phi) 
1337      !!-------------------------------------------------------------------
1338      !! Function : s11kr
1339      !!-------------------------------------------------------------------
1340      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1341   
1342      REAL(wp) ::   s11kr, p, pih
1343   
1344      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1345      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1346      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1347      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1348      REAL(wp) ::   d11, d12, d22
1349      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1350      REAL(wp) ::   Hen1t2, Hen2t1
1351      !!-------------------------------------------------------------------
1352   
1353      pih = 0.5_wp*rpi
1354      p = phi
1355   
1356      n1t2i11 = cos(z+pih-p) * cos(z+p)
1357      n1t2i12 = cos(z+pih-p) * sin(z+p)
1358      n1t2i21 = sin(z+pih-p) * cos(z+p)
1359      n1t2i22 = sin(z+pih-p) * sin(z+p)
1360      n2t1i11 = cos(z-pih+p) * cos(z-p)
1361      n2t1i12 = cos(z-pih+p) * sin(z-p)
1362      n2t1i21 = sin(z-pih+p) * cos(z-p)
1363      n2t1i22 = sin(z-pih+p) * sin(z-p)
1364      t1t2i11 = cos(z-p) * cos(z+p)
1365      t1t2i12 = cos(z-p) * sin(z+p)
1366      t1t2i21 = sin(z-p) * cos(z+p)
1367      t1t2i22 = sin(z-p) * sin(z+p)
1368      t2t1i11 = cos(z+p) * cos(z-p)
1369      t2t1i12 = cos(z+p) * sin(z-p)
1370      t2t1i21 = sin(z+p) * cos(z-p)
1371      t2t1i22 = sin(z+p) * sin(z-p)
1372   ! In expression of tensor d, with this formulatin d(x)=-d(x+pi)
1373   ! Solution, when diagonalizing always check sgn(a11-a22) if > then keep x else
1374   ! x=x-pi/2
1375      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1376      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1377      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1378      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1379      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1380      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1381   
1382      IF (-IIn1t2>=rsmall) then
1383      Hen1t2 = 1._wp
1384      ELSE
1385      Hen1t2 = 0._wp
1386      ENDIF
1387   
1388      IF (-IIn2t1>=rsmall) then
1389      Hen2t1 = 1._wp
1390      ELSE
1391      Hen2t1 = 0._wp
1392      ENDIF
1393   
1394      s11kr = (- Hen1t2 * n1t2i11 - Hen2t1 * n2t1i11)
1395
1396   END FUNCTION s11kr
1397
1398   FUNCTION s12kr(x,y,z,phi)
1399      !!-------------------------------------------------------------------
1400      !! Function : s12kr
1401      !!-------------------------------------------------------------------
1402      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1403
1404      REAL(wp) ::   s12kr, s12r0, s21r0, p, pih
1405
1406      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1407      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1408      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1409      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1410      REAL(wp) ::   d11, d12, d22
1411      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1412      REAL(wp) ::   Hen1t2, Hen2t1
1413      !!-------------------------------------------------------------------
1414      pih = 0.5_wp*rpi
1415      p = phi
1416   
1417      n1t2i11 = cos(z+pih-p) * cos(z+p)
1418      n1t2i12 = cos(z+pih-p) * sin(z+p)
1419      n1t2i21 = sin(z+pih-p) * cos(z+p)
1420      n1t2i22 = sin(z+pih-p) * sin(z+p)
1421      n2t1i11 = cos(z-pih+p) * cos(z-p)
1422      n2t1i12 = cos(z-pih+p) * sin(z-p)
1423      n2t1i21 = sin(z-pih+p) * cos(z-p)
1424      n2t1i22 = sin(z-pih+p) * sin(z-p)
1425      t1t2i11 = cos(z-p) * cos(z+p)
1426      t1t2i12 = cos(z-p) * sin(z+p)
1427      t1t2i21 = sin(z-p) * cos(z+p)
1428      t1t2i22 = sin(z-p) * sin(z+p)
1429      t2t1i11 = cos(z+p) * cos(z-p)
1430      t2t1i12 = cos(z+p) * sin(z-p)
1431      t2t1i21 = sin(z+p) * cos(z-p)
1432      t2t1i22 = sin(z+p) * sin(z-p)
1433      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1434      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1435      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1436      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1437      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1438      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1439   
1440      IF (-IIn1t2>=rsmall) then
1441      Hen1t2 = 1._wp
1442      ELSE
1443      Hen1t2 = 0._wp
1444      ENDIF
1445   
1446      IF (-IIn2t1>=rsmall) then
1447      Hen2t1 = 1._wp
1448      ELSE
1449      Hen2t1 = 0._wp
1450      ENDIF
1451   
1452      s12r0 = (- Hen1t2 * n1t2i12 - Hen2t1 * n2t1i12)
1453      s21r0 = (- Hen1t2 * n1t2i21 - Hen2t1 * n2t1i21)
1454      s12kr=0.5_wp*(s12r0+s21r0)
1455   
1456   END FUNCTION s12kr
1457
1458   FUNCTION s22kr(x,y,z,phi)
1459      !!-------------------------------------------------------------------
1460      !! Function : s22kr
1461      !!-------------------------------------------------------------------
1462      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1463
1464      REAL(wp) ::   s22kr, p, pih
1465
1466      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1467      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1468      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1469      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1470      REAL(wp) ::   d11, d12, d22
1471      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1472      REAL(wp) ::   Hen1t2, Hen2t1
1473      !!-------------------------------------------------------------------
1474      pih = 0.5_wp*rpi
1475      p = phi
1476
1477      n1t2i11 = cos(z+pih-p) * cos(z+p)
1478      n1t2i12 = cos(z+pih-p) * sin(z+p)
1479      n1t2i21 = sin(z+pih-p) * cos(z+p)
1480      n1t2i22 = sin(z+pih-p) * sin(z+p)
1481      n2t1i11 = cos(z-pih+p) * cos(z-p)
1482      n2t1i12 = cos(z-pih+p) * sin(z-p)
1483      n2t1i21 = sin(z-pih+p) * cos(z-p)
1484      n2t1i22 = sin(z-pih+p) * sin(z-p)
1485      t1t2i11 = cos(z-p) * cos(z+p)
1486      t1t2i12 = cos(z-p) * sin(z+p)
1487      t1t2i21 = sin(z-p) * cos(z+p)
1488      t1t2i22 = sin(z-p) * sin(z+p)
1489      t2t1i11 = cos(z+p) * cos(z-p)
1490      t2t1i12 = cos(z+p) * sin(z-p)
1491      t2t1i21 = sin(z+p) * cos(z-p)
1492      t2t1i22 = sin(z+p) * sin(z-p)
1493      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1494      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1495      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1496      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1497      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1498      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1499
1500      IF (-IIn1t2>=rsmall) then
1501      Hen1t2 = 1._wp
1502      ELSE
1503      Hen1t2 = 0._wp
1504      ENDIF
1505
1506      IF (-IIn2t1>=rsmall) then
1507      Hen2t1 = 1._wp
1508      ELSE
1509      Hen2t1 = 0._wp
1510      ENDIF
1511
1512      s22kr = (- Hen1t2 * n1t2i22 - Hen2t1 * n2t1i22)
1513
1514   END FUNCTION s22kr
1515
1516   FUNCTION s11ks(x,y,z,phi)
1517      !!-------------------------------------------------------------------
1518      !! Function : s11ks
1519      !!-------------------------------------------------------------------
1520      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1521
1522      REAL(wp) ::   s11ks, p, pih
1523
1524      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1525      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1526      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1527      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1528      REAL(wp) ::   d11, d12, d22
1529      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1530      REAL(wp) ::   Hen1t2, Hen2t1
1531      !!-------------------------------------------------------------------
1532      pih = 0.5_wp*rpi
1533      p = phi
1534
1535      n1t2i11 = cos(z+pih-p) * cos(z+p)
1536      n1t2i12 = cos(z+pih-p) * sin(z+p)
1537      n1t2i21 = sin(z+pih-p) * cos(z+p)
1538      n1t2i22 = sin(z+pih-p) * sin(z+p)
1539      n2t1i11 = cos(z-pih+p) * cos(z-p)
1540      n2t1i12 = cos(z-pih+p) * sin(z-p)
1541      n2t1i21 = sin(z-pih+p) * cos(z-p)
1542      n2t1i22 = sin(z-pih+p) * sin(z-p)
1543      t1t2i11 = cos(z-p) * cos(z+p)
1544      t1t2i12 = cos(z-p) * sin(z+p)
1545      t1t2i21 = sin(z-p) * cos(z+p)
1546      t1t2i22 = sin(z-p) * sin(z+p)
1547      t2t1i11 = cos(z+p) * cos(z-p)
1548      t2t1i12 = cos(z+p) * sin(z-p)
1549      t2t1i21 = sin(z+p) * cos(z-p)
1550      t2t1i22 = sin(z+p) * sin(z-p)
1551      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1552      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1553      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1554      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1555      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1556      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1557
1558      IF (-IIn1t2>=rsmall) then
1559      Hen1t2 = 1._wp
1560      ELSE
1561      Hen1t2 = 0._wp
1562      ENDIF
1563
1564      IF (-IIn2t1>=rsmall) then
1565      Hen2t1 = 1._wp
1566      ELSE
1567      Hen2t1 = 0._wp
1568      ENDIF
1569
1570      s11ks = sign(1._wp,IIt1t2+rsmall)*(Hen1t2 * t1t2i11 + Hen2t1 * t2t1i11)
1571
1572   END FUNCTION s11ks
1573
1574   FUNCTION s12ks(x,y,z,phi)
1575      !!-------------------------------------------------------------------
1576      !! Function : s12ks
1577      !!-------------------------------------------------------------------
1578      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1579
1580      REAL(wp) ::   s12ks, s12s0, s21s0, p, pih
1581
1582      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1583      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1584      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1585      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1586      REAL(wp) ::   d11, d12, d22
1587      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1588      REAL(wp) ::   Hen1t2, Hen2t1
1589      !!-------------------------------------------------------------------
1590      pih = 0.5_wp*rpi
1591      p =phi
1592
1593      n1t2i11 = cos(z+pih-p) * cos(z+p)
1594      n1t2i12 = cos(z+pih-p) * sin(z+p)
1595      n1t2i21 = sin(z+pih-p) * cos(z+p)
1596      n1t2i22 = sin(z+pih-p) * sin(z+p)
1597      n2t1i11 = cos(z-pih+p) * cos(z-p)
1598      n2t1i12 = cos(z-pih+p) * sin(z-p)
1599      n2t1i21 = sin(z-pih+p) * cos(z-p)
1600      n2t1i22 = sin(z-pih+p) * sin(z-p)
1601      t1t2i11 = cos(z-p) * cos(z+p)
1602      t1t2i12 = cos(z-p) * sin(z+p)
1603      t1t2i21 = sin(z-p) * cos(z+p)
1604      t1t2i22 = sin(z-p) * sin(z+p)
1605      t2t1i11 = cos(z+p) * cos(z-p)
1606      t2t1i12 = cos(z+p) * sin(z-p)
1607      t2t1i21 = sin(z+p) * cos(z-p)
1608      t2t1i22 = sin(z+p) * sin(z-p)
1609      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1610      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1611      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1612      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1613      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1614      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1615
1616      IF (-IIn1t2>=rsmall) then
1617      Hen1t2 = 1._wp
1618      ELSE
1619      Hen1t2 = 0._wp
1620      ENDIF
1621
1622      IF (-IIn2t1>=rsmall) then
1623      Hen2t1 = 1._wp
1624      ELSE
1625      Hen2t1 = 0._wp
1626      ENDIF
1627
1628      s12s0 = sign(1._wp,IIt1t2+rsmall)*(Hen1t2 * t1t2i12 + Hen2t1 * t2t1i12)
1629      s21s0 = sign(1._wp,IIt1t2+rsmall)*(Hen1t2 * t1t2i21 + Hen2t1 * t2t1i21)
1630      s12ks=0.5_wp*(s12s0+s21s0)
1631
1632   END FUNCTION s12ks
1633
1634   FUNCTION s22ks(x,y,z,phi) 
1635      !!-------------------------------------------------------------------
1636      !! Function : s22ks
1637      !!-------------------------------------------------------------------
1638      REAL(wp), INTENT(in   ) ::   x,y,z,phi
1639
1640      REAL(wp) ::   s22ks, p, pih
1641
1642      REAL(wp) ::   n1t2i11, n1t2i12, n1t2i21, n1t2i22
1643      REAL(wp) ::   n2t1i11, n2t1i12, n2t1i21, n2t1i22
1644      REAL(wp) ::   t1t2i11, t1t2i12, t1t2i21, t1t2i22
1645      REAL(wp) ::   t2t1i11, t2t1i12, t2t1i21, t2t1i22
1646      REAL(wp) ::   d11, d12, d22
1647      REAL(wp) ::   IIn1t2, IIn2t1, IIt1t2
1648      REAL(wp) ::   Hen1t2, Hen2t1
1649      !!-------------------------------------------------------------------
1650      pih = 0.5_wp*rpi
1651      p =phi
1652
1653      n1t2i11 = cos(z+pih-p) * cos(z+p)
1654      n1t2i12 = cos(z+pih-p) * sin(z+p)
1655      n1t2i21 = sin(z+pih-p) * cos(z+p)
1656      n1t2i22 = sin(z+pih-p) * sin(z+p)
1657      n2t1i11 = cos(z-pih+p) * cos(z-p)
1658      n2t1i12 = cos(z-pih+p) * sin(z-p)
1659      n2t1i21 = sin(z-pih+p) * cos(z-p)
1660      n2t1i22 = sin(z-pih+p) * sin(z-p)
1661      t1t2i11 = cos(z-p) * cos(z+p)
1662      t1t2i12 = cos(z-p) * sin(z+p)
1663      t1t2i21 = sin(z-p) * cos(z+p)
1664      t1t2i22 = sin(z-p) * sin(z+p)
1665      t2t1i11 = cos(z+p) * cos(z-p)
1666      t2t1i12 = cos(z+p) * sin(z-p)
1667      t2t1i21 = sin(z+p) * cos(z-p)
1668      t2t1i22 = sin(z+p) * sin(z-p)
1669      d11 = cos(y)*cos(y)*(cos(x)+sin(x)*tan(y)*tan(y))
1670      d12 = cos(y)*cos(y)*tan(y)*(-cos(x)+sin(x))
1671      d22 = cos(y)*cos(y)*(sin(x)+cos(x)*tan(y)*tan(y))
1672      IIn1t2 = n1t2i11 * d11 + (n1t2i12 + n1t2i21) * d12 + n1t2i22 * d22
1673      IIn2t1 = n2t1i11 * d11 + (n2t1i12 + n2t1i21) * d12 + n2t1i22 * d22
1674      IIt1t2 = t1t2i11 * d11 + (t1t2i12 + t1t2i21) * d12 + t1t2i22 * d22
1675
1676      IF (-IIn1t2>=rsmall) THEN
1677      Hen1t2 = 1._wp
1678      ELSE
1679      Hen1t2 = 0._wp
1680      ENDIF
1681
1682      IF (-IIn2t1>=rsmall) THEN
1683      Hen2t1 = 1._wp
1684      ELSE
1685      Hen2t1 = 0._wp
1686      ENDIF
1687
1688      s22ks = sign(1._wp,IIt1t2+rsmall)*(Hen1t2 * t1t2i22 + Hen2t1 * t2t1i22)
1689
1690   END FUNCTION s22ks
1691
1692#else
1693   !!----------------------------------------------------------------------
1694   !!   Default option         Empty module           NO SI3 sea-ice model
1695   !!----------------------------------------------------------------------
1696
1697#endif
1698   !!==============================================================================
1699END MODULE icedyn_rhg_eap
Note: See TracBrowser for help on using the repository browser.