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

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

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

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

fixed bug in stress calculation

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