bdyice_lim.F90 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408
  1. MODULE bdyice_lim
  2. !!======================================================================
  3. !! *** MODULE bdyice_lim ***
  4. !! Unstructured Open Boundary Cond. : Open boundary conditions for sea-ice (LIM2 and LIM3)
  5. !!======================================================================
  6. !! History : 3.3 ! 2010-09 (D. Storkey) Original code
  7. !! 3.4 ! 2011 (D. Storkey) rewrite in preparation for OBC-BDY merge
  8. !! - ! 2012-01 (C. Rousset) add lim3 and remove useless jk loop
  9. !!----------------------------------------------------------------------
  10. #if defined key_bdy && ( defined key_lim2 || defined key_lim3 )
  11. !!----------------------------------------------------------------------
  12. !! 'key_bdy' and Unstructured Open Boundary Conditions
  13. !! 'key_lim2' LIM-2 sea ice model
  14. !! 'key_lim3' LIM-3 sea ice model
  15. !!----------------------------------------------------------------------
  16. !! bdy_ice_lim : Application of open boundaries to ice
  17. !! bdy_ice_frs : Application of Flow Relaxation Scheme
  18. !!----------------------------------------------------------------------
  19. USE timing ! Timing
  20. USE phycst ! physical constant
  21. USE eosbn2 ! equation of state
  22. USE oce ! ocean dynamics and tracers variables
  23. #if defined key_lim2
  24. USE par_ice_2
  25. USE ice_2 ! LIM_2 ice variables
  26. USE dom_ice_2 ! sea-ice domain
  27. #elif defined key_lim3
  28. USE ice ! LIM_3 ice variables
  29. USE dom_ice ! sea-ice domain
  30. USE limvar
  31. #endif
  32. USE par_oce ! ocean parameters
  33. USE dom_oce ! ocean space and time domain variables
  34. USE sbc_oce ! Surface boundary condition: ocean fields
  35. USE bdy_oce ! ocean open boundary conditions
  36. USE lbclnk ! ocean lateral boundary conditions (or mpp link)
  37. USE in_out_manager ! write to numout file
  38. USE lib_mpp ! distributed memory computing
  39. USE lib_fortran ! to use key_nosignedzero
  40. IMPLICIT NONE
  41. PRIVATE
  42. PUBLIC bdy_ice_lim ! routine called in sbcmod
  43. PUBLIC bdy_ice_lim_dyn ! routine called in limrhg
  44. !!----------------------------------------------------------------------
  45. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  46. !! $Id: bdyice_lim.F90 2750 2016-01-12 10:42:05Z ufla $
  47. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  48. !!----------------------------------------------------------------------
  49. CONTAINS
  50. SUBROUTINE bdy_ice_lim( kt )
  51. !!----------------------------------------------------------------------
  52. !! *** SUBROUTINE bdy_ice_lim ***
  53. !!
  54. !! ** Purpose : - Apply open boundary conditions for ice (LIM2 and LIM3)
  55. !!
  56. !!----------------------------------------------------------------------
  57. INTEGER, INTENT( in ) :: kt ! Main time step counter
  58. INTEGER :: ib_bdy ! Loop index
  59. #if defined key_lim3
  60. CALL lim_var_glo2eqv
  61. #endif
  62. DO ib_bdy=1, nb_bdy
  63. SELECT CASE( cn_ice_lim(ib_bdy) )
  64. CASE('none')
  65. CYCLE
  66. CASE('frs')
  67. CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
  68. CASE DEFAULT
  69. CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' )
  70. END SELECT
  71. END DO
  72. #if defined key_lim3
  73. CALL lim_var_zapsmall
  74. CALL lim_var_agg(1)
  75. #endif
  76. END SUBROUTINE bdy_ice_lim
  77. SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy )
  78. !!------------------------------------------------------------------------------
  79. !! *** SUBROUTINE bdy_ice_frs ***
  80. !!
  81. !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields in the case
  82. !! of unstructured open boundaries.
  83. !!
  84. !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three-
  85. !! dimensional baroclinic ocean model with realistic topography. Tellus, 365-382.
  86. !!------------------------------------------------------------------------------
  87. TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices
  88. TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data
  89. INTEGER, INTENT(in) :: kt ! main time-step counter
  90. INTEGER, INTENT(in) :: ib_bdy ! BDY set index
  91. INTEGER :: jpbound ! 0 = incoming ice
  92. ! 1 = outgoing ice
  93. INTEGER :: jb, jk, jgrd, jl ! dummy loop indices
  94. INTEGER :: ji, jj, ii, ij ! local scalar
  95. REAL(wp) :: zwgt, zwgt1 ! local scalar
  96. REAL(wp) :: ztmelts, zdh
  97. #if defined key_lim2 && ! defined key_lim2_vp && defined key_agrif
  98. USE ice_2, vt_s => hsnm
  99. USE ice_2, vt_i => hicm
  100. #endif
  101. !!------------------------------------------------------------------------------
  102. !
  103. IF( nn_timing == 1 ) CALL timing_start('bdy_ice_frs')
  104. !
  105. jgrd = 1 ! Everything is at T-points here
  106. !
  107. #if defined key_lim2
  108. DO jb = 1, idx%nblenrim(jgrd)
  109. ji = idx%nbi(jb,jgrd)
  110. jj = idx%nbj(jb,jgrd)
  111. zwgt = idx%nbw(jb,jgrd)
  112. zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
  113. frld (ji,jj) = ( frld (ji,jj) * zwgt1 + dta%frld (jb) * zwgt ) * tmask(ji,jj,1) ! Leads fraction
  114. hicif(ji,jj) = ( hicif(ji,jj) * zwgt1 + dta%hicif(jb) * zwgt ) * tmask(ji,jj,1) ! Ice depth
  115. hsnif(ji,jj) = ( hsnif(ji,jj) * zwgt1 + dta%hsnif(jb) * zwgt ) * tmask(ji,jj,1) ! Snow depth
  116. END DO
  117. CALL lbc_bdy_lnk( frld, 'T', 1., ib_bdy ) ! lateral boundary conditions
  118. CALL lbc_bdy_lnk( hicif, 'T', 1., ib_bdy )
  119. CALL lbc_bdy_lnk( hsnif, 'T', 1., ib_bdy )
  120. vt_i(:,:) = hicif(:,:) * frld(:,:)
  121. vt_s(:,:) = hsnif(:,:) * frld(:,:)
  122. !
  123. #elif defined key_lim3
  124. DO jl = 1, jpl
  125. DO jb = 1, idx%nblenrim(jgrd)
  126. ji = idx%nbi(jb,jgrd)
  127. jj = idx%nbj(jb,jgrd)
  128. zwgt = idx%nbw(jb,jgrd)
  129. zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
  130. a_i (ji,jj,jl) = ( a_i (ji,jj,jl) * zwgt1 + dta%a_i (jb,jl) * zwgt ) * tmask(ji,jj,1) ! Leads fraction
  131. ht_i(ji,jj,jl) = ( ht_i(ji,jj,jl) * zwgt1 + dta%ht_i(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Ice depth
  132. ht_s(ji,jj,jl) = ( ht_s(ji,jj,jl) * zwgt1 + dta%ht_s(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Snow depth
  133. ! -----------------
  134. ! Pathological case
  135. ! -----------------
  136. ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to
  137. ! very large transformation from snow to ice (see limthd_dh.F90)
  138. ! Then, a) transfer the snow excess into the ice (different from limthd_dh)
  139. zdh = MAX( 0._wp, ( rhosn * ht_s(ji,jj,jl) + ( rhoic - rau0 ) * ht_i(ji,jj,jl) ) * r1_rau0 )
  140. ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment)
  141. !zdh = MAX( 0._wp, ht_s(ji,jj,jl) * rhosn / rhoic )
  142. ! recompute ht_i, ht_s
  143. ht_i(ji,jj,jl) = MIN( hi_max(jl), ht_i(ji,jj,jl) + zdh )
  144. ht_s(ji,jj,jl) = MAX( 0._wp, ht_s(ji,jj,jl) - zdh * rhoic / rhosn )
  145. ENDDO
  146. CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy )
  147. CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
  148. CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
  149. ENDDO
  150. ! retrieve at_i
  151. at_i(:,:) = 0._wp
  152. DO jl = 1, jpl
  153. at_i(:,:) = a_i(:,:,jl) + at_i(:,:)
  154. END DO
  155. DO jl = 1, jpl
  156. DO jb = 1, idx%nblenrim(jgrd)
  157. ji = idx%nbi(jb,jgrd)
  158. jj = idx%nbj(jb,jgrd)
  159. ! condition on ice thickness depends on the ice velocity
  160. ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values
  161. jpbound = 0; ii = ji; ij = jj;
  162. IF( u_ice(ji+1,jj ) < 0. .AND. umask(ji-1,jj ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj
  163. IF( u_ice(ji-1,jj ) > 0. .AND. umask(ji+1,jj ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj
  164. IF( v_ice(ji ,jj+1) < 0. .AND. vmask(ji ,jj-1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj+1
  165. IF( v_ice(ji ,jj-1) > 0. .AND. vmask(ji ,jj+1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj-1
  166. IF( nn_ice_lim_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj ! case ice boundaries = initial conditions
  167. ! do not make state variables dependent on velocity
  168. rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ii,ij) - 0.01 ) ) ! 0 if no ice
  169. ! concentration and thickness
  170. a_i (ji,jj,jl) = a_i (ii,ij,jl) * rswitch
  171. ht_i(ji,jj,jl) = ht_i(ii,ij,jl) * rswitch
  172. ht_s(ji,jj,jl) = ht_s(ii,ij,jl) * rswitch
  173. ! Ice and snow volumes
  174. v_i(ji,jj,jl) = ht_i(ji,jj,jl) * a_i(ji,jj,jl)
  175. v_s(ji,jj,jl) = ht_s(ji,jj,jl) * a_i(ji,jj,jl)
  176. SELECT CASE( jpbound )
  177. CASE( 0 ) ! velocity is inward
  178. ! Ice salinity, age, temperature
  179. sm_i(ji,jj,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin
  180. oa_i(ji,jj,jl) = rswitch * rn_ice_age(ib_bdy) * a_i(ji,jj,jl)
  181. t_su(ji,jj,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rn_ice_tem(ib_bdy)
  182. DO jk = 1, nlay_s
  183. t_s(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0
  184. END DO
  185. DO jk = 1, nlay_i
  186. t_i(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0
  187. s_i(ji,jj,jk,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin
  188. END DO
  189. CASE( 1 ) ! velocity is outward
  190. ! Ice salinity, age, temperature
  191. sm_i(ji,jj,jl) = rswitch * sm_i(ii,ij,jl) + ( 1.0 - rswitch ) * rn_simin
  192. oa_i(ji,jj,jl) = rswitch * oa_i(ii,ij,jl)
  193. t_su(ji,jj,jl) = rswitch * t_su(ii,ij,jl) + ( 1.0 - rswitch ) * rt0
  194. DO jk = 1, nlay_s
  195. t_s(ji,jj,jk,jl) = rswitch * t_s(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
  196. END DO
  197. DO jk = 1, nlay_i
  198. t_i(ji,jj,jk,jl) = rswitch * t_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
  199. s_i(ji,jj,jk,jl) = rswitch * s_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rn_simin
  200. END DO
  201. END SELECT
  202. ! if salinity is constant, then overwrite rn_ice_sal
  203. IF( nn_icesal == 1 ) THEN
  204. sm_i(ji,jj,jl) = rn_icesal
  205. s_i (ji,jj,:,jl) = rn_icesal
  206. ENDIF
  207. ! contents
  208. smv_i(ji,jj,jl) = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl)
  209. DO jk = 1, nlay_s
  210. ! Snow energy of melting
  211. e_s(ji,jj,jk,jl) = rswitch * rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus )
  212. ! Multiply by volume, so that heat content in J/m2
  213. e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s
  214. END DO
  215. DO jk = 1, nlay_i
  216. ztmelts = - tmut * s_i(ji,jj,jk,jl) + rt0 !Melting temperature in K
  217. ! heat content per unit volume
  218. e_i(ji,jj,jk,jl) = rswitch * rhoic * &
  219. ( cpic * ( ztmelts - t_i(ji,jj,jk,jl) ) &
  220. + lfus * ( 1.0 - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) &
  221. - rcp * ( ztmelts - rt0 ) )
  222. ! Mutliply by ice volume, and divide by number of layers to get heat content in J/m2
  223. e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * a_i(ji,jj,jl) * ht_i(ji,jj,jl) * r1_nlay_i
  224. END DO
  225. END DO
  226. CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy )
  227. CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
  228. CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
  229. CALL lbc_bdy_lnk( v_i(:,:,jl), 'T', 1., ib_bdy )
  230. CALL lbc_bdy_lnk( v_s(:,:,jl), 'T', 1., ib_bdy )
  231. CALL lbc_bdy_lnk( smv_i(:,:,jl), 'T', 1., ib_bdy )
  232. CALL lbc_bdy_lnk( sm_i(:,:,jl), 'T', 1., ib_bdy )
  233. CALL lbc_bdy_lnk( oa_i(:,:,jl), 'T', 1., ib_bdy )
  234. CALL lbc_bdy_lnk( t_su(:,:,jl), 'T', 1., ib_bdy )
  235. DO jk = 1, nlay_s
  236. CALL lbc_bdy_lnk(t_s(:,:,jk,jl), 'T', 1., ib_bdy )
  237. CALL lbc_bdy_lnk(e_s(:,:,jk,jl), 'T', 1., ib_bdy )
  238. END DO
  239. DO jk = 1, nlay_i
  240. CALL lbc_bdy_lnk(t_i(:,:,jk,jl), 'T', 1., ib_bdy )
  241. CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy )
  242. END DO
  243. END DO !jl
  244. #endif
  245. !
  246. IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_frs')
  247. !
  248. END SUBROUTINE bdy_ice_frs
  249. SUBROUTINE bdy_ice_lim_dyn( cd_type )
  250. !!------------------------------------------------------------------------------
  251. !! *** SUBROUTINE bdy_ice_lim_dyn ***
  252. !!
  253. !! ** Purpose : Apply dynamics boundary conditions for sea-ice in the cas of unstructured open boundaries.
  254. !! u_ice and v_ice are equal to the value of the adjacent grid point if this latter is not ice free
  255. !! if adjacent grid point is ice free, then u_ice and v_ice are equal to ocean velocities
  256. !!
  257. !! 2013-06 : C. Rousset
  258. !!------------------------------------------------------------------------------
  259. !!
  260. CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points
  261. INTEGER :: jb, jgrd ! dummy loop indices
  262. INTEGER :: ji, jj ! local scalar
  263. INTEGER :: ib_bdy ! Loop index
  264. REAL(wp) :: zmsk1, zmsk2, zflag
  265. !!------------------------------------------------------------------------------
  266. !
  267. IF( nn_timing == 1 ) CALL timing_start('bdy_ice_lim_dyn')
  268. !
  269. DO ib_bdy=1, nb_bdy
  270. !
  271. SELECT CASE( cn_ice_lim(ib_bdy) )
  272. CASE('none')
  273. CYCLE
  274. CASE('frs')
  275. IF( nn_ice_lim_dta(ib_bdy) == 0 ) CYCLE ! case ice boundaries = initial conditions
  276. ! do not change ice velocity (it is only computed by rheology)
  277. SELECT CASE ( cd_type )
  278. CASE ( 'U' )
  279. jgrd = 2 ! u velocity
  280. DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd)
  281. ji = idx_bdy(ib_bdy)%nbi(jb,jgrd)
  282. jj = idx_bdy(ib_bdy)%nbj(jb,jgrd)
  283. zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd)
  284. IF ( ABS( zflag ) == 1. ) THEN ! eastern and western boundaries
  285. ! one of the two zmsk is always 0 (because of zflag)
  286. zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice
  287. zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice
  288. ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
  289. u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
  290. & u_ice(ji-1,jj) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
  291. & u_oce(ji ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
  292. ELSE ! everywhere else
  293. !u_ice(ji,jj) = u_oce(ji,jj)
  294. u_ice(ji,jj) = 0._wp
  295. ENDIF
  296. ! mask ice velocities
  297. rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01_wp ) ) ! 0 if no ice
  298. u_ice(ji,jj) = rswitch * u_ice(ji,jj)
  299. ENDDO
  300. CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy )
  301. CASE ( 'V' )
  302. jgrd = 3 ! v velocity
  303. DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd)
  304. ji = idx_bdy(ib_bdy)%nbi(jb,jgrd)
  305. jj = idx_bdy(ib_bdy)%nbj(jb,jgrd)
  306. zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd)
  307. IF ( ABS( zflag ) == 1. ) THEN ! northern and southern boundaries
  308. ! one of the two zmsk is always 0 (because of zflag)
  309. zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice
  310. zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice
  311. ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
  312. v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
  313. & v_ice(ji,jj-1) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
  314. & v_oce(ji,jj ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
  315. ELSE ! everywhere else
  316. !v_ice(ji,jj) = v_oce(ji,jj)
  317. v_ice(ji,jj) = 0._wp
  318. ENDIF
  319. ! mask ice velocities
  320. rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01 ) ) ! 0 if no ice
  321. v_ice(ji,jj) = rswitch * v_ice(ji,jj)
  322. ENDDO
  323. CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy )
  324. END SELECT
  325. CASE DEFAULT
  326. CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' )
  327. END SELECT
  328. ENDDO
  329. IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn')
  330. END SUBROUTINE bdy_ice_lim_dyn
  331. #else
  332. !!---------------------------------------------------------------------------------
  333. !! Default option Empty module
  334. !!---------------------------------------------------------------------------------
  335. CONTAINS
  336. SUBROUTINE bdy_ice_lim( kt ) ! Empty routine
  337. WRITE(*,*) 'bdy_ice_lim: You should not have seen this print! error?', kt
  338. END SUBROUTINE bdy_ice_lim
  339. #endif
  340. !!=================================================================================
  341. END MODULE bdyice_lim