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/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/tests/ICE_RHEO/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r13648_ASINTER-04_laurent_bulk_ice/tests/ICE_RHEO/MY_SRC/icedyn_rhg_eap.F90 @ 14021

Last change on this file since 14021 was 14021, checked in by laurent, 3 years ago

Caught up with trunk rev 14020...

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