[8586] | 1 | MODULE 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 | |
---|
[15548] | 50 | !! for convergence tests |
---|
| 51 | INTEGER :: ncvgid ! netcdf file id |
---|
| 52 | INTEGER :: nvarid ! netcdf variable id |
---|
| 53 | REAL(wp), DIMENSION(:,:), ALLOCATABLE :: fimask ! mask at F points for the ice |
---|
| 54 | |
---|
[8586] | 55 | !! * Substitutions |
---|
[12377] | 56 | # include "do_loop_substitute.h90" |
---|
[13237] | 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 | !!---------------------------------------------------------------------- |
---|
| 63 | CONTAINS |
---|
| 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 | ! |
---|
[8586] | 136 | REAL(wp) :: zfac_x, zfac_y |
---|
[8813] | 137 | ! |
---|
[13550] | 138 | REAL(wp), DIMENSION(jpi,jpj) :: zdelta, zp_delt ! delta and P/delta at T points |
---|
[8813] | 139 | REAL(wp), DIMENSION(jpi,jpj) :: zbeta ! beta coef from Kimmritz 2017 |
---|
[8586] | 140 | ! |
---|
[8813] | 141 | REAL(wp), DIMENSION(jpi,jpj) :: zdt_m ! (dt / ice-snow_mass) on T points |
---|
[11536] | 142 | REAL(wp), DIMENSION(jpi,jpj) :: zaU , zaV ! ice fraction on U/V points |
---|
[8813] | 143 | REAL(wp), DIMENSION(jpi,jpj) :: zmU_t, zmV_t ! (ice-snow_mass / dt) on U/V points |
---|
[8586] | 144 | REAL(wp), DIMENSION(jpi,jpj) :: zmf ! coriolis parameter at T points |
---|
[13550] | 145 | REAL(wp), DIMENSION(jpi,jpj) :: v_oceU, u_oceV, v_iceU, u_iceV ! ocean/ice u/v component on V/U points |
---|
[8813] | 146 | ! |
---|
[8586] | 147 | REAL(wp), DIMENSION(jpi,jpj) :: zds ! shear |
---|
[13612] | 148 | REAL(wp), DIMENSION(jpi,jpj) :: zten_i ! tension |
---|
[8586] | 149 | REAL(wp), DIMENSION(jpi,jpj) :: zs1, zs2, zs12 ! stress tensor components |
---|
[10415] | 150 | REAL(wp), DIMENSION(jpi,jpj) :: zsshdyn ! array used for the calculation of ice surface slope: |
---|
[9049] | 151 | ! ! ocean surface (ssh_m) if ice is not embedded |
---|
[14072] | 152 | ! ! ice bottom surface if ice is embedded |
---|
[11536] | 153 | REAL(wp), DIMENSION(jpi,jpj) :: zfU , zfV ! internal stresses |
---|
| 154 | REAL(wp), DIMENSION(jpi,jpj) :: zspgU, zspgV ! surface pressure gradient at U/V points |
---|
| 155 | REAL(wp), DIMENSION(jpi,jpj) :: zCorU, zCorV ! Coriolis stress array |
---|
| 156 | REAL(wp), DIMENSION(jpi,jpj) :: ztaux_ai, ztauy_ai ! ice-atm. stress at U-V points |
---|
| 157 | REAL(wp), DIMENSION(jpi,jpj) :: ztaux_oi, ztauy_oi ! ice-ocean stress at U-V points |
---|
| 158 | REAL(wp), DIMENSION(jpi,jpj) :: ztaux_bi, ztauy_bi ! ice-OceanBottom stress at U-V points (landfast) |
---|
| 159 | REAL(wp), DIMENSION(jpi,jpj) :: ztaux_base, ztauy_base ! ice-bottom stress at U-V points (landfast) |
---|
[8813] | 160 | ! |
---|
[15548] | 161 | REAL(wp), DIMENSION(jpi,jpj) :: zmsk00, zmsk15 |
---|
[11536] | 162 | REAL(wp), DIMENSION(jpi,jpj) :: zmsk01x, zmsk01y ! dummy arrays |
---|
| 163 | REAL(wp), DIMENSION(jpi,jpj) :: zmsk00x, zmsk00y ! mask for ice presence |
---|
[8586] | 164 | |
---|
| 165 | REAL(wp), PARAMETER :: zepsi = 1.0e-20_wp ! tolerance parameter |
---|
[10891] | 166 | REAL(wp), PARAMETER :: zmmin = 1._wp ! ice mass (kg/m2) below which ice velocity becomes very small |
---|
| 167 | REAL(wp), PARAMETER :: zamin = 0.001_wp ! ice concentration below which ice velocity becomes very small |
---|
[13472] | 168 | !! --- check convergence |
---|
| 169 | REAL(wp), DIMENSION(jpi,jpj) :: zu_ice, zv_ice |
---|
[8586] | 170 | !! --- diags |
---|
[13601] | 171 | REAL(wp) :: zsig1, zsig2, zsig12, zfac, z1_strength |
---|
[14072] | 172 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zsig_I, zsig_II, zsig1_p, zsig2_p |
---|
[8586] | 173 | !! --- SIMIP diags |
---|
| 174 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_xmtrp_ice ! X-component of ice mass transport (kg/s) |
---|
| 175 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_ymtrp_ice ! Y-component of ice mass transport (kg/s) |
---|
| 176 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_xmtrp_snw ! X-component of snow mass transport (kg/s) |
---|
| 177 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_ymtrp_snw ! Y-component of snow mass transport (kg/s) |
---|
| 178 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_xatrp ! X-component of area transport (m2/s) |
---|
[14072] | 179 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zdiag_yatrp ! Y-component of area transport (m2/s) |
---|
[15548] | 180 | !! -- advect fields at the rheology time step for the calculation of strength |
---|
[15560] | 181 | !! it seems that convergence is worse when ll_advups=true. So it is not really a good idea |
---|
[15548] | 182 | LOGICAL :: ll_advups = .FALSE. |
---|
| 183 | REAL(wp) :: zdt_ups |
---|
| 184 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: za_i_ups, zv_i_ups ! tracers advected upstream |
---|
[8586] | 185 | !!------------------------------------------------------------------- |
---|
| 186 | |
---|
| 187 | IF( kt == nit000 .AND. lwp ) WRITE(numout,*) '-- ice_dyn_rhg_evp: EVP sea-ice rheology' |
---|
| 188 | ! |
---|
[13472] | 189 | ! for diagnostics and convergence tests |
---|
[15548] | 190 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
[13472] | 191 | zmsk00(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) ! 1 if ice , 0 if no ice |
---|
| 192 | END_2D |
---|
[15548] | 193 | IF( nn_rhg_chkcvg > 0 ) THEN |
---|
| 194 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
| 195 | zmsk15(ji,jj) = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - 0.15_wp ) ) ! 1 if 15% ice, 0 if less |
---|
| 196 | END_2D |
---|
| 197 | ENDIF |
---|
[13472] | 198 | ! |
---|
[8586] | 199 | !------------------------------------------------------------------------------! |
---|
| 200 | ! 0) mask at F points for the ice |
---|
| 201 | !------------------------------------------------------------------------------! |
---|
[15548] | 202 | IF( kt == nit000 ) THEN |
---|
| 203 | ! ocean/land mask |
---|
| 204 | ALLOCATE( fimask(jpi,jpj) ) |
---|
| 205 | IF( rn_ishlat == 0._wp ) THEN |
---|
| 206 | DO_2D( 0, 0, 0, 0 ) |
---|
| 207 | fimask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1) |
---|
| 208 | END_2D |
---|
| 209 | ELSE |
---|
| 210 | DO_2D( 0, 0, 0, 0 ) |
---|
| 211 | fimask(ji,jj) = tmask(ji,jj,1) * tmask(ji+1,jj,1) * tmask(ji,jj+1,1) * tmask(ji+1,jj+1,1) |
---|
| 212 | ! Lateral boundary conditions on velocity (modify fimask) |
---|
| 213 | IF( fimask(ji,jj) == 0._wp ) THEN |
---|
| 214 | fimask(ji,jj) = rn_ishlat * MIN( 1._wp , MAX( umask(ji,jj,1), umask(ji,jj+1,1), & |
---|
| 215 | & vmask(ji,jj,1), vmask(ji+1,jj,1) ) ) |
---|
| 216 | ENDIF |
---|
| 217 | END_2D |
---|
[12377] | 218 | ENDIF |
---|
[15548] | 219 | CALL lbc_lnk( 'icedyn_rhg_evp', fimask, 'F', 1._wp ) |
---|
| 220 | ENDIF |
---|
[8586] | 221 | !------------------------------------------------------------------------------! |
---|
| 222 | ! 1) define some variables and initialize arrays |
---|
| 223 | !------------------------------------------------------------------------------! |
---|
[14072] | 224 | zrhoco = rho0 * rn_cio |
---|
[8586] | 225 | |
---|
| 226 | ! ecc2: square of yield ellipse eccenticrity |
---|
| 227 | ecc2 = rn_ecc * rn_ecc |
---|
| 228 | z1_ecc2 = 1._wp / ecc2 |
---|
| 229 | |
---|
| 230 | ! alpha parameters (Bouillon 2009) |
---|
[8813] | 231 | IF( .NOT. ln_aEVP ) THEN |
---|
[13472] | 232 | zdtevp = rDt_ice / REAL( nn_nevp ) |
---|
| 233 | zalph1 = 2._wp * rn_relast * REAL( nn_nevp ) |
---|
[8813] | 234 | zalph2 = zalph1 * z1_ecc2 |
---|
[8586] | 235 | |
---|
[8813] | 236 | z1_alph1 = 1._wp / ( zalph1 + 1._wp ) |
---|
| 237 | z1_alph2 = 1._wp / ( zalph2 + 1._wp ) |
---|
[13472] | 238 | ELSE |
---|
[15548] | 239 | zdtevp = rDt_ice |
---|
[13472] | 240 | ! zalpha parameters set later on adaptatively |
---|
[8813] | 241 | ENDIF |
---|
[13472] | 242 | z1_dtevp = 1._wp / zdtevp |
---|
[14072] | 243 | |
---|
| 244 | ! Initialise stress tensor |
---|
| 245 | zs1 (:,:) = pstress1_i (:,:) |
---|
[8586] | 246 | zs2 (:,:) = pstress2_i (:,:) |
---|
| 247 | zs12(:,:) = pstress12_i(:,:) |
---|
| 248 | |
---|
| 249 | ! Ice strength |
---|
| 250 | CALL ice_strength |
---|
| 251 | |
---|
[10413] | 252 | ! landfast param from Lemieux(2016): add isotropic tensile strength (following Konig Beatty and Holland, 2010) |
---|
[13472] | 253 | IF( ln_landfast_L16 ) THEN ; zkt = rn_lf_tensile |
---|
[11536] | 254 | ELSE ; zkt = 0._wp |
---|
[10413] | 255 | ENDIF |
---|
[8586] | 256 | ! |
---|
| 257 | !------------------------------------------------------------------------------! |
---|
| 258 | ! 2) Wind / ocean stress, mass terms, coriolis terms |
---|
| 259 | !------------------------------------------------------------------------------! |
---|
[10415] | 260 | ! sea surface height |
---|
| 261 | ! embedded sea ice: compute representative ice top surface |
---|
| 262 | ! non-embedded sea ice: use ocean surface for slope calculation |
---|
| 263 | zsshdyn(:,:) = ice_var_sshdyn( ssh_m, snwice_mass, snwice_mass_b) |
---|
[8586] | 264 | |
---|
[15548] | 265 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
| 266 | zm1 = ( rhos * vt_s(ji,jj) + rhoi * vt_i(ji,jj) ) ! Ice/snow mass at U-V points |
---|
| 267 | zmf (ji,jj) = zm1 * ff_t(ji,jj) ! Coriolis at T points (m*f) |
---|
| 268 | zdt_m(ji,jj) = zdtevp / MAX( zm1, zmmin ) ! dt/m at T points (for alpha and beta coefficients) |
---|
| 269 | END_2D |
---|
| 270 | |
---|
| 271 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[8586] | 272 | |
---|
[12377] | 273 | ! ice fraction at U-V points |
---|
| 274 | zaU(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e1e2t(ji,jj) + at_i(ji+1,jj) * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1) |
---|
| 275 | zaV(ji,jj) = 0.5_wp * ( at_i(ji,jj) * e1e2t(ji,jj) + at_i(ji,jj+1) * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1) |
---|
[8586] | 276 | |
---|
[12377] | 277 | ! Ice/snow mass at U-V points |
---|
| 278 | zm1 = ( rhos * vt_s(ji ,jj ) + rhoi * vt_i(ji ,jj ) ) |
---|
| 279 | zm2 = ( rhos * vt_s(ji+1,jj ) + rhoi * vt_i(ji+1,jj ) ) |
---|
| 280 | zm3 = ( rhos * vt_s(ji ,jj+1) + rhoi * vt_i(ji ,jj+1) ) |
---|
| 281 | zmassU = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm2 * e1e2t(ji+1,jj) ) * r1_e1e2u(ji,jj) * umask(ji,jj,1) |
---|
| 282 | zmassV = 0.5_wp * ( zm1 * e1e2t(ji,jj) + zm3 * e1e2t(ji,jj+1) ) * r1_e1e2v(ji,jj) * vmask(ji,jj,1) |
---|
[8586] | 283 | |
---|
[12377] | 284 | ! Ocean currents at U-V points |
---|
[15548] | 285 | ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility) |
---|
| 286 | 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) |
---|
| 287 | 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] | 288 | |
---|
[12377] | 289 | ! m/dt |
---|
| 290 | zmU_t(ji,jj) = zmassU * z1_dtevp |
---|
| 291 | zmV_t(ji,jj) = zmassV * z1_dtevp |
---|
[14072] | 292 | |
---|
[12377] | 293 | ! Drag ice-atm. |
---|
| 294 | ztaux_ai(ji,jj) = zaU(ji,jj) * utau_ice(ji,jj) |
---|
| 295 | ztauy_ai(ji,jj) = zaV(ji,jj) * vtau_ice(ji,jj) |
---|
[8586] | 296 | |
---|
[12377] | 297 | ! Surface pressure gradient (- m*g*GRAD(ssh)) at U-V points |
---|
| 298 | zspgU(ji,jj) = - zmassU * grav * ( zsshdyn(ji+1,jj) - zsshdyn(ji,jj) ) * r1_e1u(ji,jj) |
---|
| 299 | zspgV(ji,jj) = - zmassV * grav * ( zsshdyn(ji,jj+1) - zsshdyn(ji,jj) ) * r1_e2v(ji,jj) |
---|
[8586] | 300 | |
---|
[12377] | 301 | ! masks |
---|
| 302 | zmsk00x(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassU ) ) ! 0 if no ice |
---|
| 303 | zmsk00y(ji,jj) = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zmassV ) ) ! 0 if no ice |
---|
[8586] | 304 | |
---|
[12377] | 305 | ! switches |
---|
| 306 | IF( zmassU <= zmmin .AND. zaU(ji,jj) <= zamin ) THEN ; zmsk01x(ji,jj) = 0._wp |
---|
| 307 | ELSE ; zmsk01x(ji,jj) = 1._wp ; ENDIF |
---|
| 308 | IF( zmassV <= zmmin .AND. zaV(ji,jj) <= zamin ) THEN ; zmsk01y(ji,jj) = 0._wp |
---|
| 309 | ELSE ; zmsk01y(ji,jj) = 1._wp ; ENDIF |
---|
[8586] | 310 | |
---|
[12377] | 311 | END_2D |
---|
[8586] | 312 | ! |
---|
[10413] | 313 | ! !== Landfast ice parameterization ==! |
---|
| 314 | ! |
---|
| 315 | IF( ln_landfast_L16 ) THEN !-- Lemieux 2016 |
---|
[15548] | 316 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[12377] | 317 | ! ice thickness at U-V points |
---|
| 318 | 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) |
---|
| 319 | 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) |
---|
| 320 | ! ice-bottom stress at U points |
---|
[14005] | 321 | 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] | 322 | ztaux_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvU - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaU(ji,jj) ) ) |
---|
[12377] | 323 | ! ice-bottom stress at V points |
---|
[14005] | 324 | 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] | 325 | ztauy_base(ji,jj) = - rn_lf_bfr * MAX( 0._wp, zvV - zvCr ) * EXP( -rn_crhg * ( 1._wp - zaV(ji,jj) ) ) |
---|
[12377] | 326 | ! ice_bottom stress at T points |
---|
[14005] | 327 | 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] | 328 | 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] | 329 | END_2D |
---|
[13226] | 330 | CALL lbc_lnk( 'icedyn_rhg_evp', tau_icebfr(:,:), 'T', 1.0_wp ) |
---|
[10413] | 331 | ! |
---|
[11536] | 332 | ELSE !-- no landfast |
---|
[15548] | 333 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[12377] | 334 | ztaux_base(ji,jj) = 0._wp |
---|
| 335 | ztauy_base(ji,jj) = 0._wp |
---|
| 336 | END_2D |
---|
[10413] | 337 | ENDIF |
---|
| 338 | |
---|
[8586] | 339 | !------------------------------------------------------------------------------! |
---|
| 340 | ! 3) Solution of the momentum equation, iterative procedure |
---|
| 341 | !------------------------------------------------------------------------------! |
---|
| 342 | ! |
---|
[11536] | 343 | ! ! ==================== ! |
---|
[8586] | 344 | DO jter = 1 , nn_nevp ! loop over jter ! |
---|
[14072] | 345 | ! ! ==================== ! |
---|
[10425] | 346 | l_full_nf_update = jter == nn_nevp ! false: disable full North fold update (performances) for iter = 1 to nn_nevp-1 |
---|
| 347 | ! |
---|
[13472] | 348 | ! convergence test |
---|
| 349 | IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2 ) THEN |
---|
| 350 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
| 351 | zu_ice(ji,jj) = u_ice(ji,jj) * umask(ji,jj,1) ! velocity at previous time step |
---|
| 352 | zv_ice(ji,jj) = v_ice(ji,jj) * vmask(ji,jj,1) |
---|
| 353 | END_2D |
---|
| 354 | ENDIF |
---|
[8586] | 355 | |
---|
| 356 | ! --- divergence, tension & shear (Appendix B of Hunke & Dukowicz, 2002) --- ! |
---|
[15548] | 357 | DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 ) |
---|
[8586] | 358 | |
---|
[12377] | 359 | ! shear at F points |
---|
| 360 | 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) & |
---|
| 361 | & + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj) & |
---|
[15548] | 362 | & ) * r1_e1e2f(ji,jj) * fimask(ji,jj) |
---|
[8586] | 363 | |
---|
[12377] | 364 | END_2D |
---|
[8586] | 365 | |
---|
[13550] | 366 | DO_2D( 0, 0, 0, 0 ) |
---|
[8586] | 367 | |
---|
[12377] | 368 | ! shear**2 at T points (doc eq. A16) |
---|
| 369 | zds2 = ( zds(ji,jj ) * zds(ji,jj ) * e1e2f(ji,jj ) + zds(ji-1,jj ) * zds(ji-1,jj ) * e1e2f(ji-1,jj ) & |
---|
| 370 | & + 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) & |
---|
| 371 | & ) * 0.25_wp * r1_e1e2t(ji,jj) |
---|
[14072] | 372 | |
---|
[12377] | 373 | ! divergence at T points |
---|
| 374 | zdiv = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) & |
---|
| 375 | & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) & |
---|
| 376 | & ) * r1_e1e2t(ji,jj) |
---|
| 377 | zdiv2 = zdiv * zdiv |
---|
[14072] | 378 | |
---|
[12377] | 379 | ! tension at T points |
---|
| 380 | 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) & |
---|
| 381 | & - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj) & |
---|
| 382 | & ) * r1_e1e2t(ji,jj) |
---|
| 383 | zdt2 = zdt * zdt |
---|
[14072] | 384 | |
---|
[12377] | 385 | ! delta at T points |
---|
[14072] | 386 | zdelta(ji,jj) = SQRT( zdiv2 + ( zdt2 + zds2 ) * z1_ecc2 ) |
---|
[8586] | 387 | |
---|
[15548] | 388 | ! P/delta at T points |
---|
| 389 | zp_delt(ji,jj) = strength(ji,jj) / ( zdelta(ji,jj) + rn_creepl ) |
---|
[8813] | 390 | |
---|
[13550] | 391 | END_2D |
---|
[15548] | 392 | CALL lbc_lnk( 'icedyn_rhg_evp', zdelta, 'T', 1.0_wp, zp_delt, 'T', 1.0_wp ) |
---|
[13550] | 393 | |
---|
[15548] | 394 | ! |
---|
| 395 | 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] | 396 | |
---|
| 397 | ! divergence at T points (duplication to avoid communications) |
---|
[15548] | 398 | ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility) |
---|
| 399 | zdiv = ( (e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj)) & |
---|
| 400 | & + (e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1)) & |
---|
[13550] | 401 | & ) * r1_e1e2t(ji,jj) |
---|
[14072] | 402 | |
---|
[13550] | 403 | ! tension at T points (duplication to avoid communications) |
---|
| 404 | zdt = ( ( u_ice(ji,jj) * r1_e2u(ji,jj) - u_ice(ji-1,jj) * r1_e2u(ji-1,jj) ) * e2t(ji,jj) * e2t(ji,jj) & |
---|
| 405 | & - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj) & |
---|
| 406 | & ) * r1_e1e2t(ji,jj) |
---|
[14072] | 407 | |
---|
[13472] | 408 | ! alpha for aEVP |
---|
[12377] | 409 | ! gamma = 0.5*P/(delta+creepl) * (c*pi)**2/Area * dt/m |
---|
| 410 | ! alpha = beta = sqrt(4*gamma) |
---|
| 411 | IF( ln_aEVP ) THEN |
---|
| 412 | zalph1 = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) ) |
---|
| 413 | z1_alph1 = 1._wp / ( zalph1 + 1._wp ) |
---|
| 414 | zalph2 = zalph1 |
---|
| 415 | z1_alph2 = z1_alph1 |
---|
[13472] | 416 | ! explicit: |
---|
| 417 | ! z1_alph1 = 1._wp / zalph1 |
---|
| 418 | ! z1_alph2 = 1._wp / zalph1 |
---|
| 419 | ! zalph1 = zalph1 - 1._wp |
---|
| 420 | ! zalph2 = zalph1 |
---|
[12377] | 421 | ENDIF |
---|
[14072] | 422 | |
---|
[12377] | 423 | ! stress at T points (zkt/=0 if landfast) |
---|
[13550] | 424 | zs1(ji,jj) = ( zs1(ji,jj)*zalph1 + zp_delt(ji,jj) * ( zdiv*(1._wp + zkt) - zdelta(ji,jj)*(1._wp - zkt) ) ) * z1_alph1 |
---|
| 425 | zs2(ji,jj) = ( zs2(ji,jj)*zalph2 + zp_delt(ji,jj) * ( zdt * z1_ecc2 * (1._wp + zkt) ) ) * z1_alph2 |
---|
[14072] | 426 | |
---|
[12377] | 427 | END_2D |
---|
[8586] | 428 | |
---|
[13472] | 429 | ! Save beta at T-points for further computations |
---|
| 430 | IF( ln_aEVP ) THEN |
---|
[15548] | 431 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
[13472] | 432 | zbeta(ji,jj) = MAX( 50._wp, rpi * SQRT( 0.5_wp * zp_delt(ji,jj) * r1_e1e2t(ji,jj) * zdt_m(ji,jj) ) ) |
---|
| 433 | END_2D |
---|
| 434 | ENDIF |
---|
[14072] | 435 | |
---|
[15548] | 436 | DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 ) |
---|
[8586] | 437 | |
---|
[13472] | 438 | ! alpha for aEVP |
---|
[12377] | 439 | IF( ln_aEVP ) THEN |
---|
[13472] | 440 | zalph2 = MAX( zbeta(ji,jj), zbeta(ji+1,jj), zbeta(ji,jj+1), zbeta(ji+1,jj+1) ) |
---|
[12377] | 441 | z1_alph2 = 1._wp / ( zalph2 + 1._wp ) |
---|
[13472] | 442 | ! explicit: |
---|
| 443 | ! z1_alph2 = 1._wp / zalph2 |
---|
| 444 | ! zalph2 = zalph2 - 1._wp |
---|
[12377] | 445 | ENDIF |
---|
[14072] | 446 | |
---|
[12377] | 447 | ! P/delta at F points |
---|
[15548] | 448 | ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility) |
---|
| 449 | 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] | 450 | |
---|
[12377] | 451 | ! stress at F points (zkt/=0 if landfast) |
---|
| 452 | zs12(ji,jj)= ( zs12(ji,jj) * zalph2 + zp_delf * ( zds(ji,jj) * z1_ecc2 * (1._wp + zkt) ) * 0.5_wp ) * z1_alph2 |
---|
[8586] | 453 | |
---|
[12377] | 454 | END_2D |
---|
[8586] | 455 | |
---|
| 456 | ! --- Ice internal stresses (Appendix C of Hunke and Dukowicz, 2002) --- ! |
---|
[15548] | 457 | ! (brackets added to fix the order of floating point operations for halo 1 - halo 2 compatibility) |
---|
| 458 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[12377] | 459 | ! !--- U points |
---|
[15548] | 460 | zfU(ji,jj) = 0.5_wp * ( (( zs1(ji+1,jj) - zs1(ji,jj) ) * e2u(ji,jj) & |
---|
[12377] | 461 | & + ( zs2(ji+1,jj) * e2t(ji+1,jj) * e2t(ji+1,jj) - zs2(ji,jj) * e2t(ji,jj) * e2t(ji,jj) & |
---|
[15548] | 462 | & ) * r1_e2u(ji,jj)) & |
---|
[12377] | 463 | & + ( zs12(ji,jj) * e1f(ji,jj) * e1f(ji,jj) - zs12(ji,jj-1) * e1f(ji,jj-1) * e1f(ji,jj-1) & |
---|
| 464 | & ) * 2._wp * r1_e1u(ji,jj) & |
---|
| 465 | & ) * r1_e1e2u(ji,jj) |
---|
| 466 | ! |
---|
| 467 | ! !--- V points |
---|
[15548] | 468 | zfV(ji,jj) = 0.5_wp * ( (( zs1(ji,jj+1) - zs1(ji,jj) ) * e1v(ji,jj) & |
---|
[12377] | 469 | & - ( zs2(ji,jj+1) * e1t(ji,jj+1) * e1t(ji,jj+1) - zs2(ji,jj) * e1t(ji,jj) * e1t(ji,jj) & |
---|
[15548] | 470 | & ) * r1_e1v(ji,jj)) & |
---|
[12377] | 471 | & + ( zs12(ji,jj) * e2f(ji,jj) * e2f(ji,jj) - zs12(ji-1,jj) * e2f(ji-1,jj) * e2f(ji-1,jj) & |
---|
| 472 | & ) * 2._wp * r1_e2v(ji,jj) & |
---|
| 473 | & ) * r1_e1e2v(ji,jj) |
---|
| 474 | ! |
---|
| 475 | ! !--- ice currents at U-V point |
---|
[15548] | 476 | 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) |
---|
| 477 | 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] | 478 | ! |
---|
| 479 | END_2D |
---|
[8586] | 480 | ! |
---|
| 481 | ! --- Computation of ice velocity --- ! |
---|
[8813] | 482 | ! Bouillon et al. 2013 (eq 47-48) => unstable unless alpha, beta vary as in Kimmritz 2016 & 2017 |
---|
[8586] | 483 | ! Bouillon et al. 2009 (eq 34-35) => stable |
---|
| 484 | IF( MOD(jter,2) == 0 ) THEN ! even iterations |
---|
| 485 | ! |
---|
[15548] | 486 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[12377] | 487 | ! !--- tau_io/(v_oce - v_ice) |
---|
| 488 | zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) ) & |
---|
| 489 | & + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) ) |
---|
| 490 | ! !--- Ocean-to-Ice stress |
---|
| 491 | ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) ) |
---|
| 492 | ! |
---|
| 493 | ! !--- tau_bottom/v_ice |
---|
| 494 | zvel = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) ) |
---|
| 495 | zTauB = ztauy_base(ji,jj) / zvel |
---|
| 496 | ! !--- OceanBottom-to-Ice stress |
---|
| 497 | ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj) |
---|
| 498 | ! |
---|
| 499 | ! !--- Coriolis at V-points (energy conserving formulation) |
---|
| 500 | zCorV(ji,jj) = - 0.25_wp * r1_e2v(ji,jj) * & |
---|
| 501 | & ( zmf(ji,jj ) * ( e2u(ji,jj ) * u_ice(ji,jj ) + e2u(ji-1,jj ) * u_ice(ji-1,jj ) ) & |
---|
| 502 | & + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) ) |
---|
| 503 | ! |
---|
| 504 | ! !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io |
---|
| 505 | zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj) |
---|
| 506 | ! |
---|
| 507 | ! !--- landfast switch => 0 = static friction : TauB > RHS & sign(TauB) /= sign(RHS) |
---|
| 508 | ! 1 = sliding friction : TauB < RHS |
---|
| 509 | rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) ) |
---|
| 510 | ! |
---|
| 511 | IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017) |
---|
[13472] | 512 | zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) ) |
---|
| 513 | v_ice(ji,jj) = ( ( rswitch * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) ) & ! previous velocity |
---|
| 514 | & + zRHS + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 515 | & ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
[14072] | 516 | & + ( 1._wp - rswitch ) * ( v_ice_b(ji,jj) & |
---|
[13472] | 517 | & + v_ice (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 518 | & ) / ( zbetav + 1._wp ) & |
---|
| 519 | & ) * 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] | 520 | & ) * zmsk00y(ji,jj) |
---|
| 521 | ELSE !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009) |
---|
[13472] | 522 | v_ice(ji,jj) = ( ( rswitch * ( zmV_t(ji,jj) * v_ice(ji,jj) & ! previous velocity |
---|
| 523 | & + zRHS + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 524 | & ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 525 | & + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 526 | & ) * 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 |
---|
| 527 | & ) * zmsk00y(ji,jj) |
---|
[12377] | 528 | ENDIF |
---|
| 529 | END_2D |
---|
[15548] | 530 | IF( nn_hls == 1 ) CALL lbc_lnk( 'icedyn_rhg_evp', v_ice, 'V', -1.0_wp ) |
---|
[8586] | 531 | ! |
---|
[13295] | 532 | DO_2D( 0, 0, 0, 0 ) |
---|
[12377] | 533 | ! !--- tau_io/(u_oce - u_ice) |
---|
| 534 | zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) ) & |
---|
| 535 | & + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) ) |
---|
| 536 | ! !--- Ocean-to-Ice stress |
---|
| 537 | ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) ) |
---|
| 538 | ! |
---|
| 539 | ! !--- tau_bottom/u_ice |
---|
| 540 | zvel = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) ) |
---|
| 541 | zTauB = ztaux_base(ji,jj) / zvel |
---|
| 542 | ! !--- OceanBottom-to-Ice stress |
---|
| 543 | ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj) |
---|
| 544 | ! |
---|
| 545 | ! !--- Coriolis at U-points (energy conserving formulation) |
---|
| 546 | zCorU(ji,jj) = 0.25_wp * r1_e1u(ji,jj) * & |
---|
| 547 | & ( zmf(ji ,jj) * ( e1v(ji ,jj) * v_ice(ji ,jj) + e1v(ji ,jj-1) * v_ice(ji ,jj-1) ) & |
---|
| 548 | & + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) ) |
---|
| 549 | ! |
---|
| 550 | ! !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io |
---|
| 551 | zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj) |
---|
| 552 | ! |
---|
| 553 | ! !--- landfast switch => 0 = static friction : TauB > RHS & sign(TauB) /= sign(RHS) |
---|
| 554 | ! 1 = sliding friction : TauB < RHS |
---|
| 555 | rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) ) |
---|
| 556 | ! |
---|
| 557 | IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017) |
---|
[13472] | 558 | zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) ) |
---|
| 559 | u_ice(ji,jj) = ( ( rswitch * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) ) & ! previous velocity |
---|
| 560 | & + zRHS + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 561 | & ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 562 | & + ( 1._wp - rswitch ) * ( u_ice_b(ji,jj) & |
---|
| 563 | & + u_ice (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 564 | & ) / ( zbetau + 1._wp ) & |
---|
[14072] | 565 | & ) * 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] | 566 | & ) * zmsk00x(ji,jj) |
---|
| 567 | ELSE !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009) |
---|
[13472] | 568 | u_ice(ji,jj) = ( ( rswitch * ( zmU_t(ji,jj) * u_ice(ji,jj) & ! previous velocity |
---|
| 569 | & + zRHS + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 570 | & ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 571 | & + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
[14072] | 572 | & ) * 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] | 573 | & ) * zmsk00x(ji,jj) |
---|
[12377] | 574 | ENDIF |
---|
| 575 | END_2D |
---|
[15548] | 576 | IF( nn_hls == 1 ) THEN ; CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp ) |
---|
| 577 | ELSE ; CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp, v_ice, 'V', -1.0_wp ) |
---|
| 578 | ENDIF |
---|
[8586] | 579 | ! |
---|
| 580 | ELSE ! odd iterations |
---|
| 581 | ! |
---|
[15548] | 582 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
[12377] | 583 | ! !--- tau_io/(u_oce - u_ice) |
---|
| 584 | zTauO = zaU(ji,jj) * zrhoco * SQRT( ( u_ice (ji,jj) - u_oce (ji,jj) ) * ( u_ice (ji,jj) - u_oce (ji,jj) ) & |
---|
| 585 | & + ( v_iceU(ji,jj) - v_oceU(ji,jj) ) * ( v_iceU(ji,jj) - v_oceU(ji,jj) ) ) |
---|
| 586 | ! !--- Ocean-to-Ice stress |
---|
| 587 | ztaux_oi(ji,jj) = zTauO * ( u_oce(ji,jj) - u_ice(ji,jj) ) |
---|
| 588 | ! |
---|
| 589 | ! !--- tau_bottom/u_ice |
---|
| 590 | zvel = 5.e-05_wp + SQRT( v_iceU(ji,jj) * v_iceU(ji,jj) + u_ice(ji,jj) * u_ice(ji,jj) ) |
---|
| 591 | zTauB = ztaux_base(ji,jj) / zvel |
---|
| 592 | ! !--- OceanBottom-to-Ice stress |
---|
| 593 | ztaux_bi(ji,jj) = zTauB * u_ice(ji,jj) |
---|
| 594 | ! |
---|
| 595 | ! !--- Coriolis at U-points (energy conserving formulation) |
---|
| 596 | zCorU(ji,jj) = 0.25_wp * r1_e1u(ji,jj) * & |
---|
| 597 | & ( zmf(ji ,jj) * ( e1v(ji ,jj) * v_ice(ji ,jj) + e1v(ji ,jj-1) * v_ice(ji ,jj-1) ) & |
---|
| 598 | & + zmf(ji+1,jj) * ( e1v(ji+1,jj) * v_ice(ji+1,jj) + e1v(ji+1,jj-1) * v_ice(ji+1,jj-1) ) ) |
---|
| 599 | ! |
---|
| 600 | ! !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io |
---|
| 601 | zRHS = zfU(ji,jj) + ztaux_ai(ji,jj) + zCorU(ji,jj) + zspgU(ji,jj) + ztaux_oi(ji,jj) |
---|
| 602 | ! |
---|
| 603 | ! !--- landfast switch => 0 = static friction : TauB > RHS & sign(TauB) /= sign(RHS) |
---|
| 604 | ! 1 = sliding friction : TauB < RHS |
---|
| 605 | rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztaux_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) ) |
---|
| 606 | ! |
---|
| 607 | IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017) |
---|
[13472] | 608 | zbetau = MAX( zbeta(ji,jj), zbeta(ji+1,jj) ) |
---|
| 609 | u_ice(ji,jj) = ( ( rswitch * ( zmU_t(ji,jj) * ( zbetau * u_ice(ji,jj) + u_ice_b(ji,jj) ) & ! previous velocity |
---|
| 610 | & + zRHS + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 611 | & ) / MAX( zepsi, zmU_t(ji,jj) * ( zbetau + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 612 | & + ( 1._wp - rswitch ) * ( u_ice_b(ji,jj) & |
---|
| 613 | & + u_ice (ji,jj) * MAX( 0._wp, zbetau - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 614 | & ) / ( zbetau + 1._wp ) & |
---|
[14072] | 615 | & ) * 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] | 616 | & ) * zmsk00x(ji,jj) |
---|
| 617 | ELSE !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009) |
---|
[13472] | 618 | u_ice(ji,jj) = ( ( rswitch * ( zmU_t(ji,jj) * u_ice(ji,jj) & ! previous velocity |
---|
| 619 | & + zRHS + zTauO * u_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 620 | & ) / MAX( zepsi, zmU_t(ji,jj) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 621 | & + ( 1._wp - rswitch ) * u_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 622 | & ) * 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 |
---|
| 623 | & ) * zmsk00x(ji,jj) |
---|
[12377] | 624 | ENDIF |
---|
| 625 | END_2D |
---|
[15548] | 626 | IF( nn_hls == 1 ) CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp ) |
---|
[8586] | 627 | ! |
---|
[13295] | 628 | DO_2D( 0, 0, 0, 0 ) |
---|
[12377] | 629 | ! !--- tau_io/(v_oce - v_ice) |
---|
| 630 | zTauO = zaV(ji,jj) * zrhoco * SQRT( ( v_ice (ji,jj) - v_oce (ji,jj) ) * ( v_ice (ji,jj) - v_oce (ji,jj) ) & |
---|
| 631 | & + ( u_iceV(ji,jj) - u_oceV(ji,jj) ) * ( u_iceV(ji,jj) - u_oceV(ji,jj) ) ) |
---|
| 632 | ! !--- Ocean-to-Ice stress |
---|
| 633 | ztauy_oi(ji,jj) = zTauO * ( v_oce(ji,jj) - v_ice(ji,jj) ) |
---|
| 634 | ! |
---|
| 635 | ! !--- tau_bottom/v_ice |
---|
| 636 | zvel = 5.e-05_wp + SQRT( v_ice(ji,jj) * v_ice(ji,jj) + u_iceV(ji,jj) * u_iceV(ji,jj) ) |
---|
| 637 | zTauB = ztauy_base(ji,jj) / zvel |
---|
| 638 | ! !--- OceanBottom-to-Ice stress |
---|
| 639 | ztauy_bi(ji,jj) = zTauB * v_ice(ji,jj) |
---|
| 640 | ! |
---|
| 641 | ! !--- Coriolis at v-points (energy conserving formulation) |
---|
| 642 | zCorV(ji,jj) = - 0.25_wp * r1_e2v(ji,jj) * & |
---|
| 643 | & ( zmf(ji,jj ) * ( e2u(ji,jj ) * u_ice(ji,jj ) + e2u(ji-1,jj ) * u_ice(ji-1,jj ) ) & |
---|
| 644 | & + zmf(ji,jj+1) * ( e2u(ji,jj+1) * u_ice(ji,jj+1) + e2u(ji-1,jj+1) * u_ice(ji-1,jj+1) ) ) |
---|
| 645 | ! |
---|
| 646 | ! !--- Sum of external forces (explicit solution) = F + tau_ia + Coriolis + spg + tau_io |
---|
| 647 | zRHS = zfV(ji,jj) + ztauy_ai(ji,jj) + zCorV(ji,jj) + zspgV(ji,jj) + ztauy_oi(ji,jj) |
---|
| 648 | ! |
---|
| 649 | ! !--- landfast switch => 0 = static friction : TauB > RHS & sign(TauB) /= sign(RHS) |
---|
| 650 | ! 1 = sliding friction : TauB < RHS |
---|
| 651 | rswitch = 1._wp - MIN( 1._wp, ABS( SIGN( 1._wp, zRHS + ztauy_base(ji,jj) ) - SIGN( 1._wp, zRHS ) ) ) |
---|
| 652 | ! |
---|
| 653 | IF( ln_aEVP ) THEN !--- ice velocity using aEVP (Kimmritz et al 2016 & 2017) |
---|
[13472] | 654 | zbetav = MAX( zbeta(ji,jj), zbeta(ji,jj+1) ) |
---|
| 655 | v_ice(ji,jj) = ( ( rswitch * ( zmV_t(ji,jj) * ( zbetav * v_ice(ji,jj) + v_ice_b(ji,jj) ) & ! previous velocity |
---|
| 656 | & + zRHS + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 657 | & ) / MAX( zepsi, zmV_t(ji,jj) * ( zbetav + 1._wp ) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 658 | & + ( 1._wp - rswitch ) * ( v_ice_b(ji,jj) & |
---|
| 659 | & + v_ice (ji,jj) * MAX( 0._wp, zbetav - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
[14072] | 660 | & ) / ( zbetav + 1._wp ) & |
---|
[13472] | 661 | & ) * 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] | 662 | & ) * zmsk00y(ji,jj) |
---|
| 663 | ELSE !--- ice velocity using EVP implicit formulation (cf Madec doc & Bouillon 2009) |
---|
[13472] | 664 | v_ice(ji,jj) = ( ( rswitch * ( zmV_t(ji,jj) * v_ice(ji,jj) & ! previous velocity |
---|
| 665 | & + zRHS + zTauO * v_ice(ji,jj) & ! F + tau_ia + Coriolis + spg + tau_io(only ocean part) |
---|
| 666 | & ) / MAX( zepsi, zmV_t(ji,jj) + zTauO - zTauB ) & ! m/dt + tau_io(only ice part) + landfast |
---|
| 667 | & + ( 1._wp - rswitch ) * v_ice(ji,jj) * MAX( 0._wp, 1._wp - zdtevp * rn_lf_relax ) & ! static friction => slow decrease to v=0 |
---|
| 668 | & ) * 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 |
---|
| 669 | & ) * zmsk00y(ji,jj) |
---|
[12377] | 670 | ENDIF |
---|
| 671 | END_2D |
---|
[15548] | 672 | IF( nn_hls == 1 ) THEN ; CALL lbc_lnk( 'icedyn_rhg_evp', v_ice, 'V', -1.0_wp ) |
---|
| 673 | ELSE ; CALL lbc_lnk( 'icedyn_rhg_evp', u_ice, 'U', -1.0_wp, v_ice, 'V', -1.0_wp ) |
---|
| 674 | ENDIF |
---|
[8586] | 675 | ! |
---|
[15548] | 676 | ENDIF |
---|
| 677 | ! |
---|
[8586] | 678 | #if defined key_agrif |
---|
[15548] | 679 | !! CALL agrif_interp_ice( 'U', jter, nn_nevp ) |
---|
| 680 | !! CALL agrif_interp_ice( 'V', jter, nn_nevp ) |
---|
| 681 | CALL agrif_interp_ice( 'U' ) |
---|
| 682 | CALL agrif_interp_ice( 'V' ) |
---|
[8586] | 683 | #endif |
---|
[15548] | 684 | IF( ln_bdy ) CALL bdy_ice_dyn( 'U' ) |
---|
| 685 | IF( ln_bdy ) CALL bdy_ice_dyn( 'V' ) |
---|
| 686 | ! |
---|
| 687 | ! convergence test |
---|
| 688 | IF( nn_rhg_chkcvg == 2 ) CALL rhg_cvg( kt, jter, nn_nevp, u_ice, v_ice, zu_ice, zv_ice, zmsk15 ) |
---|
| 689 | ! |
---|
| 690 | ! |
---|
| 691 | ! --- change strength according to advected a_i and v_i (upstream for now) --- ! |
---|
| 692 | IF( ll_advups .AND. ln_str_H79 ) THEN |
---|
[8586] | 693 | ! |
---|
[15548] | 694 | IF( jter == 1 ) THEN ! init |
---|
| 695 | ALLOCATE( za_i_ups(jpi,jpj,jpl), zv_i_ups(jpi,jpj,jpl) ) |
---|
| 696 | zdt_ups = rDt_ice / REAL( nn_nevp ) |
---|
| 697 | za_i_ups(:,:,:) = a_i(:,:,:) |
---|
| 698 | zv_i_ups(:,:,:) = v_i(:,:,:) |
---|
| 699 | ELSE |
---|
| 700 | CALL lbc_lnk( 'icedyn_rhg_evp', za_i_ups, 'T', 1.0_wp, zv_i_ups, 'T', 1.0_wp ) |
---|
| 701 | ENDIF |
---|
| 702 | ! |
---|
| 703 | CALL rhg_upstream( jter, zdt_ups, u_ice, v_ice, za_i_ups ) ! upstream advection: a_i |
---|
| 704 | CALL rhg_upstream( jter, zdt_ups, u_ice, v_ice, zv_i_ups ) ! upstream advection: v_i |
---|
| 705 | ! |
---|
[15560] | 706 | DO_2D( 0, 0, 0, 0 ) ! strength |
---|
[15548] | 707 | strength(ji,jj) = rn_pstar * SUM( zv_i_ups(ji,jj,:) ) * EXP( -rn_crhg * ( 1._wp - SUM( za_i_ups(ji,jj,:) ) ) ) |
---|
| 708 | END_2D |
---|
| 709 | ! |
---|
| 710 | IF( jter == nn_nevp ) THEN |
---|
| 711 | DEALLOCATE( za_i_ups, zv_i_ups ) |
---|
| 712 | ENDIF |
---|
[8586] | 713 | ENDIF |
---|
| 714 | ! ! ==================== ! |
---|
| 715 | END DO ! end loop over jter ! |
---|
| 716 | ! ! ==================== ! |
---|
[13472] | 717 | IF( ln_aEVP ) CALL iom_put( 'beta_evp' , zbeta ) |
---|
[8586] | 718 | ! |
---|
[15548] | 719 | IF( ll_advups .AND. ln_str_H79 ) CALL lbc_lnk( 'icedyn_rhg_evp', strength, 'T', 1.0_wp ) |
---|
| 720 | ! |
---|
[8586] | 721 | !------------------------------------------------------------------------------! |
---|
[14072] | 722 | ! 4) Recompute delta, shear and div (inputs for mechanical redistribution) |
---|
[8586] | 723 | !------------------------------------------------------------------------------! |
---|
[15548] | 724 | DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 ) |
---|
[8586] | 725 | |
---|
[12377] | 726 | ! shear at F points |
---|
| 727 | 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) & |
---|
| 728 | & + ( v_ice(ji+1,jj) * r1_e2v(ji+1,jj) - v_ice(ji,jj) * r1_e2v(ji,jj) ) * e2f(ji,jj) * e2f(ji,jj) & |
---|
[15548] | 729 | & ) * r1_e1e2f(ji,jj) * fimask(ji,jj) |
---|
[8586] | 730 | |
---|
[12377] | 731 | END_2D |
---|
[14072] | 732 | |
---|
[13497] | 733 | DO_2D( 0, 0, 0, 0 ) ! no vector loop |
---|
[14072] | 734 | |
---|
[12377] | 735 | ! tension**2 at T points |
---|
| 736 | 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) & |
---|
| 737 | & - ( v_ice(ji,jj) * r1_e1v(ji,jj) - v_ice(ji,jj-1) * r1_e1v(ji,jj-1) ) * e1t(ji,jj) * e1t(ji,jj) & |
---|
| 738 | & ) * r1_e1e2t(ji,jj) |
---|
| 739 | zdt2 = zdt * zdt |
---|
[13601] | 740 | |
---|
| 741 | zten_i(ji,jj) = zdt |
---|
[14072] | 742 | |
---|
[12377] | 743 | ! shear**2 at T points (doc eq. A16) |
---|
| 744 | zds2 = ( zds(ji,jj ) * zds(ji,jj ) * e1e2f(ji,jj ) + zds(ji-1,jj ) * zds(ji-1,jj ) * e1e2f(ji-1,jj ) & |
---|
| 745 | & + 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) & |
---|
| 746 | & ) * 0.25_wp * r1_e1e2t(ji,jj) |
---|
[14072] | 747 | |
---|
[12377] | 748 | ! shear at T points |
---|
| 749 | pshear_i(ji,jj) = SQRT( zdt2 + zds2 ) |
---|
[8586] | 750 | |
---|
[12377] | 751 | ! divergence at T points |
---|
| 752 | pdivu_i(ji,jj) = ( e2u(ji,jj) * u_ice(ji,jj) - e2u(ji-1,jj) * u_ice(ji-1,jj) & |
---|
| 753 | & + e1v(ji,jj) * v_ice(ji,jj) - e1v(ji,jj-1) * v_ice(ji,jj-1) & |
---|
| 754 | & ) * r1_e1e2t(ji,jj) |
---|
[14072] | 755 | |
---|
[12377] | 756 | ! delta at T points |
---|
[14072] | 757 | zfac = SQRT( pdivu_i(ji,jj) * pdivu_i(ji,jj) + ( zdt2 + zds2 ) * z1_ecc2 ) ! delta |
---|
[13601] | 758 | rswitch = 1._wp - MAX( 0._wp, SIGN( 1._wp, -zfac ) ) ! 0 if delta=0 |
---|
| 759 | pdelta_i(ji,jj) = zfac + rn_creepl * rswitch ! delta+creepl |
---|
[8586] | 760 | |
---|
[12377] | 761 | END_2D |
---|
[14433] | 762 | 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, & |
---|
| 763 | & zs1 , 'T', 1._wp, zs2 , 'T', 1._wp, zs12 , 'F', 1._wp ) |
---|
[14072] | 764 | |
---|
[8586] | 765 | ! --- Store the stress tensor for the next time step --- ! |
---|
| 766 | pstress1_i (:,:) = zs1 (:,:) |
---|
| 767 | pstress2_i (:,:) = zs2 (:,:) |
---|
| 768 | pstress12_i(:,:) = zs12(:,:) |
---|
| 769 | ! |
---|
| 770 | |
---|
| 771 | !------------------------------------------------------------------------------! |
---|
| 772 | ! 5) diagnostics |
---|
| 773 | !------------------------------------------------------------------------------! |
---|
[11536] | 774 | ! --- ice-ocean, ice-atm. & ice-oceanbottom(landfast) stresses --- ! |
---|
| 775 | IF( iom_use('utau_oi') .OR. iom_use('vtau_oi') .OR. iom_use('utau_ai') .OR. iom_use('vtau_ai') .OR. & |
---|
| 776 | & iom_use('utau_bi') .OR. iom_use('vtau_bi') ) THEN |
---|
| 777 | ! |
---|
[14433] | 778 | CALL lbc_lnk( 'icedyn_rhg_evp', ztaux_oi, 'U', -1.0_wp, ztauy_oi, 'V', -1.0_wp, & |
---|
| 779 | & ztaux_ai, 'U', -1.0_wp, ztauy_ai, 'V', -1.0_wp, & |
---|
| 780 | & ztaux_bi, 'U', -1.0_wp, ztauy_bi, 'V', -1.0_wp ) |
---|
[11536] | 781 | ! |
---|
| 782 | CALL iom_put( 'utau_oi' , ztaux_oi * zmsk00 ) |
---|
| 783 | CALL iom_put( 'vtau_oi' , ztauy_oi * zmsk00 ) |
---|
| 784 | CALL iom_put( 'utau_ai' , ztaux_ai * zmsk00 ) |
---|
| 785 | CALL iom_put( 'vtau_ai' , ztauy_ai * zmsk00 ) |
---|
| 786 | CALL iom_put( 'utau_bi' , ztaux_bi * zmsk00 ) |
---|
| 787 | CALL iom_put( 'vtau_bi' , ztauy_bi * zmsk00 ) |
---|
| 788 | ENDIF |
---|
[14072] | 789 | |
---|
[8586] | 790 | ! --- divergence, shear and strength --- ! |
---|
[11536] | 791 | IF( iom_use('icediv') ) CALL iom_put( 'icediv' , pdivu_i * zmsk00 ) ! divergence |
---|
| 792 | IF( iom_use('iceshe') ) CALL iom_put( 'iceshe' , pshear_i * zmsk00 ) ! shear |
---|
| 793 | IF( iom_use('icestr') ) CALL iom_put( 'icestr' , strength * zmsk00 ) ! strength |
---|
[15548] | 794 | IF( iom_use('icedlt') ) CALL iom_put( 'icedlt' , pdelta_i * zmsk00 ) ! delta |
---|
[8586] | 795 | |
---|
[13601] | 796 | ! --- Stress tensor invariants (SIMIP diags) --- ! |
---|
| 797 | IF( iom_use('normstr') .OR. iom_use('sheastr') ) THEN |
---|
[8586] | 798 | ! |
---|
[13601] | 799 | ALLOCATE( zsig_I(jpi,jpj) , zsig_II(jpi,jpj) ) |
---|
[14072] | 800 | ! |
---|
[15548] | 801 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
[14072] | 802 | |
---|
[13601] | 803 | ! Ice stresses |
---|
| 804 | ! sigma1, sigma2, sigma12 are some useful recombination of the stresses (Hunke and Dukowicz MWR 2002, Bouillon et al., OM2013) |
---|
| 805 | ! These are NOT stress tensor components, neither stress invariants, neither stress principal components |
---|
| 806 | ! I know, this can be confusing... |
---|
[14072] | 807 | zfac = strength(ji,jj) / ( pdelta_i(ji,jj) + rn_creepl ) |
---|
[13601] | 808 | zsig1 = zfac * ( pdivu_i(ji,jj) - pdelta_i(ji,jj) ) |
---|
| 809 | zsig2 = zfac * z1_ecc2 * zten_i(ji,jj) |
---|
| 810 | zsig12 = zfac * z1_ecc2 * pshear_i(ji,jj) |
---|
[14072] | 811 | |
---|
[13601] | 812 | ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008) |
---|
[15548] | 813 | zsig_I (ji,jj) = zsig1 * 0.5_wp ! 1st stress invariant, aka average normal stress, aka negative pressure |
---|
| 814 | zsig_II(ji,jj) = SQRT ( zsig2 * zsig2 * 0.25_wp + zsig12 * zsig12 ) ! 2nd '' '' , aka maximum shear stress |
---|
[14072] | 815 | |
---|
| 816 | END_2D |
---|
[13601] | 817 | ! |
---|
| 818 | ! Stress tensor invariants (normal and shear stress N/m) - SIMIP diags - definitions following Coon (1974) and Feltham (2008) |
---|
| 819 | IF( iom_use('normstr') ) CALL iom_put( 'normstr', zsig_I (:,:) * zmsk00(:,:) ) ! Normal stress |
---|
| 820 | IF( iom_use('sheastr') ) CALL iom_put( 'sheastr', zsig_II(:,:) * zmsk00(:,:) ) ! Maximum shear stress |
---|
[14072] | 821 | |
---|
[13601] | 822 | DEALLOCATE ( zsig_I, zsig_II ) |
---|
[14072] | 823 | |
---|
[13601] | 824 | ENDIF |
---|
[8586] | 825 | |
---|
[13601] | 826 | ! --- Normalized stress tensor principal components --- ! |
---|
| 827 | ! This are used to plot the normalized yield curve, see Lemieux & Dupont, 2020 |
---|
| 828 | ! Recommendation 1 : we use ice strength, not replacement pressure |
---|
| 829 | ! Recommendation 2 : need to use deformations at PREVIOUS iterate for viscosities |
---|
| 830 | IF( iom_use('sig1_pnorm') .OR. iom_use('sig2_pnorm') ) THEN |
---|
[8586] | 831 | ! |
---|
[14072] | 832 | ALLOCATE( zsig1_p(jpi,jpj) , zsig2_p(jpi,jpj) , zsig_I(jpi,jpj) , zsig_II(jpi,jpj) ) |
---|
| 833 | ! |
---|
[15548] | 834 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) |
---|
[14072] | 835 | |
---|
| 836 | ! Ice stresses computed with **viscosities** (delta, p/delta) at **previous** iterates |
---|
[13601] | 837 | ! and **deformations** at current iterates |
---|
| 838 | ! following Lemieux & Dupont (2020) |
---|
| 839 | zfac = zp_delt(ji,jj) |
---|
| 840 | zsig1 = zfac * ( pdivu_i(ji,jj) - ( zdelta(ji,jj) + rn_creepl ) ) |
---|
| 841 | zsig2 = zfac * z1_ecc2 * zten_i(ji,jj) |
---|
| 842 | zsig12 = zfac * z1_ecc2 * pshear_i(ji,jj) |
---|
[14072] | 843 | |
---|
[13601] | 844 | ! Stress invariants (sigma_I, sigma_II, Coon 1974, Feltham 2008), T-point |
---|
[15548] | 845 | zsig_I(ji,jj) = zsig1 * 0.5_wp ! 1st stress invariant, aka average normal stress, aka negative pressure |
---|
| 846 | zsig_II(ji,jj) = SQRT ( zsig2 * zsig2 * 0.25_wp + zsig12 * zsig12 ) ! 2nd '' '' , aka maximum shear stress |
---|
[14072] | 847 | |
---|
[13601] | 848 | ! Normalized principal stresses (used to display the ellipse) |
---|
| 849 | z1_strength = 1._wp / MAX( 1._wp, strength(ji,jj) ) |
---|
| 850 | zsig1_p(ji,jj) = ( zsig_I(ji,jj) + zsig_II(ji,jj) ) * z1_strength |
---|
| 851 | zsig2_p(ji,jj) = ( zsig_I(ji,jj) - zsig_II(ji,jj) ) * z1_strength |
---|
[14072] | 852 | END_2D |
---|
[8586] | 853 | ! |
---|
[14072] | 854 | CALL iom_put( 'sig1_pnorm' , zsig1_p ) |
---|
| 855 | CALL iom_put( 'sig2_pnorm' , zsig2_p ) |
---|
[11536] | 856 | |
---|
[13601] | 857 | DEALLOCATE( zsig1_p , zsig2_p , zsig_I, zsig_II ) |
---|
[14072] | 858 | |
---|
[8586] | 859 | ENDIF |
---|
[13472] | 860 | |
---|
[8586] | 861 | ! --- SIMIP --- ! |
---|
[11536] | 862 | IF( iom_use('dssh_dx') .OR. iom_use('dssh_dy') .OR. & |
---|
| 863 | & iom_use('corstrx') .OR. iom_use('corstry') .OR. iom_use('intstrx') .OR. iom_use('intstry') ) THEN |
---|
| 864 | ! |
---|
[14433] | 865 | CALL lbc_lnk( 'icedyn_rhg_evp', zspgU, 'U', -1.0_wp, zspgV, 'V', -1.0_wp, & |
---|
| 866 | & zCorU, 'U', -1.0_wp, zCorV, 'V', -1.0_wp, zfU, 'U', -1.0_wp, zfV, 'V', -1.0_wp ) |
---|
[8586] | 867 | |
---|
[11536] | 868 | CALL iom_put( 'dssh_dx' , zspgU * zmsk00 ) ! Sea-surface tilt term in force balance (x) |
---|
| 869 | CALL iom_put( 'dssh_dy' , zspgV * zmsk00 ) ! Sea-surface tilt term in force balance (y) |
---|
| 870 | CALL iom_put( 'corstrx' , zCorU * zmsk00 ) ! Coriolis force term in force balance (x) |
---|
| 871 | CALL iom_put( 'corstry' , zCorV * zmsk00 ) ! Coriolis force term in force balance (y) |
---|
| 872 | CALL iom_put( 'intstrx' , zfU * zmsk00 ) ! Internal force term in force balance (x) |
---|
| 873 | CALL iom_put( 'intstry' , zfV * zmsk00 ) ! Internal force term in force balance (y) |
---|
| 874 | ENDIF |
---|
| 875 | |
---|
| 876 | IF( iom_use('xmtrpice') .OR. iom_use('ymtrpice') .OR. & |
---|
| 877 | & iom_use('xmtrpsnw') .OR. iom_use('ymtrpsnw') .OR. iom_use('xatrp') .OR. iom_use('yatrp') ) THEN |
---|
| 878 | ! |
---|
| 879 | ALLOCATE( zdiag_xmtrp_ice(jpi,jpj) , zdiag_ymtrp_ice(jpi,jpj) , & |
---|
| 880 | & zdiag_xmtrp_snw(jpi,jpj) , zdiag_ymtrp_snw(jpi,jpj) , zdiag_xatrp(jpi,jpj) , zdiag_yatrp(jpi,jpj) ) |
---|
| 881 | ! |
---|
[13295] | 882 | DO_2D( 0, 0, 0, 0 ) |
---|
[12377] | 883 | ! 2D ice mass, snow mass, area transport arrays (X, Y) |
---|
| 884 | zfac_x = 0.5 * u_ice(ji,jj) * e2u(ji,jj) * zmsk00(ji,jj) |
---|
| 885 | zfac_y = 0.5 * v_ice(ji,jj) * e1v(ji,jj) * zmsk00(ji,jj) |
---|
[11536] | 886 | |
---|
[12377] | 887 | zdiag_xmtrp_ice(ji,jj) = rhoi * zfac_x * ( vt_i(ji+1,jj) + vt_i(ji,jj) ) ! ice mass transport, X-component |
---|
| 888 | zdiag_ymtrp_ice(ji,jj) = rhoi * zfac_y * ( vt_i(ji,jj+1) + vt_i(ji,jj) ) ! '' Y- '' |
---|
[11536] | 889 | |
---|
[12377] | 890 | zdiag_xmtrp_snw(ji,jj) = rhos * zfac_x * ( vt_s(ji+1,jj) + vt_s(ji,jj) ) ! snow mass transport, X-component |
---|
| 891 | zdiag_ymtrp_snw(ji,jj) = rhos * zfac_y * ( vt_s(ji,jj+1) + vt_s(ji,jj) ) ! '' Y- '' |
---|
[11536] | 892 | |
---|
[12377] | 893 | zdiag_xatrp(ji,jj) = zfac_x * ( at_i(ji+1,jj) + at_i(ji,jj) ) ! area transport, X-component |
---|
| 894 | zdiag_yatrp(ji,jj) = zfac_y * ( at_i(ji,jj+1) + at_i(ji,jj) ) ! '' Y- '' |
---|
[11536] | 895 | |
---|
[12377] | 896 | END_2D |
---|
[8586] | 897 | |
---|
[14433] | 898 | CALL lbc_lnk( 'icedyn_rhg_evp', zdiag_xmtrp_ice, 'U', -1.0_wp, zdiag_ymtrp_ice, 'V', -1.0_wp, & |
---|
| 899 | & zdiag_xmtrp_snw, 'U', -1.0_wp, zdiag_ymtrp_snw, 'V', -1.0_wp, & |
---|
| 900 | & zdiag_xatrp , 'U', -1.0_wp, zdiag_yatrp , 'V', -1.0_wp ) |
---|
[8586] | 901 | |
---|
[11536] | 902 | CALL iom_put( 'xmtrpice' , zdiag_xmtrp_ice ) ! X-component of sea-ice mass transport (kg/s) |
---|
[14072] | 903 | CALL iom_put( 'ymtrpice' , zdiag_ymtrp_ice ) ! Y-component of sea-ice mass transport |
---|
[11536] | 904 | CALL iom_put( 'xmtrpsnw' , zdiag_xmtrp_snw ) ! X-component of snow mass transport (kg/s) |
---|
| 905 | CALL iom_put( 'ymtrpsnw' , zdiag_ymtrp_snw ) ! Y-component of snow mass transport |
---|
| 906 | CALL iom_put( 'xatrp' , zdiag_xatrp ) ! X-component of ice area transport |
---|
| 907 | CALL iom_put( 'yatrp' , zdiag_yatrp ) ! Y-component of ice area transport |
---|
| 908 | |
---|
| 909 | DEALLOCATE( zdiag_xmtrp_ice , zdiag_ymtrp_ice , & |
---|
| 910 | & zdiag_xmtrp_snw , zdiag_ymtrp_snw , zdiag_xatrp , zdiag_yatrp ) |
---|
| 911 | |
---|
[8586] | 912 | ENDIF |
---|
| 913 | ! |
---|
[13472] | 914 | ! --- convergence tests --- ! |
---|
| 915 | IF( nn_rhg_chkcvg == 1 .OR. nn_rhg_chkcvg == 2 ) THEN |
---|
| 916 | IF( iom_use('uice_cvg') ) THEN |
---|
| 917 | IF( ln_aEVP ) THEN ! output: beta * ( u(t=nn_nevp) - u(t=nn_nevp-1) ) |
---|
| 918 | CALL iom_put( 'uice_cvg', MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * zbeta(:,:) * umask(:,:,1) , & |
---|
| 919 | & ABS( v_ice(:,:) - zv_ice(:,:) ) * zbeta(:,:) * vmask(:,:,1) ) * zmsk15(:,:) ) |
---|
| 920 | ELSE ! output: nn_nevp * ( u(t=nn_nevp) - u(t=nn_nevp-1) ) |
---|
| 921 | CALL iom_put( 'uice_cvg', REAL( nn_nevp ) * MAX( ABS( u_ice(:,:) - zu_ice(:,:) ) * umask(:,:,1) , & |
---|
| 922 | & ABS( v_ice(:,:) - zv_ice(:,:) ) * vmask(:,:,1) ) * zmsk15(:,:) ) |
---|
| 923 | ENDIF |
---|
| 924 | ENDIF |
---|
[14072] | 925 | ENDIF |
---|
[13472] | 926 | ! |
---|
[8586] | 927 | END SUBROUTINE ice_dyn_rhg_evp |
---|
| 928 | |
---|
[8813] | 929 | |
---|
[15548] | 930 | SUBROUTINE rhg_cvg( kt, kiter, kitermax, pu, pv, pub, pvb, pmsk15 ) |
---|
[13472] | 931 | !!---------------------------------------------------------------------- |
---|
| 932 | !! *** ROUTINE rhg_cvg *** |
---|
[14072] | 933 | !! |
---|
[13472] | 934 | !! ** Purpose : check convergence of oce rheology |
---|
| 935 | !! |
---|
| 936 | !! ** Method : create a file ice_cvg.nc containing the convergence of ice velocity |
---|
| 937 | !! during the sub timestepping of rheology so as: |
---|
| 938 | !! uice_cvg = MAX( u(t+1) - u(t) , v(t+1) - v(t) ) |
---|
| 939 | !! This routine is called every sub-iteration, so it is cpu expensive |
---|
| 940 | !! |
---|
[14072] | 941 | !! ** Note : for the first sub-iteration, uice_cvg is set to 0 (too large otherwise) |
---|
[13472] | 942 | !!---------------------------------------------------------------------- |
---|
| 943 | INTEGER , INTENT(in) :: kt, kiter, kitermax ! ocean time-step index |
---|
| 944 | REAL(wp), DIMENSION(:,:), INTENT(in) :: pu, pv, pub, pvb ! now and before velocities |
---|
[15548] | 945 | REAL(wp), DIMENSION(:,:), INTENT(in) :: pmsk15 |
---|
[13472] | 946 | !! |
---|
| 947 | INTEGER :: it, idtime, istatus |
---|
| 948 | INTEGER :: ji, jj ! dummy loop indices |
---|
[14072] | 949 | REAL(wp) :: zresm ! local real |
---|
[13472] | 950 | CHARACTER(len=20) :: clname |
---|
[15548] | 951 | LOGICAL :: ll_maxcvg |
---|
| 952 | REAL(wp), DIMENSION(jpi,jpj,2) :: zres |
---|
| 953 | REAL(wp), DIMENSION(2) :: ztmp |
---|
[13472] | 954 | !!---------------------------------------------------------------------- |
---|
[15548] | 955 | ll_maxcvg = .FALSE. |
---|
| 956 | ! |
---|
[13472] | 957 | ! create file |
---|
| 958 | IF( kt == nit000 .AND. kiter == 1 ) THEN |
---|
| 959 | ! |
---|
| 960 | IF( lwp ) THEN |
---|
| 961 | WRITE(numout,*) |
---|
| 962 | WRITE(numout,*) 'rhg_cvg : ice rheology convergence control' |
---|
| 963 | WRITE(numout,*) '~~~~~~~' |
---|
| 964 | ENDIF |
---|
| 965 | ! |
---|
| 966 | IF( lwm ) THEN |
---|
| 967 | clname = 'ice_cvg.nc' |
---|
| 968 | IF( .NOT. Agrif_Root() ) clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname) |
---|
| 969 | istatus = NF90_CREATE( TRIM(clname), NF90_CLOBBER, ncvgid ) |
---|
| 970 | istatus = NF90_DEF_DIM( ncvgid, 'time' , NF90_UNLIMITED, idtime ) |
---|
[15548] | 971 | istatus = NF90_DEF_VAR( ncvgid, 'uice_cvg', NF90_DOUBLE , (/ idtime /), nvarid ) |
---|
[13472] | 972 | istatus = NF90_ENDDEF(ncvgid) |
---|
| 973 | ENDIF |
---|
| 974 | ! |
---|
| 975 | ENDIF |
---|
| 976 | |
---|
| 977 | ! time |
---|
[15548] | 978 | it = ( kt - nit000 ) * kitermax + kiter |
---|
[14072] | 979 | |
---|
[13472] | 980 | ! convergence |
---|
| 981 | IF( kiter == 1 ) THEN ! remove the first iteration for calculations of convergence (always very large) |
---|
| 982 | zresm = 0._wp |
---|
| 983 | ELSE |
---|
[15548] | 984 | zresm = 0._wp |
---|
| 985 | IF( ll_maxcvg ) THEN ! error max over the domain |
---|
| 986 | DO_2D( 0, 0, 0, 0 ) |
---|
| 987 | zresm = MAX( zresm, MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), & |
---|
| 988 | & ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * pmsk15(ji,jj) ) |
---|
| 989 | END_2D |
---|
| 990 | CALL mpp_max( 'icedyn_rhg_evp', zresm ) |
---|
| 991 | ELSE ! error averaged over the domain |
---|
| 992 | DO_2D( 0, 0, 0, 0 ) |
---|
| 993 | zres(ji,jj,1) = MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), & |
---|
| 994 | & ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * pmsk15(ji,jj) |
---|
| 995 | zres(ji,jj,2) = pmsk15(ji,jj) |
---|
| 996 | END_2D |
---|
| 997 | ztmp(:) = glob_sum_vec( 'icedyn_rhg_evp', zres ) |
---|
| 998 | IF( ztmp(2) /= 0._wp ) zresm = ztmp(1) / ztmp(2) |
---|
| 999 | ENDIF |
---|
[13472] | 1000 | ENDIF |
---|
| 1001 | |
---|
| 1002 | IF( lwm ) THEN |
---|
| 1003 | ! write variables |
---|
| 1004 | istatus = NF90_PUT_VAR( ncvgid, nvarid, (/zresm/), (/it/), (/1/) ) |
---|
| 1005 | ! close file |
---|
[15548] | 1006 | IF( kt == nitend - nn_fsbc + 1 .AND. kiter == kitermax ) istatus = NF90_CLOSE(ncvgid) |
---|
[13472] | 1007 | ENDIF |
---|
[14072] | 1008 | |
---|
[13472] | 1009 | END SUBROUTINE rhg_cvg |
---|
| 1010 | |
---|
| 1011 | |
---|
[8586] | 1012 | SUBROUTINE rhg_evp_rst( cdrw, kt ) |
---|
| 1013 | !!--------------------------------------------------------------------- |
---|
| 1014 | !! *** ROUTINE rhg_evp_rst *** |
---|
[14072] | 1015 | !! |
---|
[8586] | 1016 | !! ** Purpose : Read or write RHG file in restart file |
---|
| 1017 | !! |
---|
| 1018 | !! ** Method : use of IOM library |
---|
| 1019 | !!---------------------------------------------------------------------- |
---|
| 1020 | CHARACTER(len=*) , INTENT(in) :: cdrw ! "READ"/"WRITE" flag |
---|
| 1021 | INTEGER, OPTIONAL, INTENT(in) :: kt ! ice time-step |
---|
| 1022 | ! |
---|
| 1023 | INTEGER :: iter ! local integer |
---|
| 1024 | INTEGER :: id1, id2, id3 ! local integers |
---|
| 1025 | !!---------------------------------------------------------------------- |
---|
| 1026 | ! |
---|
| 1027 | IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialize |
---|
| 1028 | ! ! --------------- |
---|
| 1029 | IF( ln_rstart ) THEN !* Read the restart file |
---|
| 1030 | ! |
---|
| 1031 | id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. ) |
---|
| 1032 | id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. ) |
---|
| 1033 | id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. ) |
---|
| 1034 | ! |
---|
| 1035 | IF( MIN( id1, id2, id3 ) > 0 ) THEN ! fields exist |
---|
[13286] | 1036 | CALL iom_get( numrir, jpdom_auto, 'stress1_i' , stress1_i , cd_type = 'T' ) |
---|
| 1037 | CALL iom_get( numrir, jpdom_auto, 'stress2_i' , stress2_i , cd_type = 'T' ) |
---|
| 1038 | CALL iom_get( numrir, jpdom_auto, 'stress12_i', stress12_i, cd_type = 'F' ) |
---|
[8586] | 1039 | ELSE ! start rheology from rest |
---|
[9169] | 1040 | IF(lwp) WRITE(numout,*) |
---|
| 1041 | IF(lwp) WRITE(numout,*) ' ==>>> previous run without rheology, set stresses to 0' |
---|
[8586] | 1042 | stress1_i (:,:) = 0._wp |
---|
| 1043 | stress2_i (:,:) = 0._wp |
---|
| 1044 | stress12_i(:,:) = 0._wp |
---|
| 1045 | ENDIF |
---|
| 1046 | ELSE !* Start from rest |
---|
[9169] | 1047 | IF(lwp) WRITE(numout,*) |
---|
| 1048 | IF(lwp) WRITE(numout,*) ' ==>>> start from rest: set stresses to 0' |
---|
[8586] | 1049 | stress1_i (:,:) = 0._wp |
---|
| 1050 | stress2_i (:,:) = 0._wp |
---|
| 1051 | stress12_i(:,:) = 0._wp |
---|
| 1052 | ENDIF |
---|
| 1053 | ! |
---|
| 1054 | ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file |
---|
| 1055 | ! ! ------------------- |
---|
| 1056 | IF(lwp) WRITE(numout,*) '---- rhg-rst ----' |
---|
| 1057 | iter = kt + nn_fsbc - 1 ! ice restarts are written at kt == nitrst - nn_fsbc + 1 |
---|
| 1058 | ! |
---|
[13970] | 1059 | CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i ) |
---|
| 1060 | CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i ) |
---|
[8586] | 1061 | CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i ) |
---|
| 1062 | ! |
---|
| 1063 | ENDIF |
---|
| 1064 | ! |
---|
| 1065 | END SUBROUTINE rhg_evp_rst |
---|
| 1066 | |
---|
[15548] | 1067 | SUBROUTINE rhg_upstream( jter, pdt, pu, pv, pt ) |
---|
| 1068 | !!--------------------------------------------------------------------- |
---|
| 1069 | !! *** ROUTINE rhg_upstream *** |
---|
| 1070 | !! |
---|
| 1071 | !! ** Purpose : compute the upstream fluxes and upstream guess of tracer |
---|
| 1072 | !!---------------------------------------------------------------------- |
---|
| 1073 | INTEGER , INTENT(in ) :: jter |
---|
| 1074 | REAL(wp) , INTENT(in ) :: pdt ! tracer time-step |
---|
| 1075 | REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu, pv ! 2 ice velocity components |
---|
| 1076 | REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pt ! tracer fields |
---|
| 1077 | ! |
---|
| 1078 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
| 1079 | REAL(wp) :: ztra ! local scalar |
---|
| 1080 | LOGICAL :: ll_upsxy = .TRUE. |
---|
| 1081 | REAL(wp), DIMENSION(jpi,jpj) :: zfu_ups, zfv_ups, zpt ! upstream fluxes and tracer guess |
---|
| 1082 | !!---------------------------------------------------------------------- |
---|
| 1083 | DO jl = 1, jpl |
---|
| 1084 | IF( .NOT. ll_upsxy ) THEN !** no alternate directions **! |
---|
| 1085 | ! |
---|
| 1086 | DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 ) |
---|
| 1087 | zfu_ups(ji,jj) = MAX(pu(ji,jj)*e2u(ji,jj), 0._wp) * pt(ji,jj,jl) + MIN(pu(ji,jj)*e2u(ji,jj), 0._wp) * pt(ji+1,jj,jl) |
---|
| 1088 | zfv_ups(ji,jj) = MAX(pv(ji,jj)*e1v(ji,jj), 0._wp) * pt(ji,jj,jl) + MIN(pv(ji,jj)*e1v(ji,jj), 0._wp) * pt(ji,jj+1,jl) |
---|
| 1089 | END_2D |
---|
| 1090 | ! |
---|
| 1091 | ELSE !** alternate directions **! |
---|
| 1092 | ! |
---|
| 1093 | IF( MOD(jter,2) == 1 ) THEN !== odd ice time step: adv_x then adv_y ==! |
---|
| 1094 | ! |
---|
| 1095 | DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls ) !-- flux in x-direction |
---|
| 1096 | zfu_ups(ji,jj) = MAX( pu(ji,jj)*e2u(ji,jj), 0._wp ) * pt(ji ,jj,jl) + & |
---|
| 1097 | & MIN( pu(ji,jj)*e2u(ji,jj), 0._wp ) * pt(ji+1,jj,jl) |
---|
| 1098 | END_2D |
---|
| 1099 | ! |
---|
| 1100 | DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls ) !-- first guess of tracer from u-flux |
---|
| 1101 | ztra = - ( zfu_ups(ji,jj) - zfu_ups(ji-1,jj) ) |
---|
| 1102 | zpt(ji,jj) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1) |
---|
| 1103 | END_2D |
---|
| 1104 | ! |
---|
| 1105 | DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls-1 ) !-- flux in y-direction |
---|
| 1106 | zfv_ups(ji,jj) = MAX( pv(ji,jj)*e1v(ji,jj), 0._wp ) * zpt(ji,jj ) + & |
---|
| 1107 | & MIN( pv(ji,jj)*e1v(ji,jj), 0._wp ) * zpt(ji,jj+1) |
---|
| 1108 | END_2D |
---|
| 1109 | ! |
---|
| 1110 | ELSE !== even ice time step: adv_y then adv_x ==! |
---|
| 1111 | ! |
---|
| 1112 | DO_2D( nn_hls, nn_hls, nn_hls, nn_hls-1 ) !-- flux in y-direction |
---|
| 1113 | zfv_ups(ji,jj) = MAX( pv(ji,jj)*e1v(ji,jj), 0._wp ) * pt(ji,jj ,jl) + & |
---|
| 1114 | & MIN( pv(ji,jj)*e1v(ji,jj), 0._wp ) * pt(ji,jj+1,jl) |
---|
| 1115 | END_2D |
---|
| 1116 | ! |
---|
| 1117 | DO_2D( nn_hls, nn_hls, nn_hls-1, nn_hls-1 ) !-- first guess of tracer from v-flux |
---|
| 1118 | ztra = - ( zfv_ups(ji,jj) - zfv_ups(ji,jj-1) ) |
---|
| 1119 | zpt(ji,jj) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1) |
---|
| 1120 | END_2D |
---|
| 1121 | ! |
---|
| 1122 | DO_2D( nn_hls, nn_hls-1, nn_hls-1, nn_hls-1 ) !-- flux in x-direction |
---|
| 1123 | zfu_ups(ji,jj) = MAX( pu(ji,jj)*e2u(ji,jj), 0._wp ) * zpt(ji ,jj) + & |
---|
| 1124 | & MIN( pu(ji,jj)*e2u(ji,jj), 0._wp ) * zpt(ji+1,jj) |
---|
| 1125 | END_2D |
---|
| 1126 | ! |
---|
| 1127 | ENDIF |
---|
| 1128 | ! |
---|
| 1129 | ENDIF |
---|
| 1130 | ! |
---|
| 1131 | DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 ) |
---|
| 1132 | ztra = - ( zfu_ups(ji,jj) - zfu_ups(ji-1,jj) + zfv_ups(ji,jj) - zfv_ups(ji,jj-1) ) |
---|
| 1133 | pt(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1) |
---|
| 1134 | END_2D |
---|
| 1135 | END DO |
---|
| 1136 | ! |
---|
| 1137 | END SUBROUTINE rhg_upstream |
---|
[14072] | 1138 | |
---|
[8586] | 1139 | #else |
---|
| 1140 | !!---------------------------------------------------------------------- |
---|
[9570] | 1141 | !! Default option Empty module NO SI3 sea-ice model |
---|
[8586] | 1142 | !!---------------------------------------------------------------------- |
---|
| 1143 | #endif |
---|
| 1144 | |
---|
| 1145 | !!============================================================================== |
---|
| 1146 | END MODULE icedyn_rhg_evp |
---|