sbcice_lim.F90 32 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667
  1. MODULE sbcice_lim
  2. !!======================================================================
  3. !! *** MODULE sbcice_lim ***
  4. !! Surface module : update the ocean surface boundary condition over ice
  5. !! & covered area using LIM sea-ice model
  6. !! Sea-Ice model : LIM-3 Sea ice model time-stepping
  7. !!=====================================================================
  8. !! History : 2.0 ! 2006-12 (M. Vancoppenolle) Original code
  9. !! 3.0 ! 2008-02 (C. Talandier) Surface module from icestp.F90
  10. !! - ! 2008-04 (G. Madec) sltyle and lim_ctl routine
  11. !! 3.3 ! 2010-11 (G. Madec) ice-ocean stress always computed at each ocean time-step
  12. !! 3.4 ! 2011-01 (A Porter) dynamical allocation
  13. !! - ! 2012-10 (C. Rousset) add lim_diahsb
  14. !! 3.6 ! 2014-07 (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface
  15. !!----------------------------------------------------------------------
  16. #if defined key_lim3
  17. !!----------------------------------------------------------------------
  18. !! 'key_lim3' : LIM 3.0 sea-ice model
  19. !!----------------------------------------------------------------------
  20. !! sbc_ice_lim : sea-ice model time-stepping and update ocean sbc over ice-covered area
  21. !!----------------------------------------------------------------------
  22. USE oce ! ocean dynamics and tracers
  23. USE dom_oce ! ocean space and time domain
  24. USE ice ! LIM-3: ice variables
  25. USE thd_ice ! LIM-3: thermodynamical variables
  26. USE dom_ice ! LIM-3: ice domain
  27. USE sbc_oce ! Surface boundary condition: ocean fields
  28. USE sbc_ice ! Surface boundary condition: ice fields
  29. USE sbcblk_core ! Surface boundary condition: CORE bulk
  30. USE sbcblk_clio ! Surface boundary condition: CLIO bulk
  31. USE sbccpl ! Surface boundary condition: coupled interface
  32. USE albedo ! ocean & ice albedo
  33. USE phycst ! Define parameters for the routines
  34. USE eosbn2 ! equation of state
  35. USE limdyn ! Ice dynamics
  36. USE limtrp ! Ice transport
  37. USE limhdf ! Ice horizontal diffusion
  38. USE limthd ! Ice thermodynamics
  39. USE limitd_me ! Mechanics on ice thickness distribution
  40. USE limsbc ! sea surface boundary condition
  41. USE limdiahsb ! Ice budget diagnostics
  42. USE limwri ! Ice outputs
  43. USE limrst ! Ice restarts
  44. USE limupdate1 ! update of global variables
  45. USE limupdate2 ! update of global variables
  46. USE limvar ! Ice variables switch
  47. USE limmsh ! LIM mesh
  48. USE limistate ! LIM initial state
  49. USE limthd_sal ! LIM ice thermodynamics: salinity
  50. USE c1d ! 1D vertical configuration
  51. USE lbclnk ! lateral boundary condition - MPP link
  52. USE lib_mpp ! MPP library
  53. USE wrk_nemo ! work arrays
  54. USE timing ! Timing
  55. USE iom ! I/O manager library
  56. USE in_out_manager ! I/O manager
  57. USE prtctl ! Print control
  58. USE lib_fortran !
  59. USE limctl
  60. #if defined key_bdy
  61. USE bdyice_lim ! unstructured open boundary data (bdy_ice_lim routine)
  62. #endif
  63. IMPLICIT NONE
  64. PRIVATE
  65. PUBLIC sbc_ice_lim ! routine called by sbcmod.F90
  66. PUBLIC sbc_lim_init ! routine called by sbcmod.F90
  67. !! * Substitutions
  68. # include "domzgr_substitute.h90"
  69. # include "vectopt_loop_substitute.h90"
  70. !!----------------------------------------------------------------------
  71. !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011)
  72. !! $Id: sbcice_lim.F90 8158 2017-06-09 07:09:35Z vancop $
  73. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  74. !!----------------------------------------------------------------------
  75. CONTAINS
  76. !!======================================================================
  77. SUBROUTINE sbc_ice_lim( kt, kblk )
  78. !!---------------------------------------------------------------------
  79. !! *** ROUTINE sbc_ice_lim ***
  80. !!
  81. !! ** Purpose : update the ocean surface boundary condition via the
  82. !! Louvain la Neuve Sea Ice Model time stepping
  83. !!
  84. !! ** Method : ice model time stepping
  85. !! - call the ice dynamics routine
  86. !! - call the ice advection/diffusion routine
  87. !! - call the ice thermodynamics routine
  88. !! - call the routine that computes mass and
  89. !! heat fluxes at the ice/ocean interface
  90. !! - save the outputs
  91. !! - save the outputs for restart when necessary
  92. !!
  93. !! ** Action : - time evolution of the LIM sea-ice model
  94. !! - update all sbc variables below sea-ice:
  95. !! utau, vtau, taum, wndm, qns , qsr, emp , sfx
  96. !!---------------------------------------------------------------------
  97. INTEGER, INTENT(in) :: kt ! ocean time step
  98. INTEGER, INTENT(in) :: kblk ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED)
  99. !!
  100. INTEGER :: jl ! dummy loop index
  101. REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_os, zalb_cs ! ice albedo under overcast/clear sky
  102. REAL(wp), POINTER, DIMENSION(:,: ) :: zutau_ice, zvtau_ice
  103. !!----------------------------------------------------------------------
  104. IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim')
  105. !-----------------------!
  106. ! --- Ice time step --- !
  107. !-----------------------!
  108. IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN
  109. ! mean surface ocean current at ice velocity point (C-grid dynamics : U- & V-points as the ocean)
  110. u_oce(:,:) = ssu_m(:,:) * umask(:,:,1)
  111. v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1)
  112. ! masked sea surface freezing temperature [Kelvin] (set to rt0 over land)
  113. CALL eos_fzp( sss_m(:,:) , t_bo(:,:) )
  114. t_bo(:,:) = ( t_bo(:,:) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) )
  115. ! Mask sea ice surface temperature (set to rt0 over land)
  116. DO jl = 1, jpl
  117. t_su(:,:,jl) = t_su(:,:,jl) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) )
  118. END DO
  119. !
  120. !------------------------------------------------!
  121. ! --- Dynamical coupling with the atmosphere --- !
  122. !------------------------------------------------!
  123. ! It provides the following fields:
  124. ! utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2]
  125. !-----------------------------------------------------------------
  126. SELECT CASE( kblk )
  127. CASE( jp_clio ) ; CALL blk_ice_clio_tau ! CLIO bulk formulation
  128. CASE( jp_core ) ; CALL blk_ice_core_tau ! CORE bulk formulation
  129. CASE( jp_purecpl ) ; CALL sbc_cpl_ice_tau( utau_ice , vtau_ice ) ! Coupled formulation
  130. END SELECT
  131. IF( ln_mixcpl) THEN ! Case of a mixed Bulk/Coupled formulation
  132. CALL wrk_alloc( jpi,jpj , zutau_ice, zvtau_ice)
  133. CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice )
  134. utau_ice(:,:) = utau_ice(:,:) * xcplmask(:,:,0) + zutau_ice(:,:) * ( 1. - xcplmask(:,:,0) )
  135. vtau_ice(:,:) = vtau_ice(:,:) * xcplmask(:,:,0) + zvtau_ice(:,:) * ( 1. - xcplmask(:,:,0) )
  136. CALL wrk_dealloc( jpi,jpj , zutau_ice, zvtau_ice)
  137. ENDIF
  138. !-------------------------------------------------------!
  139. ! --- ice dynamics and transport (except in 1D case) ---!
  140. !-------------------------------------------------------!
  141. numit = numit + nn_fsbc ! Ice model time step
  142. !
  143. CALL sbc_lim_bef ! Store previous ice values
  144. CALL sbc_lim_diag0 ! set diag of mass, heat and salt fluxes to 0
  145. CALL lim_rst_opn( kt ) ! Open Ice restart file
  146. !
  147. IF( .NOT. lk_c1d ) THEN
  148. !
  149. CALL lim_dyn( kt ) ! Ice dynamics ( rheology/dynamics )
  150. !
  151. CALL lim_trp( kt ) ! Ice transport ( Advection/diffusion )
  152. !
  153. IF( nn_monocat /= 2 ) CALL lim_itd_me ! Mechanical redistribution ! (ridging/rafting)
  154. !
  155. #if defined key_bdy
  156. CALL bdy_ice_lim( kt ) ! bdy ice thermo
  157. IF( ln_icectl ) CALL lim_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' )
  158. #endif
  159. !
  160. CALL lim_update1( kt ) ! Corrections
  161. !
  162. ENDIF
  163. ! previous lead fraction and ice volume for flux calculations
  164. CALL sbc_lim_bef
  165. CALL lim_var_glo2eqv ! ht_i and ht_s for ice albedo calculation
  166. CALL lim_var_agg(1) ! at_i for coupling (via pfrld)
  167. pfrld(:,:) = 1._wp - at_i(:,:)
  168. phicif(:,:) = vt_i(:,:)
  169. !------------------------------------------------------!
  170. ! --- Thermodynamical coupling with the atmosphere --- !
  171. !------------------------------------------------------!
  172. ! It provides the following fields:
  173. ! qsr_ice , qns_ice : solar & non solar heat flux over ice (T-point) [W/m2]
  174. ! qla_ice : latent heat flux over ice (T-point) [W/m2]
  175. ! dqns_ice, dqla_ice : non solar & latent heat sensistivity (T-point) [W/m2]
  176. ! tprecip , sprecip : total & solid precipitation (T-point) [Kg/m2/s]
  177. ! fr1_i0 , fr2_i0 : 1sr & 2nd fraction of qsr penetration in ice [%]
  178. !----------------------------------------------------------------------------------------
  179. CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs )
  180. CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos
  181. SELECT CASE( kblk )
  182. CASE( jp_clio ) ! CLIO bulk formulation
  183. ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo
  184. ! (alb_ice) is computed within the bulk routine
  185. CALL blk_ice_clio_flx( t_su, zalb_cs, zalb_os, alb_ice )
  186. IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
  187. IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
  188. CASE( jp_core ) ! CORE bulk formulation
  189. ! albedo depends on cloud fraction because of non-linear spectral effects
  190. alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:)
  191. CALL blk_ice_core_flx( t_su, alb_ice )
  192. IF( ln_mixcpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
  193. IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
  194. CASE ( jp_purecpl )
  195. ! albedo depends on cloud fraction because of non-linear spectral effects
  196. alb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:)
  197. CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=alb_ice, psst=sst_m, pist=t_su )
  198. IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx )
  199. END SELECT
  200. CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs )
  201. !----------------------------!
  202. ! --- ice thermodynamics --- !
  203. !----------------------------!
  204. CALL lim_thd( kt ) ! Ice thermodynamics
  205. !
  206. CALL lim_update2( kt ) ! Corrections
  207. !
  208. CALL lim_sbc_flx( kt ) ! Update surface ocean mass, heat and salt fluxes
  209. !
  210. IF(ln_limdiaout) CALL lim_diahsb( kt ) ! Diagnostics and outputs
  211. !
  212. CALL lim_wri( 1 ) ! Ice outputs
  213. !
  214. IF( kt == nit000 .AND. ln_rstart ) &
  215. & CALL iom_close( numrir ) ! close input ice restart file
  216. !
  217. IF( lrst_ice ) CALL lim_rst_write( kt ) ! Ice restart file
  218. !
  219. IF( ln_icectl ) CALL lim_ctl( kt ) ! alerts in case of model crash
  220. !
  221. ENDIF ! End sea-ice time step only
  222. !-------------------------!
  223. ! --- Ocean time step --- !
  224. !-------------------------!
  225. ! Update surface ocean stresses (only in ice-dynamic case) otherwise the atm.-ocean stresses are used everywhere
  226. IF( ln_limdyn ) CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) ) ! using before instantaneous surf. currents
  227. !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!!
  228. !
  229. IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim')
  230. !
  231. END SUBROUTINE sbc_ice_lim
  232. SUBROUTINE sbc_lim_init
  233. !!----------------------------------------------------------------------
  234. !! *** ROUTINE sbc_lim_init ***
  235. !!
  236. !! ** purpose : Allocate all the dynamic arrays of the LIM-3 modules
  237. !!----------------------------------------------------------------------
  238. INTEGER :: ierr
  239. INTEGER :: ji, jj
  240. !!----------------------------------------------------------------------
  241. IF(lwp) WRITE(numout,*)
  242. IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition'
  243. IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ via Louvain la Neuve Ice Model (LIM-3) time stepping'
  244. !
  245. ! Open the reference and configuration namelist files and namelist output file
  246. CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp )
  247. CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg', 'OLD', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp )
  248. IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 )
  249. CALL ice_run ! set some ice run parameters
  250. !
  251. ! ! Allocate the ice arrays
  252. ierr = ice_alloc () ! ice variables
  253. ierr = ierr + dom_ice_alloc () ! domain
  254. ierr = ierr + sbc_ice_alloc () ! surface forcing
  255. ierr = ierr + thd_ice_alloc () ! thermodynamics
  256. ierr = ierr + lim_itd_me_alloc () ! ice thickness distribution - mechanics
  257. !
  258. IF( lk_mpp ) CALL mpp_sum( ierr )
  259. IF( ierr /= 0 ) CALL ctl_stop('STOP', 'sbc_lim_init : unable to allocate ice arrays')
  260. !
  261. ! ! adequation jpk versus ice/snow layers/categories
  262. IF( jpl > jpk .OR. (nlay_i+1) > jpk .OR. nlay_s > jpk ) &
  263. & CALL ctl_stop( 'STOP', &
  264. & 'sbc_lim_init: the 3rd dimension of workspace arrays is too small.', &
  265. & 'use more ocean levels or less ice/snow layers/categories.' )
  266. !
  267. CALL lim_itd_init ! ice thickness distribution initialization
  268. !
  269. CALL lim_hdf_init ! set ice horizontal diffusion computation parameters
  270. !
  271. CALL lim_thd_init ! set ice thermodynics parameters
  272. !
  273. CALL lim_thd_sal_init ! set ice salinity parameters
  274. !
  275. CALL lim_msh ! ice mesh initialization
  276. !
  277. CALL lim_itd_me_init ! ice thickness distribution initialization for mecanical deformation
  278. ! ! Initial sea-ice state
  279. IF( .NOT. ln_rstart ) THEN ! start from rest: sea-ice deduced from sst
  280. numit = 0
  281. numit = nit000 - 1
  282. CALL lim_istate
  283. ELSE ! start from a restart file
  284. CALL lim_rst_read
  285. numit = nit000 - 1
  286. ENDIF
  287. CALL lim_var_agg(2)
  288. CALL lim_var_glo2eqv
  289. !
  290. CALL lim_sbc_init ! ice surface boundary condition
  291. !
  292. IF( ln_limdiaout) CALL lim_diahsb_init ! initialization for diags
  293. !
  294. fr_i(:,:) = at_i(:,:) ! initialisation of sea-ice fraction
  295. tn_ice(:,:,:) = t_su(:,:,:) ! initialisation of surface temp for coupled simu
  296. !
  297. DO jj = 1, jpj
  298. DO ji = 1, jpi
  299. IF( gphit(ji,jj) > 0._wp ) THEN ; rn_amax_2d(ji,jj) = rn_amax_n ! NH
  300. ELSE ; rn_amax_2d(ji,jj) = rn_amax_s ! SH
  301. ENDIF
  302. ENDDO
  303. ENDDO
  304. !
  305. nstart = numit + nn_fsbc
  306. nitrun = nitend - nit000 + 1
  307. nlast = numit + nitrun
  308. !
  309. IF( nstock == 0 ) nstock = nlast + 1
  310. !
  311. END SUBROUTINE sbc_lim_init
  312. SUBROUTINE ice_run
  313. !!-------------------------------------------------------------------
  314. !! *** ROUTINE ice_run ***
  315. !!
  316. !! ** Purpose : Definition some run parameter for ice model
  317. !!
  318. !! ** Method : Read the namicerun namelist and check the parameter
  319. !! values called at the first timestep (nit000)
  320. !!
  321. !! ** input : Namelist namicerun
  322. !!-------------------------------------------------------------------
  323. INTEGER :: ios ! Local integer output status for namelist read
  324. NAMELIST/namicerun/ jpl, nlay_i, nlay_s, cn_icerst_in, cn_icerst_indir, cn_icerst_out, cn_icerst_outdir, &
  325. & ln_limdyn, rn_amax_n, rn_amax_s, ln_limdiahsb, ln_limdiaout, ln_icectl, iiceprt, jiceprt
  326. !!-------------------------------------------------------------------
  327. !
  328. REWIND( numnam_ice_ref ) ! Namelist namicerun in reference namelist : Parameters for ice
  329. READ ( numnam_ice_ref, namicerun, IOSTAT = ios, ERR = 901)
  330. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in reference namelist', lwp )
  331. REWIND( numnam_ice_cfg ) ! Namelist namicerun in configuration namelist : Parameters for ice
  332. READ ( numnam_ice_cfg, namicerun, IOSTAT = ios, ERR = 902 )
  333. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in configuration namelist', lwp )
  334. IF(lwm) WRITE ( numoni, namicerun )
  335. !
  336. !
  337. IF(lwp) THEN ! control print
  338. WRITE(numout,*)
  339. WRITE(numout,*) 'ice_run : ice share parameters for dynamics/advection/thermo of sea-ice'
  340. WRITE(numout,*) ' ~~~~~~'
  341. WRITE(numout,*) ' number of ice categories = ', jpl
  342. WRITE(numout,*) ' number of ice layers = ', nlay_i
  343. WRITE(numout,*) ' number of snow layers = ', nlay_s
  344. WRITE(numout,*) ' switch for ice dynamics (1) or not (0) ln_limdyn = ', ln_limdyn
  345. WRITE(numout,*) ' maximum ice concentration for NH = ', rn_amax_n
  346. WRITE(numout,*) ' maximum ice concentration for SH = ', rn_amax_s
  347. WRITE(numout,*) ' Diagnose heat/salt budget or not ln_limdiahsb = ', ln_limdiahsb
  348. WRITE(numout,*) ' Output heat/salt budget or not ln_limdiaout = ', ln_limdiaout
  349. WRITE(numout,*) ' control prints in ocean.out for (i,j)=(iiceprt,jiceprt) = ', ln_icectl
  350. WRITE(numout,*) ' i-index for control prints (ln_icectl=true) = ', iiceprt
  351. WRITE(numout,*) ' j-index for control prints (ln_icectl=true) = ', jiceprt
  352. ENDIF
  353. !
  354. ! sea-ice timestep and inverse
  355. rdt_ice = nn_fsbc * rdttra(1)
  356. r1_rdtice = 1._wp / rdt_ice
  357. ! inverse of nlay_i and nlay_s
  358. r1_nlay_i = 1._wp / REAL( nlay_i, wp )
  359. r1_nlay_s = 1._wp / REAL( nlay_s, wp )
  360. !
  361. #if defined key_bdy
  362. IF( lwp .AND. ln_limdiahsb ) CALL ctl_warn('online conservation check activated but it does not work with BDY')
  363. #endif
  364. !
  365. END SUBROUTINE ice_run
  366. SUBROUTINE lim_itd_init
  367. !!------------------------------------------------------------------
  368. !! *** ROUTINE lim_itd_init ***
  369. !!
  370. !! ** Purpose : Initializes the ice thickness distribution
  371. !! ** Method : ...
  372. !! ** input : Namelist namiceitd
  373. !!-------------------------------------------------------------------
  374. INTEGER :: ios ! Local integer output status for namelist read
  375. NAMELIST/namiceitd/ nn_catbnd, rn_himean
  376. !
  377. INTEGER :: jl ! dummy loop index
  378. REAL(wp) :: zc1, zc2, zc3, zx1 ! local scalars
  379. REAL(wp) :: zhmax, znum, zden, zalpha !
  380. !!------------------------------------------------------------------
  381. !
  382. REWIND( numnam_ice_ref ) ! Namelist namiceitd in reference namelist : Parameters for ice
  383. READ ( numnam_ice_ref, namiceitd, IOSTAT = ios, ERR = 903)
  384. 903 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in reference namelist', lwp )
  385. REWIND( numnam_ice_cfg ) ! Namelist namiceitd in configuration namelist : Parameters for ice
  386. READ ( numnam_ice_cfg, namiceitd, IOSTAT = ios, ERR = 904 )
  387. 904 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in configuration namelist', lwp )
  388. IF(lwm) WRITE ( numoni, namiceitd )
  389. !
  390. !
  391. IF(lwp) THEN ! control print
  392. WRITE(numout,*)
  393. WRITE(numout,*) 'ice_itd : ice cat distribution'
  394. WRITE(numout,*) ' ~~~~~~'
  395. WRITE(numout,*) ' shape of ice categories distribution nn_catbnd = ', nn_catbnd
  396. WRITE(numout,*) ' mean ice thickness in the domain (only active if nn_catbnd=2) rn_himean = ', rn_himean
  397. ENDIF
  398. !----------------------------------
  399. !- Thickness categories boundaries
  400. !----------------------------------
  401. IF(lwp) WRITE(numout,*)
  402. IF(lwp) WRITE(numout,*) 'lim_itd_init : Initialization of ice cat distribution '
  403. IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
  404. hi_max(:) = 0._wp
  405. SELECT CASE ( nn_catbnd )
  406. !----------------------
  407. CASE (1) ! tanh function (CICE)
  408. !----------------------
  409. zc1 = 3._wp / REAL( jpl, wp )
  410. zc2 = 10._wp * zc1
  411. zc3 = 3._wp
  412. DO jl = 1, jpl
  413. zx1 = REAL( jl-1, wp ) / REAL( jpl, wp )
  414. hi_max(jl) = hi_max(jl-1) + zc1 + zc2 * (1._wp + TANH( zc3 * (zx1 - 1._wp ) ) )
  415. END DO
  416. !----------------------
  417. CASE (2) ! h^(-alpha) function
  418. !----------------------
  419. zalpha = 0.05 ! exponent of the transform function
  420. zhmax = 3.*rn_himean
  421. DO jl = 1, jpl
  422. znum = jpl * ( zhmax+1 )**zalpha
  423. zden = ( jpl - jl ) * ( zhmax+1 )**zalpha + jl
  424. hi_max(jl) = ( znum / zden )**(1./zalpha) - 1
  425. END DO
  426. END SELECT
  427. DO jl = 1, jpl
  428. hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp
  429. END DO
  430. ! Set hi_max(jpl) to a big value to ensure that all ice is thinner than hi_max(jpl)
  431. hi_max(jpl) = 99._wp
  432. IF(lwp) WRITE(numout,*) ' Thickness category boundaries '
  433. IF(lwp) WRITE(numout,*) ' hi_max ', hi_max(0:jpl)
  434. !
  435. END SUBROUTINE lim_itd_init
  436. SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_limflx )
  437. !!---------------------------------------------------------------------
  438. !! *** ROUTINE ice_lim_flx ***
  439. !!
  440. !! ** Purpose : update the ice surface boundary condition by averaging and / or
  441. !! redistributing fluxes on ice categories
  442. !!
  443. !! ** Method : average then redistribute
  444. !!
  445. !! ** Action :
  446. !!---------------------------------------------------------------------
  447. INTEGER , INTENT(in ) :: k_limflx ! =-1 do nothing; =0 average ;
  448. ! =1 average and redistribute ; =2 redistribute
  449. REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: ptn_ice ! ice surface temperature
  450. REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb_ice ! ice albedo
  451. REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqns_ice ! non solar flux
  452. REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqsr_ice ! net solar flux
  453. REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdqn_ice ! non solar flux sensitivity
  454. REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pevap_ice ! sublimation
  455. REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdevap_ice ! sublimation sensitivity
  456. !
  457. INTEGER :: jl ! dummy loop index
  458. !
  459. REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m ! Mean albedo over all categories
  460. REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m ! Mean temperature over all categories
  461. !
  462. REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m ! Mean solar heat flux over all categories
  463. REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m ! Mean non solar heat flux over all categories
  464. REAL(wp), POINTER, DIMENSION(:,:) :: z_evap_m ! Mean sublimation over all categories
  465. REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m ! Mean d(qns)/dT over all categories
  466. REAL(wp), POINTER, DIMENSION(:,:) :: z_devap_m ! Mean d(evap)/dT over all categories
  467. !!----------------------------------------------------------------------
  468. IF( nn_timing == 1 ) CALL timing_start('ice_lim_flx')
  469. !
  470. !
  471. SELECT CASE( k_limflx ) !== averaged on all ice categories ==!
  472. CASE( 0 , 1 )
  473. CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m)
  474. !
  475. z_qns_m (:,:) = fice_ice_ave ( pqns_ice (:,:,:) )
  476. z_qsr_m (:,:) = fice_ice_ave ( pqsr_ice (:,:,:) )
  477. z_dqn_m (:,:) = fice_ice_ave ( pdqn_ice (:,:,:) )
  478. z_evap_m (:,:) = fice_ice_ave ( pevap_ice (:,:,:) )
  479. z_devap_m(:,:) = fice_ice_ave ( pdevap_ice (:,:,:) )
  480. DO jl = 1, jpl
  481. pdqn_ice (:,:,jl) = z_dqn_m(:,:)
  482. pdevap_ice(:,:,jl) = z_devap_m(:,:)
  483. END DO
  484. !
  485. DO jl = 1, jpl
  486. pqns_ice (:,:,jl) = z_qns_m(:,:)
  487. pqsr_ice (:,:,jl) = z_qsr_m(:,:)
  488. pevap_ice(:,:,jl) = z_evap_m(:,:)
  489. END DO
  490. !
  491. CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m)
  492. END SELECT
  493. SELECT CASE( k_limflx ) !== redistribution on all ice categories ==!
  494. CASE( 1 , 2 )
  495. CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m )
  496. !
  497. zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) )
  498. ztem_m(:,:) = fice_ice_ave ( ptn_ice (:,:,:) )
  499. DO jl = 1, jpl
  500. pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
  501. pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
  502. pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) )
  503. END DO
  504. !
  505. CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m )
  506. END SELECT
  507. !
  508. IF( nn_timing == 1 ) CALL timing_stop('ice_lim_flx')
  509. !
  510. END SUBROUTINE ice_lim_flx
  511. SUBROUTINE sbc_lim_bef
  512. !!----------------------------------------------------------------------
  513. !! *** ROUTINE sbc_lim_bef ***
  514. !!
  515. !! ** purpose : store ice variables at "before" time step
  516. !!----------------------------------------------------------------------
  517. a_i_b (:,:,:) = a_i (:,:,:) ! ice area
  518. e_i_b (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy
  519. v_i_b (:,:,:) = v_i (:,:,:) ! ice volume
  520. v_s_b (:,:,:) = v_s (:,:,:) ! snow volume
  521. e_s_b (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy
  522. smv_i_b(:,:,:) = smv_i(:,:,:) ! salt content
  523. oa_i_b (:,:,:) = oa_i (:,:,:) ! areal age content
  524. u_ice_b(:,:) = u_ice(:,:)
  525. v_ice_b(:,:) = v_ice(:,:)
  526. END SUBROUTINE sbc_lim_bef
  527. SUBROUTINE sbc_lim_diag0
  528. !!----------------------------------------------------------------------
  529. !! *** ROUTINE sbc_lim_diag0 ***
  530. !!
  531. !! ** purpose : set ice-ocean and ice-atm. fluxes to zeros at the beggining
  532. !! of the time step
  533. !!----------------------------------------------------------------------
  534. sfx (:,:) = 0._wp ;
  535. sfx_bri(:,:) = 0._wp ;
  536. sfx_sni(:,:) = 0._wp ; sfx_opw(:,:) = 0._wp
  537. sfx_bog(:,:) = 0._wp ; sfx_dyn(:,:) = 0._wp
  538. sfx_bom(:,:) = 0._wp ; sfx_sum(:,:) = 0._wp
  539. sfx_res(:,:) = 0._wp ; sfx_sub(:,:) = 0._wp
  540. wfx_snw(:,:) = 0._wp ; wfx_ice(:,:) = 0._wp
  541. wfx_sni(:,:) = 0._wp ; wfx_opw(:,:) = 0._wp
  542. wfx_bog(:,:) = 0._wp ; wfx_dyn(:,:) = 0._wp
  543. wfx_bom(:,:) = 0._wp ; wfx_sum(:,:) = 0._wp
  544. wfx_res(:,:) = 0._wp ; wfx_sub(:,:) = 0._wp
  545. wfx_spr(:,:) = 0._wp ;
  546. wfx_snw_dyn(:,:) = 0._wp ; wfx_snw_sum(:,:) = 0._wp
  547. wfx_snw_sub(:,:) = 0._wp ; wfx_ice_sub(:,:) = 0._wp
  548. hfx_thd(:,:) = 0._wp ;
  549. hfx_snw(:,:) = 0._wp ; hfx_opw(:,:) = 0._wp
  550. hfx_bog(:,:) = 0._wp ; hfx_dyn(:,:) = 0._wp
  551. hfx_bom(:,:) = 0._wp ; hfx_sum(:,:) = 0._wp
  552. hfx_res(:,:) = 0._wp ; hfx_sub(:,:) = 0._wp
  553. hfx_spr(:,:) = 0._wp ; hfx_dif(:,:) = 0._wp
  554. hfx_err(:,:) = 0._wp ; hfx_err_rem(:,:) = 0._wp
  555. hfx_err_dif(:,:) = 0._wp
  556. wfx_err_sub(:,:) = 0._wp
  557. afx_tot(:,:) = 0._wp ;
  558. afx_dyn(:,:) = 0._wp ; afx_thd(:,:) = 0._wp
  559. diag_heat(:,:) = 0._wp ; diag_smvi(:,:) = 0._wp ;
  560. diag_vice(:,:) = 0._wp ; diag_vsnw(:,:) = 0._wp ;
  561. ! SIMIP diagnostics
  562. diag_dms_dyn(:,:) = 0._wp ; diag_dmi_dyn(:,:) = 0._wp
  563. diag_fc_bo(:,:) = 0._wp ; diag_fc_su(:,:) = 0._wp
  564. END SUBROUTINE sbc_lim_diag0
  565. FUNCTION fice_cell_ave ( ptab )
  566. !!--------------------------------------------------------------------------
  567. !! * Compute average over categories, for grid cell (ice covered and free ocean)
  568. !!--------------------------------------------------------------------------
  569. REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave
  570. REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab
  571. INTEGER :: jl ! Dummy loop index
  572. fice_cell_ave (:,:) = 0.0_wp
  573. DO jl = 1, jpl
  574. fice_cell_ave (:,:) = fice_cell_ave (:,:) + a_i (:,:,jl) * ptab (:,:,jl)
  575. END DO
  576. END FUNCTION fice_cell_ave
  577. FUNCTION fice_ice_ave ( ptab )
  578. !!--------------------------------------------------------------------------
  579. !! * Compute average over categories, for ice covered part of grid cell
  580. !!--------------------------------------------------------------------------
  581. REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave
  582. REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab
  583. fice_ice_ave (:,:) = 0.0_wp
  584. WHERE ( at_i (:,:) > 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:)
  585. END FUNCTION fice_ice_ave
  586. #else
  587. !!----------------------------------------------------------------------
  588. !! Default option Dummy module NO LIM 3.0 sea-ice model
  589. !!----------------------------------------------------------------------
  590. CONTAINS
  591. SUBROUTINE sbc_ice_lim ( kt, kblk ) ! Dummy routine
  592. INTEGER, INTENT(in) :: kt, kblk
  593. WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk
  594. END SUBROUTINE sbc_ice_lim
  595. SUBROUTINE sbc_lim_init ! Dummy routine
  596. END SUBROUTINE sbc_lim_init
  597. #endif
  598. !!======================================================================
  599. END MODULE sbcice_lim