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_evp.F90 in NEMO/trunk/src/ICE – NEMO

source: NEMO/trunk/src/ICE/icedyn_rhg_evp.F90 @ 15309

Last change on this file since 15309 was 15266, checked in by clem, 3 years ago

forgotten output in evp

  • Property svn:keywords set to Id
File size: 61.9 KB
RevLine 
[8586]1MODULE icedyn_rhg_evp
2   !!======================================================================
3   !!                     ***  MODULE  icedyn_rhg_evp  ***
4   !!   Sea-Ice dynamics : rheology Elasto-Viscous-Plastic
5   !!======================================================================
6   !! History :   -   !  2007-03  (M.A. Morales Maqueda, S. Bouillon) Original code
[9656]7   !!            3.0  !  2008-03  (M. Vancoppenolle) adaptation to new model
[14072]8   !!             -   !  2008-11  (M. Vancoppenolle, S. Bouillon, Y. Aksenov) add surface tilt in ice rheolohy
[9604]9   !!            3.3  !  2009-05  (G.Garric)    addition of the evp case
[14072]10   !!            3.4  !  2011-01  (A. Porter)   dynamical allocation
[9604]11   !!            3.5  !  2012-08  (R. Benshila) AGRIF
[9656]12   !!            3.6  !  2016-06  (C. Rousset)  Rewriting + landfast ice + mEVP (Bouillon 2013)
[9604]13   !!            3.7  !  2017     (C. Rousset)  add aEVP (Kimmritz 2016-2017)
14   !!            4.0  !  2018     (many people) SI3 [aka Sea Ice cube]
[8586]15   !!----------------------------------------------------------------------
[9570]16#if defined key_si3
[8586]17   !!----------------------------------------------------------------------
[9570]18   !!   'key_si3'                                       SI3 sea-ice model
[8586]19   !!----------------------------------------------------------------------
[8813]20   !!   ice_dyn_rhg_evp : computes ice velocities from EVP rheology
21   !!   rhg_evp_rst     : read/write EVP fields in ice restart
[8586]22   !!----------------------------------------------------------------------
23   USE phycst         ! Physical constant
24   USE dom_oce        ! Ocean domain
25   USE sbc_oce , ONLY : ln_ice_embd, nn_fsbc, ssh_m
26   USE sbc_ice , ONLY : utau_ice, vtau_ice, snwice_mass, snwice_mass_b
27   USE ice            ! sea-ice: ice variables
[10332]28   USE icevar         ! ice_var_sshdyn
[8586]29   USE icedyn_rdgrft  ! sea-ice: ice strength
[14072]30   USE bdy_oce , ONLY : ln_bdy
31   USE bdyice
[8813]32#if defined key_agrif
[9596]33   USE agrif_ice_interp
[8813]34#endif
[8586]35   !
36   USE in_out_manager ! I/O manager
37   USE iom            ! I/O manager library
38   USE lib_mpp        ! MPP library
39   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
40   USE lbclnk         ! lateral boundary conditions (or mpp links)
41   USE prtctl         ! Print control
42
[13472]43   USE netcdf         ! NetCDF library for convergence test
[8586]44   IMPLICIT NONE
45   PRIVATE
46
47   PUBLIC   ice_dyn_rhg_evp   ! called by icedyn_rhg.F90
48   PUBLIC   rhg_evp_rst       ! called by icedyn_rhg.F90
49
[13472]50   !! for convergence tests
51   INTEGER ::   ncvgid   ! netcdf file id
52   INTEGER ::   nvarid   ! netcdf variable id
[15014]53   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   fimask   ! mask at F points for the ice
54
55   !! * Substitutions
56#  include "do_loop_substitute.h90"
[15062]57#  include "domzgr_substitute.h90"
[8586]58   !!----------------------------------------------------------------------
[9598]59   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
[10069]60   !! $Id$
[10068]61   !! Software governed by the CeCILL license (see ./LICENSE)
[8586]62   !!----------------------------------------------------------------------
63CONTAINS
64
[12377]65   SUBROUTINE ice_dyn_rhg_evp( kt, Kmm, pstress1_i, pstress2_i, pstress12_i, pshear_i, pdivu_i, pdelta_i )
[8586]66      !!-------------------------------------------------------------------
67      !!                 ***  SUBROUTINE ice_dyn_rhg_evp  ***
[8813]68      !!                             EVP-C-grid
[8586]69      !!
70      !! ** purpose : determines sea ice drift from wind stress, ice-ocean
[14072]71      !!  stress and sea-surface slope. Ice-ice interaction is described by
72      !!  a non-linear elasto-viscous-plastic (EVP) law including shear
73      !!  strength and a bulk rheology (Hunke and Dukowicz, 2002).
[8586]74      !!
75      !!  The points in the C-grid look like this, dear reader
76      !!
77      !!                              (ji,jj)
78      !!                                 |
79      !!                                 |
80      !!                      (ji-1,jj)  |  (ji,jj)
[14072]81      !!                             ---------
[8586]82      !!                            |         |
83      !!                            | (ji,jj) |------(ji,jj)
84      !!                            |         |
[14072]85      !!                             ---------
[8586]86      !!                     (ji-1,jj-1)     (ji,jj-1)
87      !!
88      !! ** Inputs  : - wind forcing (stress), oceanic currents
89      !!                ice total volume (vt_i) per unit area
90      !!                snow total volume (vt_s) per unit area
91      !!
[14072]92      !! ** Action  : - compute u_ice, v_ice : the components of the
[8586]93      !!                sea-ice velocity vector
94      !!              - compute delta_i, shear_i, divu_i, which are inputs
95      !!                of the ice thickness distribution
96      !!
97      !! ** Steps   : 0) compute mask at F point
[14072]98      !!              1) Compute ice snow mass, ice strength
[8586]99      !!              2) Compute wind, oceanic stresses, mass terms and
100      !!                 coriolis terms of the momentum equation
101      !!              3) Solve the momentum equation (iterative procedure)
102      !!              4) Recompute delta, shear and divergence
103      !!                 (which are inputs of the ITD) & store stress
104      !!                 for the next time step
105      !!              5) Diagnostics including charge ellipse
106      !!
[8813]107      !! ** Notes   : There is the possibility to use aEVP from the nice work of Kimmritz et al. (2016 & 2017)
108      !!              by setting up ln_aEVP=T (i.e. changing alpha and beta parameters).
109      !!              This is an upgraded version of mEVP from Bouillon et al. 2013
110      !!              (i.e. more stable and better convergence)
[8586]111      !!
112      !! References : Hunke and Dukowicz, JPO97
113      !!              Bouillon et al., Ocean Modelling 2009
114      !!              Bouillon et al., Ocean Modelling 2013
[8813]115      !!              Kimmritz et al., Ocean Modelling 2016 & 2017
[8586]116      !!-------------------------------------------------------------------
[8813]117      INTEGER                 , INTENT(in   ) ::   kt                                    ! time step
[12377]118      INTEGER                 , INTENT(in   ) ::   Kmm                                   ! ocean time level index
[8813]119      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pstress1_i, pstress2_i, pstress12_i   !
120      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   pshear_i  , pdivu_i   , pdelta_i      !
[8586]121      !!
122      INTEGER ::   ji, jj       ! dummy loop indices
123      INTEGER ::   jter         ! local integers
[8813]124      !
[12489]125      REAL(wp) ::   zrhoco                                              ! rho0 * rn_cio
[9049]126      REAL(wp) ::   zdtevp, z1_dtevp                                    ! time step for subcycling
127      REAL(wp) ::   ecc2, z1_ecc2                                       ! square of yield ellipse eccenticity
128      REAL(wp) ::   zalph1, z1_alph1, zalph2, z1_alph2                  ! alpha coef from Bouillon 2009 or Kimmritz 2017
[13472]129      REAl(wp) ::   zbetau, zbetav
[10413]130      REAL(wp) ::   zm1, zm2, zm3, zmassU, zmassV, zvU, zvV             ! ice/snow mass and volume
[13550]131      REAL(wp) ::   zp_delf, zds2, zdt, zdt2, zdiv, zdiv2               ! temporary scalars
[11536]132      REAL(wp) ::   zTauO, zTauB, zRHS, zvel                            ! temporary scalars
[10413]133      REAL(wp) ::   zkt                                                 ! isotropic tensile strength for landfast ice
134      REAL(wp) ::   zvCr                                                ! critical ice volume above which ice is landfast
[8813]135      !
[9049]136      REAL(wp) ::   zintb, zintn                                        ! dummy argument
[8586]137      REAL(wp) ::   zfac_x, zfac_y
138      REAL(wp) ::   zshear, zdum1, zdum2
[8813]139      !
[13550]140      REAL(wp), DIMENSION(jpi,jpj) ::   zdelta, zp_delt                 ! delta and P/delta at T points
[8813]141      REAL(wp), DIMENSION(jpi,jpj) ::   zbeta                           ! beta coef from Kimmritz 2017
[8586]142      !
[8813]143      REAL(wp), DIMENSION(jpi,jpj) ::   zdt_m                           ! (dt / ice-snow_mass) on T points
[11536]144      REAL(wp), DIMENSION(jpi,jpj) ::   zaU  , zaV                      ! ice fraction on U/V points
[8813]145      REAL(wp), DIMENSION(jpi,jpj) ::   zmU_t, zmV_t                    ! (ice-snow_mass / dt) on U/V points
[8586]146      REAL(wp), DIMENSION(jpi,jpj) ::   zmf                             ! coriolis parameter at T points
[13550]147      REAL(wp), DIMENSION(jpi,jpj) ::   v_oceU, u_oceV, v_iceU, u_iceV  ! ocean/ice u/v component on V/U points
[8813]148      !
[8586]149      REAL(wp), DIMENSION(jpi,jpj) ::   zds                             ! shear
[13612]150      REAL(wp), DIMENSION(jpi,jpj) ::   zten_i                          ! tension
[8586]151      REAL(wp), DIMENSION(jpi,jpj) ::   zs1, zs2, zs12                  ! stress tensor components
[10415]152      REAL(wp), DIMENSION(jpi,jpj) ::   zsshdyn                         ! array used for the calculation of ice surface slope:
[9049]153      !                                                                 !    ocean surface (ssh_m) if ice is not embedded
[14072]154      !                                                                 !    ice bottom surface if ice is embedded
[11536]155      REAL(wp), DIMENSION(jpi,jpj) ::   zfU  , zfV                      ! internal stresses
156      REAL(wp), DIMENSION(jpi,jpj) ::   zspgU, zspgV                    ! surface pressure gradient at U/V points
157      REAL(wp), DIMENSION(jpi,jpj) ::   zCorU, zCorV                    ! Coriolis stress array
158      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_ai, ztauy_ai              ! ice-atm. stress at U-V points
159      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_oi, ztauy_oi              ! ice-ocean stress at U-V points
160      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_bi, ztauy_bi              ! ice-OceanBottom stress at U-V points (landfast)
161      REAL(wp), DIMENSION(jpi,jpj) ::   ztaux_base, ztauy_base          ! ice-bottom stress at U-V points (landfast)
[8813]162      !
[15014]163      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk00, zmsk15
[11536]164      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk01x, zmsk01y                ! dummy arrays
165      REAL(wp), DIMENSION(jpi,jpj) ::   zmsk00x, zmsk00y                ! mask for ice presence
[8586]166
167      REAL(wp), PARAMETER          ::   zepsi  = 1.0e-20_wp             ! tolerance parameter
[10891]168      REAL(wp), PARAMETER          ::   zmmin  = 1._wp                  ! ice mass (kg/m2)  below which ice velocity becomes very small
169      REAL(wp), PARAMETER          ::   zamin  = 0.001_wp               ! ice concentration below which ice velocity becomes very small
[13472]170      !! --- check convergence
171      REAL(wp), DIMENSION(jpi,jpj) ::   zu_ice, zv_ice
[8586]172      !! --- diags
[13601]173      REAL(wp) ::   zsig1, zsig2, zsig12, zfac, z1_strength
[14072]174      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zsig_I, zsig_II, zsig1_p, zsig2_p
[8586]175      !! --- SIMIP diags
176      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s)
177      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s)
178      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s)
179      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s)
180      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_xatrp     ! X-component of area transport (m2/s)
[14072]181      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zdiag_yatrp     ! Y-component of area transport (m2/s)
[8586]182      !!-------------------------------------------------------------------
183
184      IF( kt == nit000 .AND. lwp )   WRITE(numout,*) '-- ice_dyn_rhg_evp: EVP sea-ice rheology'
185      !
[13472]186      ! for diagnostics and convergence tests
[15049]187      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[13472]188         zmsk00(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06  ) ) ! 1 if ice    , 0 if no ice
189      END_2D
[15014]190      IF( nn_rhg_chkcvg > 0 ) THEN
[15049]191         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[15014]192            zmsk15(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - 0.15_wp ) ) ! 1 if 15% ice, 0 if less
193         END_2D
194      ENDIF
[13472]195      !
[8586]196      !------------------------------------------------------------------------------!
197      ! 0) mask at F points for the ice
198      !------------------------------------------------------------------------------!
[15014]199      IF( kt == nit000 ) THEN
200         ! ocean/land mask
201         ALLOCATE( fimask(jpi,jpj) )
202         IF( rn_ishlat == 0._wp ) THEN
203            DO_2D( 0, 0, 0, 0 )
204               fimask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
205            END_2D
206         ELSE
207            DO_2D( 0, 0, 0, 0 )
208               fimask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1)
209               ! Lateral boundary conditions on velocity (modify fimask)
210               IF( fimask(ji,jj) == 0._wp ) THEN
211                  fimask(ji,jj) = rn_ishlat * MIN( 1._wp , MAX( umask(ji,jj,1), umask(ji,jj+1,1), &
212                     &                                          vmask(ji,jj,1), vmask(ji+1,jj,1) ) )
213               ENDIF
214            END_2D
[12377]215         ENDIF
[15014]216         CALL lbc_lnk( 'icedyn_rhg_evp', fimask, 'F', 1._wp )
217      ENDIF
[8586]218      !------------------------------------------------------------------------------!
219      ! 1) define some variables and initialize arrays
220      !------------------------------------------------------------------------------!
[14072]221      zrhoco = rho0 * rn_cio
[8586]222
223      ! ecc2: square of yield ellipse eccenticrity
224      ecc2    = rn_ecc * rn_ecc
225      z1_ecc2 = 1._wp / ecc2
226
227      ! alpha parameters (Bouillon 2009)
[8813]228      IF( .NOT. ln_aEVP ) THEN
[13472]229         zdtevp   = rDt_ice / REAL( nn_nevp )
230         zalph1 =   2._wp * rn_relast * REAL( nn_nevp )
[8813]231         zalph2 = zalph1 * z1_ecc2
[8586]232
[8813]233         z1_alph1 = 1._wp / ( zalph1 + 1._wp )
234         z1_alph2 = 1._wp / ( zalph2 + 1._wp )
[13472]235      ELSE
236         zdtevp   = rdt_ice
237         ! zalpha parameters set later on adaptatively
[8813]238      ENDIF
[13472]239      z1_dtevp = 1._wp / zdtevp
[14072]240
241      ! Initialise stress tensor
242      zs1 (:,:) = pstress1_i (:,:)
[8586]243      zs2 (:,:) = pstress2_i (:,:)
244      zs12(:,:) = pstress12_i(:,:)
245
246      ! Ice strength
247      CALL ice_strength
248
[10413]249      ! landfast param from Lemieux(2016): add isotropic tensile strength (following Konig Beatty and Holland, 2010)
[13472]250      IF( ln_landfast_L16 ) THEN   ;   zkt = rn_lf_tensile
[11536]251      ELSE                         ;   zkt = 0._wp
[10413]252      ENDIF
[8586]253      !
254      !------------------------------------------------------------------------------!
255      ! 2) Wind / ocean stress, mass terms, coriolis terms
256      !------------------------------------------------------------------------------!
[10415]257      ! sea surface height
258      !    embedded sea ice: compute representative ice top surface
259      !    non-embedded sea ice: use ocean surface for slope calculation
260      zsshdyn(:,:) = ice_var_sshdyn( ssh_m, snwice_mass, snwice_mass_b)
[8586]261
[15049]262      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
263         zm1          = ( rhos * vt_s(ji,jj) + rhoi * vt_i(ji,jj) )  ! Ice/snow mass at U-V points
264         zmf  (ji,jj) = zm1 * ff_t(ji,jj)                            ! Coriolis at T points (m*f)
265         zdt_m(ji,jj) = zdtevp / MAX( zm1, zmmin )                   ! dt/m at T points (for alpha and beta coefficients)
266      END_2D
267     
268      DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[8586]269
[12377]270         ! ice fraction at U-V points
271         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)
272         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)
[8586]273
[12377]274         ! Ice/snow mass at U-V points
275         zm1 = ( rhos * vt_s(ji  ,jj  ) + rhoi * vt_i(ji  ,jj  ) )
276         zm2 = ( rhos * vt_s(ji+1,jj  ) + rhoi * vt_i(ji+1,jj  ) )
277         zm3 = ( rhos * vt_s(ji  ,jj+1) + rhoi * vt_i(ji  ,jj+1) )
278         zmassU = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm2 * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1)
279         zmassV = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm3 * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1)
[8586]280
[12377]281         ! Ocean currents at U-V points
[15263]282         ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility)
283         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)
284         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)
[8586]285
[12377]286         ! m/dt
287         zmU_t(ji,jj)    = zmassU * z1_dtevp
288         zmV_t(ji,jj)    = zmassV * z1_dtevp
[14072]289
[12377]290         ! Drag ice-atm.
291         ztaux_ai(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj)
292         ztauy_ai(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj)
[8586]293
[12377]294         ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points
295         zspgU(ji,jj)    = - zmassU * grav * ( zsshdyn(ji+1,jj) - zsshdyn(ji,jj) ) * r1_e1u(ji,jj)
296         zspgV(ji,jj)    = - zmassV * grav * ( zsshdyn(ji,jj+1) - zsshdyn(ji,jj) ) * r1_e2v(ji,jj)
[8586]297
[12377]298         ! masks
299         zmsk00x(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) )  ! 0 if no ice
300         zmsk00y(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) )  ! 0 if no ice
[8586]301
[12377]302         ! switches
303         IF( zmassU <= zmmin .AND. zaU(ji,jj) <= zamin ) THEN   ;   zmsk01x(ji,jj) = 0._wp
304         ELSE                                                   ;   zmsk01x(ji,jj) = 1._wp   ;   ENDIF
305         IF( zmassV <= zmmin .AND. zaV(ji,jj) <= zamin ) THEN   ;   zmsk01y(ji,jj) = 0._wp
306         ELSE                                                   ;   zmsk01y(ji,jj) = 1._wp   ;   ENDIF
[8586]307
[12377]308      END_2D
[8586]309      !
[10413]310      !                                  !== Landfast ice parameterization ==!
311      !
312      IF( ln_landfast_L16 ) THEN         !-- Lemieux 2016
[15049]313         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]314            ! ice thickness at U-V points
315            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)
316            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)
317            ! ice-bottom stress at U points
[14005]318            zvCr = zaU(ji,jj) * rn_lf_depfra * hu(ji,jj,Kmm) * ( 1._wp - icb_mask(ji,jj) ) ! if grounded icebergs are read: ocean depth = 0
[13472]319            ztaux_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvU - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaU(ji,jj) ) )
[12377]320            ! ice-bottom stress at V points
[14005]321            zvCr = zaV(ji,jj) * rn_lf_depfra * hv(ji,jj,Kmm) * ( 1._wp - icb_mask(ji,jj) ) ! if grounded icebergs are read: ocean depth = 0
[13472]322            ztauy_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvV - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaV(ji,jj) ) )
[12377]323            ! ice_bottom stress at T points
[14005]324            zvCr = at_i(ji,jj) * rn_lf_depfra * ht(ji,jj) * ( 1._wp - icb_mask(ji,jj) )    ! if grounded icebergs are read: ocean depth = 0
[13472]325            tau_icebfr(ji,jj) = - rn_lf_bfr * MAX( 0._wp, vt_i(ji,jj) - zvCr ) * EXP( -rn_crhg * ( 1._wp - at_i(ji,jj) ) )
[12377]326         END_2D
[13226]327         CALL lbc_lnk( 'icedyn_rhg_evp', tau_icebfr(:,:), 'T', 1.0_wp )
[10413]328         !
[11536]329      ELSE                               !-- no landfast
[15049]330         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]331            ztaux_base(ji,jj) = 0._wp
332            ztauy_base(ji,jj) = 0._wp
333         END_2D
[10413]334      ENDIF
335
[8586]336      !------------------------------------------------------------------------------!
337      ! 3) Solution of the momentum equation, iterative procedure
338      !------------------------------------------------------------------------------!
339      !
[11536]340      !                                               ! ==================== !
[8586]341      DO jter = 1 , nn_nevp                           !    loop over jter    !
[14072]342         !                                            ! ==================== !
[10425]343         l_full_nf_update = jter == nn_nevp   ! false: disable full North fold update (performances) for iter = 1 to nn_nevp-1
344         !
[13472]345         ! convergence test
346         IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2  ) THEN
347            DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
348               zu_ice(ji,jj) = u_ice(ji,jj) * umask(ji,jj,1) ! velocity at previous time step
349               zv_ice(ji,jj) = v_ice(ji,jj) * vmask(ji,jj,1)
350            END_2D
351         ENDIF
[8586]352
353         ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- !
[15049]354         DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 )
[8586]355
[12377]356            ! shear at F points
357            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)   &
358               &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
[15014]359               &         ) * r1_e1e2f(ji,jj) * fimask(ji,jj)
[8586]360
[12377]361         END_2D
[8586]362
[13550]363         DO_2D( 0, 0, 0, 0 )
[8586]364
[12377]365            ! shear**2 at T points (doc eq. A16)
366            zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
367               &   + 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)  &
368               &   ) * 0.25_wp * r1_e1e2t(ji,jj)
[14072]369
[12377]370            ! divergence at T points
371            zdiv  = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
372               &    + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
373               &    ) * r1_e1e2t(ji,jj)
374            zdiv2 = zdiv * zdiv
[14072]375
[12377]376            ! tension at T points
377            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)   &
378               &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
379               &   ) * r1_e1e2t(ji,jj)
380            zdt2 = zdt * zdt
[14072]381
[12377]382            ! delta at T points
[14072]383            zdelta(ji,jj) = SQRT( zdiv2 + ( zdt2 + zds2 ) * z1_ecc2 )
[8586]384
[13550]385         END_2D
386         CALL lbc_lnk( 'icedyn_rhg_evp', zdelta, 'T', 1.0_wp )
[8813]387
[13550]388         ! P/delta at T points
[15049]389         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[13550]390            zp_delt(ji,jj) = strength(ji,jj) / ( zdelta(ji,jj) + rn_creepl )
391         END_2D
392
[15049]393         DO_2D( nn_hls-1, nn_hls, nn_hls-1, nn_hls )   ! loop ends at jpi,jpj so that no lbc_lnk are needed for zs1 and zs2
[13550]394
395            ! divergence at T points (duplication to avoid communications)
[15263]396            ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility)
397            zdiv  = ( (e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj))   &
398               &    + (e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1))   &
[13550]399               &    ) * r1_e1e2t(ji,jj)
[14072]400
[13550]401            ! tension at T points (duplication to avoid communications)
402            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)   &
403               &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
404               &   ) * r1_e1e2t(ji,jj)
[14072]405
[13472]406            ! alpha for aEVP
[12377]407            !   gamma = 0.5*P/(delta+creepl) * (c*pi)**2/Area * dt/m
408            !   alpha = beta = sqrt(4*gamma)
409            IF( ln_aEVP ) THEN
410               zalph1   = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
411               z1_alph1 = 1._wp / ( zalph1 + 1._wp )
412               zalph2   = zalph1
413               z1_alph2 = z1_alph1
[13472]414               ! explicit:
415               ! z1_alph1 = 1._wp / zalph1
416               ! z1_alph2 = 1._wp / zalph1
417               ! zalph1 = zalph1 - 1._wp
418               ! zalph2 = zalph1
[12377]419            ENDIF
[14072]420
[12377]421            ! stress at T points (zkt/=0 if landfast)
[13550]422            zs1(ji,jj) = ( zs1(ji,jj)*zalph1 + zp_delt(ji,jj) * ( zdiv*(1._wp + zkt) - zdelta(ji,jj)*(1._wp - zkt) ) ) * z1_alph1
423            zs2(ji,jj) = ( zs2(ji,jj)*zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 * (1._wp + zkt) ) ) * z1_alph2
[14072]424
[12377]425         END_2D
[8586]426
[13472]427         ! Save beta at T-points for further computations
428         IF( ln_aEVP ) THEN
[15049]429            DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[13472]430               zbeta(ji,jj) = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) )
431            END_2D
432         ENDIF
[14072]433
[15049]434         DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 )
[8586]435
[13472]436            ! alpha for aEVP
[12377]437            IF( ln_aEVP ) THEN
[13472]438               zalph2   = MAX( zbeta(ji,jj), zbeta(ji+1,jj), zbeta(ji,jj+1), zbeta(ji+1,jj+1) )
[12377]439               z1_alph2 = 1._wp / ( zalph2 + 1._wp )
[13472]440               ! explicit:
441               ! z1_alph2 = 1._wp / zalph2
442               ! zalph2 = zalph2 - 1._wp
[12377]443            ENDIF
[14072]444
[12377]445            ! P/delta at F points
[15263]446            ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility)
447            zp_delf = 0.25_wp * ( (zp_delt(ji,jj) + zp_delt(ji+1,jj)) + (zp_delt(ji,jj+1) + zp_delt(ji+1,jj+1)) )
[14072]448
[12377]449            ! stress at F points (zkt/=0 if landfast)
450            zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 * (1._wp + zkt) ) * 0.5_wp ) * z1_alph2
[8586]451
[12377]452         END_2D
[8586]453
454         ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- !
[15263]455         ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility)
[15049]456         DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]457            !                   !--- U points
[15263]458            zfU(ji,jj) = 0.5_wp * ( (( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj)                                             &
[12377]459               &                  + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj)    &
[15263]460               &                    ) * r1_e2u(ji,jj))                                                                      &
[12377]461               &                  + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1)  &
462               &                    ) * 2._wp * r1_e1u(ji,jj)                                                              &
463               &                  ) * r1_e1e2u(ji,jj)
464            !
465            !                !--- V points
[15263]466            zfV(ji,jj) = 0.5_wp * ( (( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj)                                             &
[12377]467               &                  - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj)    &
[15263]468               &                    ) * r1_e1v(ji,jj))                                                                      &
[12377]469               &                  + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj)  &
470               &                    ) * 2._wp * r1_e2v(ji,jj)                                                              &
471               &                  ) * r1_e1e2v(ji,jj)
472            !
473            !                !--- ice currents at U-V point
[15263]474            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)
475            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)
[12377]476            !
477         END_2D
[8586]478         !
479         ! --- Computation of ice velocity --- !
[8813]480         !  Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta vary as in Kimmritz 2016 & 2017
[8586]481         !  Bouillon et al. 2009 (eq 34-35) => stable
482         IF( MOD(jter,2) == 0 ) THEN ! even iterations
483            !
[15049]484            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]485               !                 !--- tau_io/(v_oce - v_ice)
486               zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
487                  &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
488               !                 !--- Ocean-to-Ice stress
489               ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
490               !
491               !                 !--- tau_bottom/v_ice
492               zvel  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
493               zTauB = ztauy_base(ji,jj) / zvel
494               !                 !--- OceanBottom-to-Ice stress
495               ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
496               !
497               !                 !--- Coriolis at V-points (energy conserving formulation)
498               zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
499                  &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
500                  &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
501               !
502               !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
503               zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
504               !
505               !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
506               !                                         1 = sliding friction : TauB < RHS
507               rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
508               !
509               IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
[13472]510                  zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) )
511                  v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) )         & ! previous velocity
512                     &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
513                     &                                    ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
[14072]514                     &            + ( 1._wp - rswitch ) * (  v_ice_b(ji,jj)                                                   &
[13472]515                     &                                     + v_ice  (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
516                     &                                    ) / ( zbetav + 1._wp )                                              &
517                     &             ) * 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
[12377]518                     &           )   * zmsk00y(ji,jj)
519               ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
[13472]520                  v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * v_ice(ji,jj)                                       & ! previous velocity
521                     &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
522                     &                                    ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
523                     &            + ( 1._wp - rswitch ) *   v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
524                     &             ) * zmsk01y(ji,jj) + v_oce(ji,jj) * 0.01_wp * ( 1._wp - zmsk01y(ji,jj) )                   & ! v_ice = v_oce/100 if mass < zmmin & conc < zamin
525                     &            )  * zmsk00y(ji,jj)
[12377]526               ENDIF
527            END_2D
[15049]528            IF( nn_hls == 1 )   CALL lbc_lnk( 'icedyn_rhg_evp', v_ice, 'V', -1.0_wp )
[8586]529            !
[13295]530            DO_2D( 0, 0, 0, 0 )
[12377]531               !                 !--- tau_io/(u_oce - u_ice)
532               zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
533                  &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
534               !                 !--- Ocean-to-Ice stress
535               ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
536               !
537               !                 !--- tau_bottom/u_ice
538               zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
539               zTauB = ztaux_base(ji,jj) / zvel
540               !                 !--- OceanBottom-to-Ice stress
541               ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
542               !
543               !                 !--- Coriolis at U-points (energy conserving formulation)
544               zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
545                  &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
546                  &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
547               !
548               !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
549               zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
550               !
551               !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
552               !                                         1 = sliding friction : TauB < RHS
553               rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
554               !
555               IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
[13472]556                  zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) )
557                  u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) )         & ! previous velocity
558                     &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
559                     &                                    ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
560                     &            + ( 1._wp - rswitch ) * (  u_ice_b(ji,jj)                                                   &
561                     &                                     + u_ice  (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
562                     &                                    ) / ( zbetau + 1._wp )                                              &
[14072]563                     &             ) * 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
[12377]564                     &           )   * zmsk00x(ji,jj)
565               ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
[13472]566                  u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * u_ice(ji,jj)                                       & ! previous velocity
567                     &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
568                     &                                    ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
569                     &            + ( 1._wp - rswitch ) *   u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
[14072]570                     &             ) * 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
[13472]571                     &           )   * zmsk00x(ji,jj)
[12377]572               ENDIF
573            END_2D
[15049]574            IF( nn_hls == 1 ) THEN   ;   CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp )
575            ELSE                     ;   CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp, v_ice, 'V', -1.0_wp )
576            ENDIF
[8586]577            !
578         ELSE ! odd iterations
579            !
[15049]580            DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
[12377]581               !                 !--- tau_io/(u_oce - u_ice)
582               zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) )  &
583                  &                              + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) )
584               !                 !--- Ocean-to-Ice stress
585               ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) )
586               !
587               !                 !--- tau_bottom/u_ice
588               zvel  = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) )
589               zTauB = ztaux_base(ji,jj) / zvel
590               !                 !--- OceanBottom-to-Ice stress
591               ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj)
592               !
593               !                 !--- Coriolis at U-points (energy conserving formulation)
594               zCorU(ji,jj)  =   0.25_wp * r1_e1u(ji,jj) *  &
595                  &    ( zmf(ji  ,jj) * ( e1v(ji  ,jj) * v_ice(ji  ,jj) + e1v(ji  ,jj-1) * v_ice(ji  ,jj-1) )  &
596                  &    + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) )
597               !
598               !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
599               zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj)
600               !
601               !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
602               !                                         1 = sliding friction : TauB < RHS
603               rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
604               !
605               IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
[13472]606                  zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) )
607                  u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) )         & ! previous velocity
608                     &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
609                     &                                    ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
610                     &            + ( 1._wp - rswitch ) * (  u_ice_b(ji,jj)                                                   &
611                     &                                     + u_ice  (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
612                     &                                    ) / ( zbetau + 1._wp )                                              &
[14072]613                     &             ) * 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
[12377]614                     &           )   * zmsk00x(ji,jj)
615               ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
[13472]616                  u_ice(ji,jj) = ( (          rswitch   * ( zmU_t(ji,jj) * u_ice(ji,jj)                                       & ! previous velocity
617                     &                                    + zRHS + zTauO * u_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
618                     &                                    ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
619                     &            + ( 1._wp - rswitch ) *   u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
620                     &             ) * 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
621                     &           )   * zmsk00x(ji,jj)
[12377]622               ENDIF
623            END_2D
[15049]624            IF( nn_hls == 1 )   CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp )
[8586]625            !
[13295]626            DO_2D( 0, 0, 0, 0 )
[12377]627               !                 !--- tau_io/(v_oce - v_ice)
628               zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) )  &
629                  &                              + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) )
630               !                 !--- Ocean-to-Ice stress
631               ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) )
632               !
633               !                 !--- tau_bottom/v_ice
634               zvel  = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) )
635               zTauB = ztauy_base(ji,jj) / zvel
636               !                 !--- OceanBottom-to-Ice stress
637               ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj)
638               !
639               !                 !--- Coriolis at v-points (energy conserving formulation)
640               zCorV(ji,jj)  = - 0.25_wp * r1_e2v(ji,jj) *  &
641                  &    ( zmf(ji,jj  ) * ( e2u(ji,jj  ) * u_ice(ji,jj  ) + e2u(ji-1,jj  ) * u_ice(ji-1,jj  ) )  &
642                  &    + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) )
643               !
644               !                 !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io
645               zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj)
646               !
647               !                 !--- landfast switch => 0 = static  friction : TauB > RHS & sign(TauB) /= sign(RHS)
648               !                                         1 = sliding friction : TauB < RHS
649               rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) )
650               !
651               IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017)
[13472]652                  zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) )
653                  v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) )         & ! previous velocity
654                     &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
655                     &                                    ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast
656                     &            + ( 1._wp - rswitch ) * (  v_ice_b(ji,jj)                                                   &
657                     &                                     + v_ice  (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax )     & ! static friction => slow decrease to v=0
[14072]658                     &                                    ) / ( zbetav + 1._wp )                                              &
[13472]659                     &             ) * 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
[12377]660                     &           )   * zmsk00y(ji,jj)
661               ELSE               !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009)
[13472]662                  v_ice(ji,jj) = ( (          rswitch   * ( zmV_t(ji,jj) * v_ice(ji,jj)                                       & ! previous velocity
663                     &                                    + zRHS + zTauO * v_ice(ji,jj)                                       & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part)
664                     &                                    ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB )                      & ! m/dt + tau_io(only ice part) + landfast
665                     &            + ( 1._wp - rswitch ) *   v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax )         & ! static friction => slow decrease to v=0
666                     &             ) * 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
667                     &           )   * zmsk00y(ji,jj)
[12377]668               ENDIF
669            END_2D
[15049]670            IF( nn_hls == 1 ) THEN   ;   CALL lbc_lnk( 'icedyn_rhg_evp', v_ice, 'V', -1.0_wp )
671            ELSE                     ;   CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp, v_ice, 'V', -1.0_wp )
672            ENDIF
[8586]673            !
[15049]674         ENDIF
675         !
[8586]676#if defined key_agrif
[15049]677!!       CALL agrif_interp_ice( 'U', jter, nn_nevp )
678!!       CALL agrif_interp_ice( 'V', jter, nn_nevp )
679         CALL agrif_interp_ice( 'U' )
680         CALL agrif_interp_ice( 'V' )
[8586]681#endif
[15049]682         IF( ln_bdy )   CALL bdy_ice_dyn( 'U' )
683         IF( ln_bdy )   CALL bdy_ice_dyn( 'V' )
684         !
[13472]685         ! convergence test
[15014]686         IF( nn_rhg_chkcvg == 2 )   CALL rhg_cvg( kt, jter, nn_nevp, u_ice, v_ice, zu_ice, zv_ice, zmsk15 )
[8586]687         !
688         !                                                ! ==================== !
689      END DO                                              !  end loop over jter  !
690      !                                                   ! ==================== !
[13472]691      IF( ln_aEVP )   CALL iom_put( 'beta_evp' , zbeta )
[8586]692      !
693      !------------------------------------------------------------------------------!
[14072]694      ! 4) Recompute delta, shear and div (inputs for mechanical redistribution)
[8586]695      !------------------------------------------------------------------------------!
[15049]696      DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 )
[8586]697
[12377]698         ! shear at F points
699         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)   &
700            &         + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj)   &
[15014]701            &         ) * r1_e1e2f(ji,jj) * fimask(ji,jj)
[8586]702
[12377]703      END_2D
[14072]704
[13497]705      DO_2D( 0, 0, 0, 0 )   ! no vector loop
[14072]706
[12377]707         ! tension**2 at T points
708         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)   &
709            &   - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj)   &
710            &   ) * r1_e1e2t(ji,jj)
711         zdt2 = zdt * zdt
[13601]712
713         zten_i(ji,jj) = zdt
[14072]714
[12377]715         ! shear**2 at T points (doc eq. A16)
716         zds2 = ( zds(ji,jj  ) * zds(ji,jj  ) * e1e2f(ji,jj  ) + zds(ji-1,jj  ) * zds(ji-1,jj  ) * e1e2f(ji-1,jj  )  &
717            &   + 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)  &
718            &   ) * 0.25_wp * r1_e1e2t(ji,jj)
[14072]719
[12377]720         ! shear at T points
721         pshear_i(ji,jj) = SQRT( zdt2 + zds2 )
[8586]722
[12377]723         ! divergence at T points
724         pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)   &
725            &             + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)   &
726            &             ) * r1_e1e2t(ji,jj)
[14072]727
[12377]728         ! delta at T points
[14072]729         zfac            = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) ! delta
[13601]730         rswitch         = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zfac ) ) ! 0 if delta=0
731         pdelta_i(ji,jj) = zfac + rn_creepl * rswitch ! delta+creepl
[8586]732
[12377]733      END_2D
[14433]734      CALL lbc_lnk( 'icedyn_rhg_evp', pshear_i, 'T', 1._wp, pdivu_i, 'T', 1._wp, pdelta_i, 'T', 1._wp, zten_i, 'T', 1._wp, &
735         &                            zs1     , 'T', 1._wp, zs2    , 'T', 1._wp, zs12    , 'F', 1._wp )
[14072]736
[8586]737      ! --- Store the stress tensor for the next time step --- !
738      pstress1_i (:,:) = zs1 (:,:)
739      pstress2_i (:,:) = zs2 (:,:)
740      pstress12_i(:,:) = zs12(:,:)
741      !
742
743      !------------------------------------------------------------------------------!
744      ! 5) diagnostics
745      !------------------------------------------------------------------------------!
[11536]746      ! --- ice-ocean, ice-atm. & ice-oceanbottom(landfast) stresses --- !
747      IF(  iom_use('utau_oi') .OR. iom_use('vtau_oi') .OR. iom_use('utau_ai') .OR. iom_use('vtau_ai') .OR. &
748         & iom_use('utau_bi') .OR. iom_use('vtau_bi') ) THEN
749         !
[14433]750         CALL lbc_lnk( 'icedyn_rhg_evp', ztaux_oi, 'U', -1.0_wp, ztauy_oi, 'V', -1.0_wp, &
751            &                            ztaux_ai, 'U', -1.0_wp, ztauy_ai, 'V', -1.0_wp, &
752            &                            ztaux_bi, 'U', -1.0_wp, ztauy_bi, 'V', -1.0_wp )
[11536]753         !
754         CALL iom_put( 'utau_oi' , ztaux_oi * zmsk00 )
755         CALL iom_put( 'vtau_oi' , ztauy_oi * zmsk00 )
756         CALL iom_put( 'utau_ai' , ztaux_ai * zmsk00 )
757         CALL iom_put( 'vtau_ai' , ztauy_ai * zmsk00 )
758         CALL iom_put( 'utau_bi' , ztaux_bi * zmsk00 )
759         CALL iom_put( 'vtau_bi' , ztauy_bi * zmsk00 )
760      ENDIF
[14072]761
[8586]762      ! --- divergence, shear and strength --- !
[11536]763      IF( iom_use('icediv') )   CALL iom_put( 'icediv' , pdivu_i  * zmsk00 )   ! divergence
764      IF( iom_use('iceshe') )   CALL iom_put( 'iceshe' , pshear_i * zmsk00 )   ! shear
765      IF( iom_use('icestr') )   CALL iom_put( 'icestr' , strength * zmsk00 )   ! strength
[15266]766      IF( iom_use('icedlt') )   CALL iom_put( 'icedlt' , pdelta_i * zmsk00 )   ! delta
[8586]767
[13601]768      ! --- Stress tensor invariants (SIMIP diags) --- !
769      IF( iom_use('normstr') .OR. iom_use('sheastr') ) THEN
[8586]770         !
[13601]771         ALLOCATE( zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )
[14072]772         !
[15049]773         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[14072]774
[13601]775            ! Ice stresses
776            ! sigma1, sigma2, sigma12 are some useful recombination of the stresses (Hunke and Dukowicz MWR 2002, Bouillon et al., OM2013)
777            ! These are NOT stress tensor components, neither stress invariants, neither stress principal components
778            ! I know, this can be confusing...
[14072]779            zfac             =   strength(ji,jj) / ( pdelta_i(ji,jj) + rn_creepl )
[13601]780            zsig1            =   zfac * ( pdivu_i(ji,jj) - pdelta_i(ji,jj) )
781            zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
782            zsig12           =   zfac * z1_ecc2 * pshear_i(ji,jj)
[14072]783
[13601]784            ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008)
[13612]785            zsig_I (ji,jj)   =   zsig1 * 0.5_wp                                           ! 1st stress invariant, aka average normal stress, aka negative pressure
786            zsig_II(ji,jj)   =   SQRT ( MAX( 0._wp, zsig2 * zsig2 * 0.25_wp + zsig12 ) )  ! 2nd  ''       '', aka maximum shear stress
[14072]787
788         END_2D
[13601]789         !
790         ! Stress tensor invariants (normal and shear stress N/m) - SIMIP diags - definitions following Coon (1974) and Feltham (2008)
791         IF( iom_use('normstr') )   CALL iom_put( 'normstr', zsig_I (:,:) * zmsk00(:,:) ) ! Normal stress
792         IF( iom_use('sheastr') )   CALL iom_put( 'sheastr', zsig_II(:,:) * zmsk00(:,:) ) ! Maximum shear stress
[14072]793
[13601]794         DEALLOCATE ( zsig_I, zsig_II )
[14072]795
[13601]796      ENDIF
[8586]797
[13601]798      ! --- Normalized stress tensor principal components --- !
799      ! This are used to plot the normalized yield curve, see Lemieux & Dupont, 2020
800      ! Recommendation 1 : we use ice strength, not replacement pressure
801      ! Recommendation 2 : need to use deformations at PREVIOUS iterate for viscosities
802      IF( iom_use('sig1_pnorm') .OR. iom_use('sig2_pnorm') ) THEN
[8586]803         !
[14072]804         ALLOCATE( zsig1_p(jpi,jpj) , zsig2_p(jpi,jpj) , zsig_I(jpi,jpj) , zsig_II(jpi,jpj) )
805         !
[15049]806         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
[14072]807
808            ! Ice stresses computed with **viscosities** (delta, p/delta) at **previous** iterates
[13601]809            !                        and **deformations** at current iterates
810            !                        following Lemieux & Dupont (2020)
811            zfac             =   zp_delt(ji,jj)
812            zsig1            =   zfac * ( pdivu_i(ji,jj) - ( zdelta(ji,jj) + rn_creepl ) )
813            zsig2            =   zfac * z1_ecc2 * zten_i(ji,jj)
814            zsig12           =   zfac * z1_ecc2 * pshear_i(ji,jj)
[14072]815
[13601]816            ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008), T-point
[13612]817            zsig_I(ji,jj)    =   zsig1 * 0.5_wp                                            ! 1st stress invariant, aka average normal stress, aka negative pressure
818            zsig_II(ji,jj)   =   SQRT ( MAX( 0._wp, zsig2 * zsig2 * 0.25_wp + zsig12 ) )   ! 2nd  ''       '', aka maximum shear stress
[14072]819
[13601]820            ! Normalized  principal stresses (used to display the ellipse)
821            z1_strength      =   1._wp / MAX( 1._wp, strength(ji,jj) )
822            zsig1_p(ji,jj)   =   ( zsig_I(ji,jj) + zsig_II(ji,jj) ) * z1_strength
823            zsig2_p(ji,jj)   =   ( zsig_I(ji,jj) - zsig_II(ji,jj) ) * z1_strength
[14072]824         END_2D
[8586]825         !
[14072]826         CALL iom_put( 'sig1_pnorm' , zsig1_p )
827         CALL iom_put( 'sig2_pnorm' , zsig2_p )
[11536]828
[13601]829         DEALLOCATE( zsig1_p , zsig2_p , zsig_I, zsig_II )
[14072]830
[8586]831      ENDIF
[13472]832
[8586]833      ! --- SIMIP --- !
[11536]834      IF(  iom_use('dssh_dx') .OR. iom_use('dssh_dy') .OR. &
835         & iom_use('corstrx') .OR. iom_use('corstry') .OR. iom_use('intstrx') .OR. iom_use('intstry') ) THEN
836         !
[14433]837         CALL lbc_lnk( 'icedyn_rhg_evp', zspgU, 'U', -1.0_wp, zspgV, 'V', -1.0_wp, &
838            &                            zCorU, 'U', -1.0_wp, zCorV, 'V', -1.0_wp, zfU, 'U', -1.0_wp, zfV, 'V', -1.0_wp )
[8586]839
[11536]840         CALL iom_put( 'dssh_dx' , zspgU * zmsk00 )   ! Sea-surface tilt term in force balance (x)
841         CALL iom_put( 'dssh_dy' , zspgV * zmsk00 )   ! Sea-surface tilt term in force balance (y)
842         CALL iom_put( 'corstrx' , zCorU * zmsk00 )   ! Coriolis force term in force balance (x)
843         CALL iom_put( 'corstry' , zCorV * zmsk00 )   ! Coriolis force term in force balance (y)
844         CALL iom_put( 'intstrx' , zfU   * zmsk00 )   ! Internal force term in force balance (x)
845         CALL iom_put( 'intstry' , zfV   * zmsk00 )   ! Internal force term in force balance (y)
846      ENDIF
847
848      IF(  iom_use('xmtrpice') .OR. iom_use('ymtrpice') .OR. &
849         & iom_use('xmtrpsnw') .OR. iom_use('ymtrpsnw') .OR. iom_use('xatrp') .OR. iom_use('yatrp') ) THEN
850         !
851         ALLOCATE( zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) , &
852            &      zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp(jpi,jpj) , zdiag_yatrp(jpi,jpj) )
853         !
[13295]854         DO_2D( 0, 0, 0, 0 )
[12377]855            ! 2D ice mass, snow mass, area transport arrays (X, Y)
856            zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zmsk00(ji,jj)
857            zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zmsk00(ji,jj)
[11536]858
[12377]859            zdiag_xmtrp_ice(ji,jj) = rhoi * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component
860            zdiag_ymtrp_ice(ji,jj) = rhoi * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) !        ''           Y-   ''
[11536]861
[12377]862            zdiag_xmtrp_snw(ji,jj) = rhos * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component
863            zdiag_ymtrp_snw(ji,jj) = rhos * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) !          ''          Y-   ''
[11536]864
[12377]865            zdiag_xatrp(ji,jj)     = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) )        ! area transport,      X-component
866            zdiag_yatrp(ji,jj)     = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) )        !        ''            Y-   ''
[11536]867
[12377]868         END_2D
[8586]869
[14433]870         CALL lbc_lnk( 'icedyn_rhg_evp', zdiag_xmtrp_ice, 'U', -1.0_wp, zdiag_ymtrp_ice, 'V', -1.0_wp, &
871            &                            zdiag_xmtrp_snw, 'U', -1.0_wp, zdiag_ymtrp_snw, 'V', -1.0_wp, &
872            &                            zdiag_xatrp    , 'U', -1.0_wp, zdiag_yatrp    , 'V', -1.0_wp )
[8586]873
[11536]874         CALL iom_put( 'xmtrpice' , zdiag_xmtrp_ice )   ! X-component of sea-ice mass transport (kg/s)
[14072]875         CALL iom_put( 'ymtrpice' , zdiag_ymtrp_ice )   ! Y-component of sea-ice mass transport
[11536]876         CALL iom_put( 'xmtrpsnw' , zdiag_xmtrp_snw )   ! X-component of snow mass transport (kg/s)
877         CALL iom_put( 'ymtrpsnw' , zdiag_ymtrp_snw )   ! Y-component of snow mass transport
878         CALL iom_put( 'xatrp'    , zdiag_xatrp     )   ! X-component of ice area transport
879         CALL iom_put( 'yatrp'    , zdiag_yatrp     )   ! Y-component of ice area transport
880
881         DEALLOCATE( zdiag_xmtrp_ice , zdiag_ymtrp_ice , &
882            &        zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp , zdiag_yatrp )
883
[8586]884      ENDIF
885      !
[13472]886      ! --- convergence tests --- !
887      IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2 ) THEN
888         IF( iom_use('uice_cvg') ) THEN
889            IF( ln_aEVP ) THEN   ! output: beta * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
890               CALL iom_put( 'uice_cvg', MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * zbeta(:,:) * umask(:,:,1) , &
891                  &                           ABS( v_ice(:,:) - zv_ice(:,:) ) * zbeta(:,:) * vmask(:,:,1) ) * zmsk15(:,:) )
892            ELSE                 ! output: nn_nevp * ( u(t=nn_nevp) - u(t=nn_nevp-1) )
893               CALL iom_put( 'uice_cvg', REAL( nn_nevp ) * MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * umask(:,:,1) , &
894                  &                                             ABS( v_ice(:,:) - zv_ice(:,:) ) * vmask(:,:,1) ) * zmsk15(:,:) )
895            ENDIF
896         ENDIF
[14072]897      ENDIF
[13472]898      !
[8586]899   END SUBROUTINE ice_dyn_rhg_evp
900
[8813]901
[15014]902   SUBROUTINE rhg_cvg( kt, kiter, kitermax, pu, pv, pub, pvb, pmsk15 )
[13472]903      !!----------------------------------------------------------------------
904      !!                    ***  ROUTINE rhg_cvg  ***
[14072]905      !!
[13472]906      !! ** Purpose :   check convergence of oce rheology
907      !!
908      !! ** Method  :   create a file ice_cvg.nc containing the convergence of ice velocity
909      !!                during the sub timestepping of rheology so as:
910      !!                  uice_cvg = MAX( u(t+1) - u(t) , v(t+1) - v(t) )
911      !!                This routine is called every sub-iteration, so it is cpu expensive
912      !!
[14072]913      !! ** Note    :   for the first sub-iteration, uice_cvg is set to 0 (too large otherwise)
[13472]914      !!----------------------------------------------------------------------
915      INTEGER ,                 INTENT(in) ::   kt, kiter, kitermax       ! ocean time-step index
916      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pu, pv, pub, pvb          ! now and before velocities
[15014]917      REAL(wp), DIMENSION(:,:), INTENT(in) ::   pmsk15
[13472]918      !!
919      INTEGER           ::   it, idtime, istatus
920      INTEGER           ::   ji, jj          ! dummy loop indices
[14072]921      REAL(wp)          ::   zresm           ! local real
[13472]922      CHARACTER(len=20) ::   clname
923      !!----------------------------------------------------------------------
924
925      ! create file
926      IF( kt == nit000 .AND. kiter == 1 ) THEN
927         !
928         IF( lwp ) THEN
929            WRITE(numout,*)
930            WRITE(numout,*) 'rhg_cvg : ice rheology convergence control'
931            WRITE(numout,*) '~~~~~~~'
932         ENDIF
933         !
934         IF( lwm ) THEN
935            clname = 'ice_cvg.nc'
936            IF( .NOT. Agrif_Root() )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
937            istatus = NF90_CREATE( TRIM(clname), NF90_CLOBBER, ncvgid )
938            istatus = NF90_DEF_DIM( ncvgid, 'time'  , NF90_UNLIMITED, idtime )
939            istatus = NF90_DEF_VAR( ncvgid, 'uice_cvg', NF90_DOUBLE   , (/ idtime /), nvarid )
940            istatus = NF90_ENDDEF(ncvgid)
941         ENDIF
942         !
943      ENDIF
944
945      ! time
946      it = ( kt - 1 ) * kitermax + kiter
[14072]947
[13472]948      ! convergence
949      IF( kiter == 1 ) THEN ! remove the first iteration for calculations of convergence (always very large)
950         zresm = 0._wp
951      ELSE
[15014]952         zresm = 0._wp
953         DO_2D( 0, 0, 0, 0 )
954            zresm = MAX( zresm, MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), &
955               &                     ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * pmsk15(ji,jj) )
[13472]956         END_2D
957         CALL mpp_max( 'icedyn_rhg_evp', zresm )   ! max over the global domain
958      ENDIF
959
960      IF( lwm ) THEN
961         ! write variables
962         istatus = NF90_PUT_VAR( ncvgid, nvarid, (/zresm/), (/it/), (/1/) )
963         ! close file
[13601]964         IF( kt == nitend - nn_fsbc + 1 )   istatus = NF90_CLOSE(ncvgid)
[13472]965      ENDIF
[14072]966
[13472]967   END SUBROUTINE rhg_cvg
968
969
[8586]970   SUBROUTINE rhg_evp_rst( cdrw, kt )
971      !!---------------------------------------------------------------------
972      !!                   ***  ROUTINE rhg_evp_rst  ***
[14072]973      !!
[8586]974      !! ** Purpose :   Read or write RHG file in restart file
975      !!
976      !! ** Method  :   use of IOM library
977      !!----------------------------------------------------------------------
978      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
979      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
980      !
981      INTEGER  ::   iter            ! local integer
982      INTEGER  ::   id1, id2, id3   ! local integers
983      !!----------------------------------------------------------------------
984      !
985      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
986         !                                   ! ---------------
987         IF( ln_rstart ) THEN                   !* Read the restart file
988            !
989            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
990            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
991            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
992            !
993            IF( MIN( id1, id2, id3 ) > 0 ) THEN      ! fields exist
[13286]994               CALL iom_get( numrir, jpdom_auto, 'stress1_i' , stress1_i , cd_type = 'T' )
995               CALL iom_get( numrir, jpdom_auto, 'stress2_i' , stress2_i , cd_type = 'T' )
996               CALL iom_get( numrir, jpdom_auto, 'stress12_i', stress12_i, cd_type = 'F' )
[8586]997            ELSE                                     ! start rheology from rest
[9169]998               IF(lwp) WRITE(numout,*)
999               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without rheology, set stresses to 0'
[8586]1000               stress1_i (:,:) = 0._wp
1001               stress2_i (:,:) = 0._wp
1002               stress12_i(:,:) = 0._wp
1003            ENDIF
1004         ELSE                                   !* Start from rest
[9169]1005            IF(lwp) WRITE(numout,*)
1006            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set stresses to 0'
[8586]1007            stress1_i (:,:) = 0._wp
1008            stress2_i (:,:) = 0._wp
1009            stress12_i(:,:) = 0._wp
1010         ENDIF
1011         !
1012      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
1013         !                                   ! -------------------
1014         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
1015         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
1016         !
[13970]1017         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i )
1018         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i )
[8586]1019         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
1020         !
1021      ENDIF
1022      !
1023   END SUBROUTINE rhg_evp_rst
1024
[14072]1025
[8586]1026#else
1027   !!----------------------------------------------------------------------
[9570]1028   !!   Default option         Empty module           NO SI3 sea-ice model
[8586]1029   !!----------------------------------------------------------------------
1030#endif
1031
1032   !!==============================================================================
1033END MODULE icedyn_rhg_evp
Note: See TracBrowser for help on using the repository browser.