sbcblk_clio.F90 59 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120
  1. MODULE sbcblk_clio
  2. !!======================================================================
  3. !! *** MODULE sbcblk_clio ***
  4. !! Ocean forcing: bulk thermohaline forcing of the ocean (or ice)
  5. !!=====================================================================
  6. !! History : OPA ! 1997-06 (Louvain-La-Neuve) Original code
  7. !! ! 2001-04 (C. Ethe) add flx_blk_declin
  8. !! NEMO 2.0 ! 2002-08 (C. Ethe, G. Madec) F90: Free form and module
  9. !! 3.0 ! 2008-03 (C. Talandier, G. Madec) surface module + LIM3
  10. !! 3.2 ! 2009-04 (B. Lemaire) Introduce iom_put
  11. !!----------------------------------------------------------------------
  12. !!----------------------------------------------------------------------
  13. !! sbc_blk_clio : CLIO bulk formulation: read and update required input fields
  14. !! blk_clio_oce : ocean CLIO bulk formulea: compute momentum, heat and freswater fluxes for the ocean
  15. !! blk_ice_clio : ice CLIO bulk formulea: compute momentum, heat and freswater fluxes for the sea-ice
  16. !! blk_clio_qsr_oce : shortwave radiation for ocean computed from the cloud cover
  17. !! blk_clio_qsr_ice : shortwave radiation for ice computed from the cloud cover
  18. !! flx_blk_declin : solar declination
  19. !!----------------------------------------------------------------------
  20. USE oce ! ocean dynamics and tracers
  21. USE dom_oce ! ocean space and time domain
  22. USE phycst ! physical constants
  23. USE fldread ! read input fields
  24. USE sbc_oce ! Surface boundary condition: ocean fields
  25. USE iom ! I/O manager library
  26. USE in_out_manager ! I/O manager
  27. USE lib_mpp ! distribued memory computing library
  28. USE wrk_nemo ! work arrays
  29. USE timing ! Timing
  30. USE lbclnk ! ocean lateral boundary conditions (or mpp link)
  31. USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
  32. USE albedo
  33. USE prtctl ! Print control
  34. #if defined key_lim3
  35. USE ice
  36. USE sbc_ice ! Surface boundary condition: ice fields
  37. USE limthd_dh ! for CALL lim_thd_snwblow
  38. #elif defined key_lim2
  39. USE ice_2
  40. USE sbc_ice ! Surface boundary condition: ice fields
  41. USE par_ice_2 ! Surface boundary condition: ice fields
  42. #endif
  43. IMPLICIT NONE
  44. PRIVATE
  45. PUBLIC sbc_blk_clio ! routine called by sbcmod.F90
  46. #if defined key_lim2 || defined key_lim3
  47. PUBLIC blk_ice_clio_tau ! routine called by sbcice_lim.F90
  48. PUBLIC blk_ice_clio_flx ! routine called by sbcice_lim.F90
  49. #endif
  50. INTEGER , PARAMETER :: jpfld = 7 ! maximum number of files to read
  51. INTEGER , PARAMETER :: jp_utau = 1 ! index of wind stress (i-component) (N/m2) at U-point
  52. INTEGER , PARAMETER :: jp_vtau = 2 ! index of wind stress (j-component) (N/m2) at V-point
  53. INTEGER , PARAMETER :: jp_wndm = 3 ! index of 10m wind module (m/s) at T-point
  54. INTEGER , PARAMETER :: jp_humi = 4 ! index of specific humidity ( % )
  55. INTEGER , PARAMETER :: jp_ccov = 5 ! index of cloud cover ( % )
  56. INTEGER , PARAMETER :: jp_tair = 6 ! index of 10m air temperature (Kelvin)
  57. INTEGER , PARAMETER :: jp_prec = 7 ! index of total precipitation (rain+snow) (Kg/m2/s)
  58. TYPE(FLD),ALLOCATABLE,DIMENSION(:) :: sf ! structure of input fields (file informations, fields read)
  59. INTEGER, PARAMETER :: jpintsr = 24 ! number of time step between sunrise and sunset
  60. ! ! uses for heat flux computation
  61. LOGICAL :: lbulk_init = .TRUE. ! flag, bulk initialization done or not)
  62. REAL(wp) :: cai = 1.40e-3 ! best estimate of atm drag in order to get correct FS export in ORCA2-LIM
  63. REAL(wp) :: cao = 1.00e-3 ! chosen by default ==> should depends on many things... !!gmto be updated
  64. REAL(wp) :: rdtbs2 !:
  65. REAL(wp), DIMENSION(19) :: budyko ! BUDYKO's coefficient (cloudiness effect on LW radiation)
  66. DATA budyko / 1.00, 0.98, 0.95, 0.92, 0.89, 0.86, 0.83, 0.80, 0.78, 0.75, &
  67. & 0.72, 0.69, 0.67, 0.64, 0.61, 0.58, 0.56, 0.53, 0.50 /
  68. REAL(wp), DIMENSION(20) :: tauco ! cloud optical depth coefficient
  69. DATA tauco / 6.6, 6.6, 7.0, 7.2, 7.1, 6.8, 6.5, 6.6, 7.1, 7.6, &
  70. & 6.6, 6.1, 5.6, 5.5, 5.8, 5.8, 5.6, 5.6, 5.6, 5.6 /
  71. !!
  72. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: sbudyko ! cloudiness effect on LW radiation
  73. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: stauc ! cloud optical depth
  74. REAL(wp) :: eps20 = 1.e-20 ! constant values
  75. !! * Substitutions
  76. # include "vectopt_loop_substitute.h90"
  77. !!----------------------------------------------------------------------
  78. !! NEMO/OPA 4.0 , NEMO Consortium (2011)
  79. !! $Id: sbcblk_clio.F90 6399 2016-03-22 17:17:23Z clem $
  80. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  81. !!----------------------------------------------------------------------
  82. CONTAINS
  83. SUBROUTINE sbc_blk_clio( kt )
  84. !!---------------------------------------------------------------------
  85. !! *** ROUTINE sbc_blk_clio ***
  86. !!
  87. !! ** Purpose : provide at each time step the surface ocean fluxes
  88. !! (momentum, heat, freshwater and runoff)
  89. !!
  90. !! ** Method : (1) READ each fluxes in NetCDF files:
  91. !! the i-component of the stress (N/m2)
  92. !! the j-component of the stress (N/m2)
  93. !! the 10m wind speed module (m/s)
  94. !! the 10m air temperature (Kelvin)
  95. !! the 10m specific humidity (%)
  96. !! the cloud cover (%)
  97. !! the total precipitation (rain+snow) (Kg/m2/s)
  98. !! (2) CALL blk_oce_clio
  99. !!
  100. !! C A U T I O N : never mask the surface stress fields
  101. !! the stress is assumed to be in the (i,j) mesh referential
  102. !!
  103. !! ** Action : defined at each time-step at the air-sea interface
  104. !! - utau, vtau i- and j-component of the wind stress
  105. !! - taum wind stress module at T-point
  106. !! - wndm 10m wind module at T-point over free ocean or leads in presence of sea-ice
  107. !! - qns non-solar heat flux including latent heat of solid
  108. !! precip. melting and emp heat content
  109. !! - qsr solar heat flux
  110. !! - emp upward mass flux (evap. - precip)
  111. !! - sfx salt flux; set to zero at nit000 but possibly non-zero
  112. !! if ice is present (computed in limsbc(_2).F90)
  113. !!----------------------------------------------------------------------
  114. INTEGER, INTENT( in ) :: kt ! ocean time step
  115. !!
  116. INTEGER :: ifpr, jfpr ! dummy indices
  117. INTEGER :: ierr0, ierr1, ierr2, ierr3 ! return error code
  118. INTEGER :: ios ! Local integer output status for namelist read
  119. !!
  120. CHARACTER(len=100) :: cn_dir ! Root directory for location of CLIO files
  121. TYPE(FLD_N), DIMENSION(jpfld) :: slf_i ! array of namelist informations on the fields to read
  122. TYPE(FLD_N) :: sn_utau, sn_vtau, sn_wndm, sn_tair ! informations about the fields to be read
  123. TYPE(FLD_N) :: sn_humi, sn_ccov, sn_prec ! " "
  124. !!
  125. NAMELIST/namsbc_clio/ cn_dir, sn_utau, sn_vtau, sn_wndm, sn_humi, &
  126. & sn_ccov, sn_tair, sn_prec
  127. !!---------------------------------------------------------------------
  128. ! ! ====================== !
  129. IF( kt == nit000 ) THEN ! First call kt=nit000 !
  130. ! ! ====================== !
  131. REWIND( numnam_ref ) ! Namelist namsbc_clio in reference namelist : CLIO files
  132. READ ( numnam_ref, namsbc_clio, IOSTAT = ios, ERR = 901)
  133. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_clio in reference namelist', lwp )
  134. REWIND( numnam_cfg ) ! Namelist namsbc_clio in configuration namelist : CLIO files
  135. READ ( numnam_cfg, namsbc_clio, IOSTAT = ios, ERR = 902 )
  136. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_clio in configuration namelist', lwp )
  137. IF(lwm) WRITE ( numond, namsbc_clio )
  138. ! store namelist information in an array
  139. slf_i(jp_utau) = sn_utau ; slf_i(jp_vtau) = sn_vtau ; slf_i(jp_wndm) = sn_wndm
  140. slf_i(jp_tair) = sn_tair ; slf_i(jp_humi) = sn_humi
  141. slf_i(jp_ccov) = sn_ccov ; slf_i(jp_prec) = sn_prec
  142. ! set sf structure
  143. ALLOCATE( sf(jpfld), STAT=ierr0 )
  144. IF( ierr0 > 0 ) CALL ctl_stop( 'STOP', 'sbc_blk_clio: unable to allocate sf structure' )
  145. DO ifpr= 1, jpfld
  146. ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) , STAT=ierr1)
  147. IF( slf_i(ifpr)%ln_tint ) ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) , STAT=ierr2 )
  148. END DO
  149. IF( ierr1+ierr2 > 0 ) CALL ctl_stop( 'STOP', 'sbc_blk_clio: unable to allocate sf array structure' )
  150. ! fill sf with slf_i and control print
  151. CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_clio', 'flux formulation for ocean surface boundary condition', 'namsbc_clio' )
  152. ! allocate sbcblk clio arrays
  153. ALLOCATE( sbudyko(jpi,jpj) , stauc(jpi,jpj), STAT=ierr3 )
  154. IF( ierr3 > 0 ) CALL ctl_stop( 'STOP', 'sbc_blk_clio: unable to allocate arrays' )
  155. !
  156. sfx(:,:) = 0._wp ! salt flux; zero unless ice is present (computed in limsbc(_2).F90)
  157. !
  158. ENDIF
  159. ! ! ====================== !
  160. ! ! At each time-step !
  161. ! ! ====================== !
  162. !
  163. CALL fld_read( kt, nn_fsbc, sf ) ! input fields provided at the current time-step
  164. !
  165. IF( MOD( kt - 1, nn_fsbc ) == 0 ) CALL blk_oce_clio( sf, sst_m )
  166. !
  167. END SUBROUTINE sbc_blk_clio
  168. SUBROUTINE blk_oce_clio( sf, pst )
  169. !!---------------------------------------------------------------------------
  170. !! *** ROUTINE blk_oce_clio ***
  171. !!
  172. !! ** Purpose : Compute momentum, heat and freshwater fluxes at ocean surface
  173. !! using CLIO bulk formulea
  174. !!
  175. !! ** Method : The flux of heat at the ocean surfaces are derived
  176. !! from semi-empirical ( or bulk ) formulae which relate the flux to
  177. !! the properties of the surface and of the lower atmosphere. Here, we
  178. !! follow the work of Oberhuber, 1988
  179. !! - momentum flux (stresses) directly read in files at U- and V-points
  180. !! - compute ocean/ice albedos (call albedo_oce/albedo_ice)
  181. !! - compute shortwave radiation for ocean (call blk_clio_qsr_oce)
  182. !! - compute long-wave radiation for the ocean
  183. !! - compute the turbulent heat fluxes over the ocean
  184. !! - deduce the evaporation over the ocean
  185. !! ** Action : Fluxes over the ocean:
  186. !! - utau, vtau i- and j-component of the wind stress
  187. !! - taum wind stress module at T-point
  188. !! - wndm 10m wind module at T-point over free ocean or leads in presence of sea-ice
  189. !! - qns non-solar heat flux including latent heat of solid
  190. !! precip. melting and emp heat content
  191. !! - qsr solar heat flux
  192. !! - emp suface mass flux (evap.-precip.)
  193. !! ** Nota : sf has to be a dummy argument for AGRIF on NEC
  194. !!----------------------------------------------------------------------
  195. TYPE(fld), INTENT(in), DIMENSION(:) :: sf ! input data
  196. REAL(wp) , INTENT(in), DIMENSION(jpi,jpj) :: pst ! surface temperature [Celcius]
  197. !!
  198. INTEGER :: ji, jj ! dummy loop indices
  199. !!
  200. REAL(wp) :: zrhova, zcsho, zcleo, zcldeff ! temporary scalars
  201. REAL(wp) :: zqsato, zdteta, zdeltaq, ztvmoy, zobouks ! - -
  202. REAL(wp) :: zpsims, zpsihs, zpsils, zobouku, zxins, zpsimu ! - -
  203. REAL(wp) :: zpsihu, zpsilu, zstab,zpsim, zpsih, zpsil ! - -
  204. REAL(wp) :: zvatmg, zcmn, zchn, zcln, zcmcmn, zdenum ! - -
  205. REAL(wp) :: zdtetar, ztvmoyr, zlxins, zchcm, zclcm ! - -
  206. REAL(wp) :: zmt1, zmt2, zmt3, ztatm3, ztamr, ztaevbk ! - -
  207. REAL(wp) :: zsst, ztatm, zcco1, zpatm, zcmax, zrmax ! - -
  208. REAL(wp) :: zrhoa, zev, zes, zeso, zqatm, zevsqr ! - -
  209. REAL(wp) :: ztx2, zty2, zcevap, zcprec ! - -
  210. REAL(wp) :: ztau1 ! TECLIM change
  211. REAL(wp), DIMENSION(jpi,jpj) :: ztau2 ! TECLIM change
  212. REAL(wp), POINTER, DIMENSION(:,:) :: zqlw ! long-wave heat flux over ocean
  213. REAL(wp), POINTER, DIMENSION(:,:) :: zqla ! latent heat flux over ocean
  214. REAL(wp), POINTER, DIMENSION(:,:) :: zqsb ! sensible heat flux over ocean
  215. !!---------------------------------------------------------------------
  216. !
  217. IF( nn_timing == 1 ) CALL timing_start('blk_oce_clio')
  218. !
  219. CALL wrk_alloc( jpi,jpj, zqlw, zqla, zqsb )
  220. zpatm = 101000._wp ! atmospheric pressure (assumed constant here)
  221. !------------------------------------!
  222. ! momentum fluxes (utau, vtau ) !
  223. !------------------------------------!
  224. !CDIR COLLAPSE
  225. utau(:,:) = sf(jp_utau)%fnow(:,:,1)
  226. !CDIR COLLAPSE
  227. vtau(:,:) = sf(jp_vtau)%fnow(:,:,1)
  228. ! begin TECLIM change
  229. CALL lbc_lnk( utau, 'U', -1. )
  230. CALL lbc_lnk( vtau, 'V', -1. )
  231. ! In the way we use CLIO forcings in TECLIM, utau and vtau are not the
  232. ! wind stress components as they should be, but the wind components. In
  233. ! this case (corresponding to sf(jp_utau)%clvar == 'uwnd') a special
  234. ! treatement is needed in the next section (ie a change from wind speed
  235. ! to wind stress over the ocean). This is inspired by what is done in the
  236. ! subroutine lim_sbc_tau in limsbc.F90.
  237. !
  238. ! Notes : * 1.3 = air density
  239. ! * cao is used
  240. !------------------------------------!
  241. ! wind stress module (taum ) !
  242. !------------------------------------!
  243. IF( sf(jp_utau)%clvar == 'uwnd' ) THEN
  244. !CDIR NOVERRCHK
  245. DO jj = 2, jpjm1
  246. !CDIR NOVERRCHK
  247. DO ji = fs_2, fs_jpim1 ! vector opt.
  248. ztx2 = utau(ji-1,jj ) + utau(ji,jj)
  249. zty2 = vtau(ji ,jj-1) + vtau(ji,jj)
  250. ztau1 = 0.25 * ( ztx2 * ztx2 + zty2 * zty2 )
  251. taum(ji,jj) = 1.3 * cao * ztau1
  252. ztau2(ji,jj) = 1.3 * cao * SQRT ( ztau1 )
  253. END DO
  254. END DO
  255. taum(:,:) = taum(:,:) * tmask(:,:,1)
  256. CALL lbc_lnk( taum, 'T', 1. )
  257. CALL lbc_lnk( ztau2, 'T', 1. )
  258. !CDIR NOVERRCHK
  259. DO jj = 2, jpjm1
  260. !CDIR NOVERRCHK
  261. DO ji = fs_2, fs_jpim1 ! vector opt.
  262. utau(ji,jj) = 0.5 * ( ztau2(ji,jj) + ztau2(ji+1,jj) ) * utau(ji,jj)
  263. vtau(ji,jj) = 0.5 * ( ztau2(ji,jj) + ztau2(ji,jj+1) ) * vtau(ji,jj)
  264. END DO
  265. END DO
  266. utau(:,:) = utau(:,:) * umask(:,:,1)
  267. vtau(:,:) = vtau(:,:) * vmask(:,:,1)
  268. CALL lbc_lnk( utau, 'U', -1. )
  269. CALL lbc_lnk( vtau, 'V', -1. )
  270. ELSE
  271. !CDIR NOVERRCHK
  272. DO jj = 2, jpjm1
  273. !CDIR NOVERRCHK
  274. DO ji = fs_2, fs_jpim1 ! vector opt.
  275. ztx2 = utau(ji-1,jj ) + utau(ji,jj)
  276. zty2 = vtau(ji ,jj-1) + vtau(ji,jj)
  277. taum(ji,jj) = 0.5 * SQRT( ztx2 * ztx2 + zty2 * zty2 )
  278. END DO
  279. END DO
  280. utau(:,:) = utau(:,:) * umask(:,:,1)
  281. vtau(:,:) = vtau(:,:) * vmask(:,:,1)
  282. taum(:,:) = taum(:,:) * tmask(:,:,1)
  283. CALL lbc_lnk( taum, 'T', 1. )
  284. ENDIF
  285. ! end TECLIM change
  286. !------------------------------------!
  287. ! store the wind speed (wndm ) !
  288. !------------------------------------!
  289. !CDIR COLLAPSE
  290. wndm(:,:) = sf(jp_wndm)%fnow(:,:,1)
  291. wndm(:,:) = wndm(:,:) * tmask(:,:,1)
  292. !------------------------------------------------!
  293. ! Shortwave radiation for ocean and snow/ice !
  294. !------------------------------------------------!
  295. CALL blk_clio_qsr_oce( qsr )
  296. qsr(:,:) = qsr(:,:) * tmask(:,:,1) ! no shortwave radiation into the ocean beneath ice shelf
  297. !------------------------!
  298. ! Other ocean fluxes !
  299. !------------------------!
  300. !CDIR NOVERRCHK
  301. !CDIR COLLAPSE
  302. DO jj = 1, jpj
  303. !CDIR NOVERRCHK
  304. DO ji = 1, jpi
  305. !
  306. zsst = pst(ji,jj) + rt0 ! converte Celcius to Kelvin the SST
  307. ztatm = sf(jp_tair)%fnow(ji,jj,1) ! and set minimum value far above 0 K (=rt0 over land)
  308. zcco1 = 1.0 - sf(jp_ccov)%fnow(ji,jj,1) ! fraction of clear sky ( 1 - cloud cover)
  309. zrhoa = zpatm / ( 287.04 * ztatm ) ! air density (equation of state for dry air)
  310. ztamr = ztatm - rtt ! Saturation water vapour
  311. zmt1 = SIGN( 17.269, ztamr ) ! ||
  312. zmt2 = SIGN( 21.875, ztamr ) ! \ /
  313. zmt3 = SIGN( 28.200, -ztamr ) ! \/
  314. zes = 611.0 * EXP( ABS( ztamr ) * MIN ( zmt1, zmt2 ) / ( ztatm - 35.86 + MAX( 0.e0, zmt3 ) ) )
  315. zev = sf(jp_humi)%fnow(ji,jj,1) * zes ! vapour pressure
  316. zevsqr = SQRT( zev * 0.01 ) ! square-root of vapour pressure
  317. zqatm = 0.622 * zev / ( zpatm - 0.378 * zev ) ! specific humidity
  318. !--------------------------------------!
  319. ! long-wave radiation over the ocean ! ( Berliand 1952 ; all latitudes )
  320. !--------------------------------------!
  321. ztatm3 = ztatm * ztatm * ztatm
  322. zcldeff = 1.0 - sbudyko(ji,jj) * sf(jp_ccov)%fnow(ji,jj,1) * sf(jp_ccov)%fnow(ji,jj,1)
  323. ztaevbk = ztatm * ztatm3 * zcldeff * ( 0.39 - 0.05 * zevsqr )
  324. !
  325. zqlw(ji,jj) = - emic * stefan * ( ztaevbk + 4. * ztatm3 * ( zsst - ztatm ) )
  326. !--------------------------------------------------
  327. ! Latent and sensible heat fluxes over the ocean
  328. !--------------------------------------------------
  329. ! ! vapour pressure at saturation of ocean
  330. zeso = 611.0 * EXP ( 17.2693884 * ( zsst - rtt ) * tmask(ji,jj,1) / ( zsst - 35.86 ) )
  331. zqsato = ( 0.622 * zeso ) / ( zpatm - 0.378 * zeso ) ! humidity close to the ocean surface (at saturation)
  332. ! Drag coefficients from Large and Pond (1981,1982)
  333. ! ! Stability parameters
  334. zdteta = zsst - ztatm
  335. zdeltaq = zqatm - zqsato
  336. ztvmoy = ztatm * ( 1. + 2.2e-3 * ztatm * zqatm )
  337. zdenum = MAX( sf(jp_wndm)%fnow(ji,jj,1) * sf(jp_wndm)%fnow(ji,jj,1) * ztvmoy, eps20 )
  338. zdtetar = zdteta / zdenum
  339. ztvmoyr = ztvmoy * ztvmoy * zdeltaq / zdenum
  340. ! ! case of stable atmospheric conditions
  341. zobouks = -70.0 * 10. * ( zdtetar + 3.2e-3 * ztvmoyr )
  342. zobouks = MAX( 0.e0, zobouks )
  343. zpsims = -7.0 * zobouks
  344. zpsihs = zpsims
  345. zpsils = zpsims
  346. ! ! case of unstable atmospheric conditions
  347. zobouku = MIN( 0.e0, -100.0 * 10.0 * ( zdtetar + 2.2e-3 * ztvmoyr ) )
  348. zxins = ( 1. - 16. * zobouku )**0.25
  349. zlxins = LOG( ( 1. + zxins * zxins ) / 2. )
  350. zpsimu = 2. * LOG( ( 1 + zxins ) * 0.5 ) + zlxins - 2. * ATAN( zxins ) + rpi * 0.5
  351. zpsihu = 2. * zlxins
  352. zpsilu = zpsihu
  353. ! ! intermediate values
  354. zstab = MAX( 0.e0, SIGN( 1.e0, zdteta ) )
  355. zpsim = zstab * zpsimu + ( 1.0 - zstab ) * zpsims
  356. zpsih = zstab * zpsihu + ( 1.0 - zstab ) * zpsihs
  357. zpsil = zpsih
  358. zvatmg = MAX( 0.032 * 1.5e-3 * sf(jp_wndm)%fnow(ji,jj,1) * sf(jp_wndm)%fnow(ji,jj,1) / grav, eps20 )
  359. zcmn = vkarmn / LOG ( 10. / zvatmg )
  360. zchn = 0.0327 * zcmn
  361. zcln = 0.0346 * zcmn
  362. zcmcmn = 1. / ( 1. - zcmn * zpsim / vkarmn )
  363. ! sometimes the ratio zchn * zpsih / ( vkarmn * zcmn ) is too close to 1 and zchcm becomes very very big
  364. zcmax = 0.1 ! choice for maximum value of the heat transfer coefficient, guided by my intuition
  365. zrmax = 1 - 3.e-4 / zcmax ! maximum value of the ratio
  366. zchcm = zcmcmn / ( 1. - MIN ( zchn * zpsih / ( vkarmn * zcmn ) , zrmax ) )
  367. zclcm = zchcm
  368. ! ! transfert coef. (Large and Pond 1981,1982)
  369. zcsho = zchn * zchcm
  370. zcleo = zcln * zclcm
  371. zrhova = zrhoa * sf(jp_wndm)%fnow(ji,jj,1)
  372. ! sensible heat flux
  373. zqsb(ji,jj) = zrhova * zcsho * 1004.0 * ( zsst - ztatm )
  374. ! latent heat flux (bounded by zero)
  375. zqla(ji,jj) = MAX( 0.e0, zrhova * zcleo * 2.5e+06 * ( zqsato - zqatm ) )
  376. !
  377. END DO
  378. END DO
  379. ! ----------------------------------------------------------------------------- !
  380. ! III Total FLUXES !
  381. ! ----------------------------------------------------------------------------- !
  382. zcevap = rcp / cevap ! convert zqla ==> evap (Kg/m2/s) ==> m/s ==> W/m2
  383. zcprec = rcp / rday ! convert prec ( mm/day ==> m/s) ==> W/m2
  384. !CDIR COLLAPSE
  385. emp(:,:) = zqla(:,:) / cevap & ! freshwater flux
  386. & - sf(jp_prec)%fnow(:,:,1) / rday * tmask(:,:,1)
  387. !
  388. !CDIR COLLAPSE
  389. qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:) & ! Downward Non Solar flux
  390. & - zqla(:,:) * pst(:,:) * zcevap & ! remove evap. heat content at SST in Celcius
  391. & + sf(jp_prec)%fnow(:,:,1) * sf(jp_tair)%fnow(:,:,1) * zcprec ! add precip. heat content at Tair in Celcius
  392. qns(:,:) = qns(:,:) * tmask(:,:,1)
  393. #if defined key_lim3
  394. qns_oce(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)
  395. qsr_oce(:,:) = qsr(:,:)
  396. #endif
  397. ! NB: if sea-ice model, the snow precip are computed and the associated heat is added to qns (see blk_ice_clio)
  398. IF ( nn_ice == 0 ) THEN
  399. CALL iom_put( "qlw_oce" , zqlw ) ! output downward longwave heat over the ocean
  400. CALL iom_put( "qsb_oce" , - zqsb ) ! output downward sensible heat over the ocean
  401. CALL iom_put( "qla_oce" , - zqla ) ! output downward latent heat over the ocean
  402. CALL iom_put( "qemp_oce", qns-zqlw+zqsb+zqla ) ! output downward heat content of E-P over the ocean
  403. CALL iom_put( "qns_oce" , qns ) ! output downward non solar heat over the ocean
  404. CALL iom_put( "qsr_oce" , qsr ) ! output downward solar heat over the ocean
  405. CALL iom_put( "qt_oce" , qns+qsr ) ! output total downward heat over the ocean
  406. ENDIF
  407. IF(ln_ctl) THEN
  408. CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_clio: zqsb : ', tab2d_2=zqlw , clinfo2=' zqlw : ')
  409. CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_clio: zqla : ', tab2d_2=qsr , clinfo2=' qsr : ')
  410. CALL prt_ctl(tab2d_1=pst , clinfo1=' blk_oce_clio: pst : ', tab2d_2=emp , clinfo2=' emp : ')
  411. CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_clio: utau : ', mask1=umask, &
  412. & tab2d_2=vtau , clinfo2=' vtau : ', mask2=vmask )
  413. ENDIF
  414. CALL wrk_dealloc( jpi,jpj, zqlw, zqla, zqsb )
  415. !
  416. IF( nn_timing == 1 ) CALL timing_stop('blk_oce_clio')
  417. !
  418. END SUBROUTINE blk_oce_clio
  419. # if defined key_lim2 || defined key_lim3
  420. SUBROUTINE blk_ice_clio_tau
  421. !!---------------------------------------------------------------------------
  422. !! *** ROUTINE blk_ice_clio_tau ***
  423. !!
  424. !! ** Purpose : Computation momentum flux at the ice-atm interface
  425. !!
  426. !! ** Method : Read utau from a forcing file. Rearrange if C-grid
  427. !!
  428. !!----------------------------------------------------------------------
  429. REAL(wp) :: zcoef
  430. INTEGER :: ji, jj ! dummy loop indices
  431. !!---------------------------------------------------------------------
  432. !
  433. IF( nn_timing == 1 ) CALL timing_start('blk_ice_clio_tau')
  434. SELECT CASE( cp_ice_msh )
  435. CASE( 'C' ) ! C-grid ice dynamics
  436. zcoef = cai / cao ! Change from air-sea stress to air-ice stress
  437. utau_ice(:,:) = zcoef * utau(:,:)
  438. vtau_ice(:,:) = zcoef * vtau(:,:)
  439. CASE( 'I' ) ! I-grid ice dynamics: I-point (i.e. F-point lower-left corner)
  440. zcoef = 0.5_wp * cai / cao ! Change from air-sea stress to air-ice stress
  441. DO jj = 2, jpj ! stress from ocean U- and V-points to ice U,V point
  442. DO ji = 2, jpi ! I-grid : no vector opt.
  443. utau_ice(ji,jj) = zcoef * ( utau(ji-1,jj ) + utau(ji-1,jj-1) )
  444. vtau_ice(ji,jj) = zcoef * ( vtau(ji ,jj-1) + vtau(ji-1,jj-1) )
  445. END DO
  446. END DO
  447. CALL lbc_lnk( utau_ice(:,:), 'I', -1. ) ; CALL lbc_lnk( vtau_ice(:,:), 'I', -1. ) ! I-point
  448. END SELECT
  449. IF(ln_ctl) THEN
  450. CALL prt_ctl(tab2d_1=utau_ice , clinfo1=' blk_ice_clio: utau_ice : ', tab2d_2=vtau_ice , clinfo2=' vtau_ice : ')
  451. ENDIF
  452. IF( nn_timing == 1 ) CALL timing_stop('blk_ice_clio_tau')
  453. END SUBROUTINE blk_ice_clio_tau
  454. #endif
  455. # if defined key_lim2 || defined key_lim3
  456. SUBROUTINE blk_ice_clio_flx( ptsu , palb_cs, palb_os, palb )
  457. !!---------------------------------------------------------------------------
  458. !! *** ROUTINE blk_ice_clio_flx ***
  459. !!
  460. !! ** Purpose : Computation of the heat fluxes at ocean and snow/ice
  461. !! surface the solar heat at ocean and snow/ice surfaces and the
  462. !! sensitivity of total heat fluxes to the SST variations
  463. !!
  464. !! ** Method : The flux of heat at the ice and ocean surfaces are derived
  465. !! from semi-empirical ( or bulk ) formulae which relate the flux to
  466. !! the properties of the surface and of the lower atmosphere. Here, we
  467. !! follow the work of Oberhuber, 1988
  468. !!
  469. !! ** Action : call albedo_oce/albedo_ice to compute ocean/ice albedo
  470. !! - snow precipitation
  471. !! - solar flux at the ocean and ice surfaces
  472. !! - the long-wave radiation for the ocean and sea/ice
  473. !! - turbulent heat fluxes over water and ice
  474. !! - evaporation over water
  475. !! - total heat fluxes sensitivity over ice (dQ/dT)
  476. !! - latent heat flux sensitivity over ice (dQla/dT)
  477. !! - qns : modified the non solar heat flux over the ocean
  478. !! to take into account solid precip latent heat flux
  479. !!----------------------------------------------------------------------
  480. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ptsu ! ice surface temperature [Kelvin]
  481. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: palb_cs ! ice albedo (clear sky) (alb_ice_cs) [-]
  482. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: palb_os ! ice albedo (overcast sky) (alb_ice_os) [-]
  483. REAL(wp), INTENT( out), DIMENSION(:,:,:) :: palb ! ice albedo (actual value) [-]
  484. !!
  485. INTEGER :: ji, jj, jl ! dummy loop indices
  486. !!
  487. REAL(wp) :: zmt1, zmt2, zmt3, ztatm3 ! temporary scalars
  488. REAL(wp) :: ztaevbk, zind1, zind2, zind3, ztamr ! - -
  489. REAL(wp) :: zesi, zqsati, zdesidt ! - -
  490. REAL(wp) :: zdqla, zcldeff, zev, zes, zpatm, zrhova ! - -
  491. REAL(wp) :: zcshi, zclei, zrhovaclei, zrhovacshi ! - -
  492. REAL(wp) :: ztice3, zticemb, zticemb2, zdqlw, zdqsb ! - -
  493. REAL(wp) :: z1_lsub ! - -
  494. !!
  495. REAL(wp), DIMENSION(:,:) , POINTER :: ztatm ! Tair in Kelvin
  496. REAL(wp), DIMENSION(:,:) , POINTER :: zqatm ! specific humidity
  497. REAL(wp), DIMENSION(:,:) , POINTER :: zevsqr ! vapour pressure square-root
  498. REAL(wp), DIMENSION(:,:) , POINTER :: zrhoa ! air density
  499. REAL(wp), DIMENSION(:,:,:), POINTER :: z_qlw, z_qsb
  500. REAL(wp), DIMENSION(:,:) , POINTER :: zevap, zsnw
  501. !!---------------------------------------------------------------------
  502. !
  503. IF( nn_timing == 1 ) CALL timing_start('blk_ice_clio_flx')
  504. !
  505. CALL wrk_alloc( jpi,jpj, ztatm, zqatm, zevsqr, zrhoa )
  506. CALL wrk_alloc( jpi,jpj, jpl, z_qlw, z_qsb )
  507. zpatm = 101000. ! atmospheric pressure (assumed constant here)
  508. !--------------------------------------------------------------------------------
  509. ! Determine cloud optical depths as a function of latitude (Chou et al., 1981).
  510. ! and the correction factor for taking into account the effect of clouds
  511. !--------------------------------------------------------------------------------
  512. !CDIR NOVERRCHK
  513. !CDIR COLLAPSE
  514. DO jj = 1, jpj
  515. !CDIR NOVERRCHK
  516. DO ji = 1, jpi
  517. ztatm (ji,jj) = sf(jp_tair)%fnow(ji,jj,1) ! air temperature in Kelvins
  518. zrhoa(ji,jj) = zpatm / ( 287.04 * ztatm(ji,jj) ) ! air density (equation of state for dry air)
  519. ztamr = ztatm(ji,jj) - rtt ! Saturation water vapour
  520. zmt1 = SIGN( 17.269, ztamr )
  521. zmt2 = SIGN( 21.875, ztamr )
  522. zmt3 = SIGN( 28.200, -ztamr )
  523. zes = 611.0 * EXP( ABS( ztamr ) * MIN ( zmt1, zmt2 ) &
  524. & / ( ztatm(ji,jj) - 35.86 + MAX( 0.e0, zmt3 ) ) )
  525. zev = sf(jp_humi)%fnow(ji,jj,1) * zes ! vapour pressure
  526. zevsqr(ji,jj) = SQRT( zev * 0.01 ) ! square-root of vapour pressure
  527. zqatm(ji,jj) = 0.622 * zev / ( zpatm - 0.378 * zev ) ! specific humidity
  528. !----------------------------------------------------
  529. ! Computation of snow precipitation (Ledley, 1985) |
  530. !----------------------------------------------------
  531. zmt1 = 253.0 - ztatm(ji,jj) ; zind1 = MAX( 0.e0, SIGN( 1.e0, zmt1 ) )
  532. zmt2 = ( 272.0 - ztatm(ji,jj) ) / 38.0 ; zind2 = MAX( 0.e0, SIGN( 1.e0, zmt2 ) )
  533. zmt3 = ( 281.0 - ztatm(ji,jj) ) / 18.0 ; zind3 = MAX( 0.e0, SIGN( 1.e0, zmt3 ) )
  534. sprecip(ji,jj) = sf(jp_prec)%fnow(ji,jj,1) / rday & ! rday = converte mm/day to kg/m2/s
  535. & * ( zind1 & ! solid (snow) precipitation [kg/m2/s]
  536. & + ( 1.0 - zind1 ) * ( zind2 * ( 0.5 + zmt2 ) &
  537. & + ( 1.0 - zind2 ) * zind3 * zmt3 ) )
  538. !----------------------------------------------------!
  539. ! fraction of net penetrative shortwave radiation !
  540. !----------------------------------------------------!
  541. ! fraction of qsr_ice which is NOT absorbed in the thin surface layer
  542. ! and thus which penetrates inside the ice cover ( Maykut and Untersteiner, 1971 ; Elbert anbd Curry, 1993 )
  543. fr1_i0(ji,jj) = 0.18 * ( 1.e0 - sf(jp_ccov)%fnow(ji,jj,1) ) + 0.35 * sf(jp_ccov)%fnow(ji,jj,1)
  544. fr2_i0(ji,jj) = 0.82 * ( 1.e0 - sf(jp_ccov)%fnow(ji,jj,1) ) + 0.65 * sf(jp_ccov)%fnow(ji,jj,1)
  545. END DO
  546. END DO
  547. CALL iom_put( 'snowpre', sprecip ) ! Snow precipitation
  548. !-----------------------------------------------------------!
  549. ! snow/ice Shortwave radiation (abedo already computed) !
  550. !-----------------------------------------------------------!
  551. CALL blk_clio_qsr_ice( palb_cs, palb_os, qsr_ice )
  552. DO jl = 1, jpl
  553. palb(:,:,jl) = ( palb_cs(:,:,jl) * ( 1.e0 - sf(jp_ccov)%fnow(:,:,1) ) &
  554. & + palb_os(:,:,jl) * sf(jp_ccov)%fnow(:,:,1) )
  555. END DO
  556. ! ! ========================== !
  557. DO jl = 1, jpl ! Loop over ice categories !
  558. ! ! ========================== !
  559. !CDIR NOVERRCHK
  560. !CDIR COLLAPSE
  561. DO jj = 1 , jpj
  562. !CDIR NOVERRCHK
  563. DO ji = 1, jpi
  564. !-------------------------------------------!
  565. ! long-wave radiation over ice categories ! ( Berliand 1952 ; all latitudes )
  566. !-------------------------------------------!
  567. ztatm3 = ztatm(ji,jj) * ztatm(ji,jj) * ztatm(ji,jj)
  568. zcldeff = 1.0 - sbudyko(ji,jj) * sf(jp_ccov)%fnow(ji,jj,1) * sf(jp_ccov)%fnow(ji,jj,1)
  569. ztaevbk = ztatm3 * ztatm(ji,jj) * zcldeff * ( 0.39 - 0.05 * zevsqr(ji,jj) )
  570. !
  571. z_qlw(ji,jj,jl) = - emic * stefan * ( ztaevbk + 4. * ztatm3 * ( ptsu(ji,jj,jl) - ztatm(ji,jj) ) )
  572. !----------------------------------------
  573. ! Turbulent heat fluxes over snow/ice ( Latent and sensible )
  574. !----------------------------------------
  575. ! vapour pressure at saturation of ice (tmask to avoid overflow in the exponential)
  576. zesi = 611.0 * EXP( 21.8745587 * tmask(ji,jj,1) * ( ptsu(ji,jj,jl) - rtt )/ ( ptsu(ji,jj,jl) - 7.66 ) )
  577. ! humidity close to the ice surface (at saturation)
  578. zqsati = ( 0.622 * zesi ) / ( zpatm - 0.378 * zesi )
  579. ! computation of intermediate values
  580. zticemb = ptsu(ji,jj,jl) - 7.66
  581. zticemb2 = zticemb * zticemb
  582. ztice3 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
  583. zdesidt = zesi * ( 9.5 * LOG( 10.0 ) * ( rtt - 7.66 ) / zticemb2 )
  584. ! Transfer cofficients assumed to be constant (Parkinson 1979 ; Maykut 1982)
  585. zcshi = 1.75e-03
  586. zclei = zcshi
  587. ! sensible and latent fluxes over ice
  588. zrhova = zrhoa(ji,jj) * sf(jp_wndm)%fnow(ji,jj,1) ! computation of intermediate values
  589. zrhovaclei = zrhova * zcshi * 2.834e+06
  590. zrhovacshi = zrhova * zclei * 1004.0
  591. ! sensible heat flux
  592. z_qsb(ji,jj,jl) = zrhovacshi * ( ptsu(ji,jj,jl) - ztatm(ji,jj) )
  593. ! latent heat flux
  594. qla_ice(ji,jj,jl) = MAX( 0.e0, zrhovaclei * ( zqsati - zqatm(ji,jj) ) )
  595. ! sensitivity of non solar fluxes (dQ/dT) (long-wave, sensible and latent fluxes)
  596. zdqlw = 4.0 * emic * stefan * ztice3
  597. zdqsb = zrhovacshi
  598. zdqla = zrhovaclei * ( zdesidt * ( zqsati * zqsati / ( zesi * zesi ) ) * ( zpatm / 0.622 ) )
  599. !
  600. dqla_ice(ji,jj,jl) = zdqla ! latent flux sensitivity
  601. dqns_ice(ji,jj,jl) = -( zdqlw + zdqsb + zdqla ) ! total non solar sensitivity
  602. END DO
  603. !
  604. END DO
  605. !
  606. END DO
  607. !
  608. ! ----------------------------------------------------------------------------- !
  609. ! Total FLUXES !
  610. ! ----------------------------------------------------------------------------- !
  611. !
  612. !CDIR COLLAPSE
  613. qns_ice(:,:,:) = z_qlw (:,:,:) - z_qsb (:,:,:) - qla_ice (:,:,:) ! Downward Non Solar flux
  614. !CDIR COLLAPSE
  615. tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) / rday ! total precipitation [kg/m2/s]
  616. !
  617. ! ----------------------------------------------------------------------------- !
  618. ! Correct the OCEAN non solar flux with the existence of solid precipitation !
  619. ! ---------------=====--------------------------------------------------------- !
  620. !CDIR COLLAPSE
  621. qns(:,:) = qns(:,:) & ! update the non-solar heat flux with:
  622. & - sprecip(:,:) * lfus & ! remove melting solid precip
  623. & + sprecip(:,:) * MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow - rt0 ) * cpic & ! add solid P at least below melting
  624. & - sprecip(:,:) * sf(jp_tair)%fnow(:,:,1) * rcp ! remove solid precip. at Tair
  625. #if defined key_lim3
  626. ! ----------------------------------------------------------------------------- !
  627. ! Distribute evapo, precip & associated heat over ice and ocean
  628. ! ---------------=====--------------------------------------------------------- !
  629. CALL wrk_alloc( jpi,jpj, zevap, zsnw )
  630. ! --- evaporation --- !
  631. z1_lsub = 1._wp / Lsub
  632. evap_ice (:,:,:) = qla_ice (:,:,:) * z1_lsub ! sublimation
  633. devap_ice(:,:,:) = dqla_ice(:,:,:) * z1_lsub
  634. zevap (:,:) = emp(:,:) + tprecip(:,:) ! evaporation over ocean
  635. ! --- evaporation minus precipitation --- !
  636. zsnw(:,:) = 0._wp
  637. CALL lim_thd_snwblow( pfrld, zsnw ) ! snow redistribution by wind
  638. emp_oce(:,:) = pfrld(:,:) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * ( 1._wp - zsnw )
  639. emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
  640. emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
  641. ! --- heat flux associated with emp --- !
  642. qemp_oce(:,:) = - pfrld(:,:) * zevap(:,:) * sst_m(:,:) * rcp & ! evap
  643. & + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp & ! liquid precip
  644. & + sprecip(:,:) * ( 1._wp - zsnw ) * & ! solid precip
  645. & ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
  646. qemp_ice(:,:) = sprecip(:,:) * zsnw * & ! solid precip (only)
  647. & ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
  648. ! --- total solar and non solar fluxes --- !
  649. qns_tot(:,:) = pfrld(:,:) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 ) + qemp_ice(:,:) + qemp_oce(:,:)
  650. qsr_tot(:,:) = pfrld(:,:) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
  651. ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
  652. qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
  653. ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) --- !
  654. DO jl = 1, jpl
  655. qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * cpic * tmask(:,:,1) - lfus )
  656. ! but then qemp_ice should also include sublimation
  657. END DO
  658. CALL wrk_dealloc( jpi,jpj, zevap, zsnw )
  659. #endif
  660. !!gm : not necessary as all input data are lbc_lnk...
  661. CALL lbc_lnk( fr1_i0 (:,:) , 'T', 1. )
  662. CALL lbc_lnk( fr2_i0 (:,:) , 'T', 1. )
  663. DO jl = 1, jpl
  664. CALL lbc_lnk( qns_ice (:,:,jl) , 'T', 1. )
  665. CALL lbc_lnk( dqns_ice(:,:,jl) , 'T', 1. )
  666. CALL lbc_lnk( qla_ice (:,:,jl) , 'T', 1. )
  667. CALL lbc_lnk( dqla_ice(:,:,jl) , 'T', 1. )
  668. END DO
  669. !!gm : mask is not required on forcing
  670. DO jl = 1, jpl
  671. qns_ice (:,:,jl) = qns_ice (:,:,jl) * tmask(:,:,1)
  672. qla_ice (:,:,jl) = qla_ice (:,:,jl) * tmask(:,:,1)
  673. dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * tmask(:,:,1)
  674. dqla_ice(:,:,jl) = dqla_ice(:,:,jl) * tmask(:,:,1)
  675. END DO
  676. CALL wrk_dealloc( jpi,jpj, ztatm, zqatm, zevsqr, zrhoa )
  677. CALL wrk_dealloc( jpi,jpj, jpl , z_qlw, z_qsb )
  678. IF(ln_ctl) THEN
  679. CALL prt_ctl(tab3d_1=z_qsb , clinfo1=' blk_ice_clio: z_qsb : ', tab3d_2=z_qlw , clinfo2=' z_qlw : ', kdim=jpl)
  680. CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice_clio: z_qla : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice : ', kdim=jpl)
  681. CALL prt_ctl(tab3d_1=dqns_ice , clinfo1=' blk_ice_clio: dqns_ice : ', tab3d_2=qns_ice , clinfo2=' qns_ice : ', kdim=jpl)
  682. CALL prt_ctl(tab3d_1=dqla_ice , clinfo1=' blk_ice_clio: dqla_ice : ', tab3d_2=ptsu , clinfo2=' ptsu : ', kdim=jpl)
  683. CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice_clio: tprecip : ', tab2d_2=sprecip , clinfo2=' sprecip : ')
  684. ENDIF
  685. IF( nn_timing == 1 ) CALL timing_stop('blk_ice_clio_flx')
  686. !
  687. END SUBROUTINE blk_ice_clio_flx
  688. #endif
  689. SUBROUTINE blk_clio_qsr_oce( pqsr_oce )
  690. !!---------------------------------------------------------------------------
  691. !! *** ROUTINE blk_clio_qsr_oce ***
  692. !!
  693. !! ** Purpose : Computation of the shortwave radiation at the ocean and the
  694. !! snow/ice surfaces.
  695. !!
  696. !! ** Method : - computed qsr from the cloud cover for both ice and ocean
  697. !! - also initialise sbudyko and stauc once for all
  698. !!----------------------------------------------------------------------
  699. REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: pqsr_oce ! shortwave radiation over the ocean
  700. !!
  701. INTEGER, PARAMETER :: jp24 = 24 ! sampling of the daylight period (sunrise to sunset) into 24 equal parts
  702. !!
  703. INTEGER :: ji, jj, jt ! dummy loop indices
  704. INTEGER :: indaet ! = -1, 0, 1 for odd, normal and leap years resp.
  705. INTEGER :: iday ! integer part of day
  706. INTEGER :: indxb, indxc ! index for cloud depth coefficient
  707. REAL(wp) :: zalat , zclat, zcmue, zcmue2 ! local scalars
  708. REAL(wp) :: zmt1, zmt2, zmt3 !
  709. REAL(wp) :: zdecl, zsdecl , zcdecl !
  710. REAL(wp) :: za_oce, ztamr !
  711. REAL(wp) :: zdl, zlha ! local scalars
  712. REAL(wp) :: zlmunoon, zcldcor, zdaycor !
  713. REAL(wp) :: zxday, zdist, zcoef, zcoef1 !
  714. REAL(wp) :: zes
  715. REAL(wp), DIMENSION(:,:), POINTER :: zev ! vapour pressure
  716. REAL(wp), DIMENSION(:,:), POINTER :: zdlha, zlsrise, zlsset ! 2D workspace
  717. REAL(wp), DIMENSION(:,:), POINTER :: zps, zpc ! sine (cosine) of latitude per sine (cosine) of solar declination
  718. !!---------------------------------------------------------------------
  719. !
  720. IF( nn_timing == 1 ) CALL timing_start('blk_clio_qsr_oce')
  721. !
  722. CALL wrk_alloc( jpi,jpj, zev, zdlha, zlsrise, zlsset, zps, zpc )
  723. IF( lbulk_init ) THEN ! Initilization at first time step only
  724. rdtbs2 = nn_fsbc * rdt * 0.5
  725. ! cloud optical depths as a function of latitude (Chou et al., 1981).
  726. ! and the correction factor for taking into account the effect of clouds
  727. DO jj = 1, jpj
  728. DO ji = 1 , jpi
  729. zalat = ( 90.e0 - ABS( gphit(ji,jj) ) ) / 5.e0
  730. zclat = ( 95.e0 - gphit(ji,jj) ) / 10.e0
  731. indxb = 1 + INT( zalat )
  732. indxc = 1 + INT( zclat )
  733. zdl = zclat - INT( zclat )
  734. ! correction factor to account for the effect of clouds
  735. sbudyko(ji,jj) = budyko(indxb)
  736. stauc (ji,jj) = ( 1.e0 - zdl ) * tauco( indxc ) + zdl * tauco( indxc + 1 )
  737. END DO
  738. END DO
  739. lbulk_init = .FALSE.
  740. ENDIF
  741. ! Saturated water vapour and vapour pressure
  742. ! ------------------------------------------
  743. !CDIR NOVERRCHK
  744. !CDIR COLLAPSE
  745. DO jj = 1, jpj
  746. !CDIR NOVERRCHK
  747. DO ji = 1, jpi
  748. ztamr = sf(jp_tair)%fnow(ji,jj,1) - rtt
  749. zmt1 = SIGN( 17.269, ztamr )
  750. zmt2 = SIGN( 21.875, ztamr )
  751. zmt3 = SIGN( 28.200, -ztamr )
  752. zes = 611.0 * EXP( ABS( ztamr ) * MIN ( zmt1, zmt2 ) & ! Saturation water vapour
  753. & / ( sf(jp_tair)%fnow(ji,jj,1) - 35.86 + MAX( 0.e0, zmt3 ) ) )
  754. zev(ji,jj) = sf(jp_humi)%fnow(ji,jj,1) * zes * 1.0e-05 ! vapour pressure
  755. END DO
  756. END DO
  757. !-----------------------------------!
  758. ! Computation of solar irradiance !
  759. !-----------------------------------!
  760. !!gm : hard coded leap year ???
  761. indaet = 1 ! = -1, 0, 1 for odd, normal and leap years resp.
  762. zxday = nday_year + rdtbs2 / rday ! day of the year at which the fluxes are calculated
  763. iday = INT( zxday ) ! (centred at the middle of the ice time step)
  764. CALL flx_blk_declin( indaet, iday, zdecl ) ! solar declination of the current day
  765. zsdecl = SIN( zdecl * rad ) ! its sine
  766. zcdecl = COS( zdecl * rad ) ! its cosine
  767. ! correction factor added for computation of shortwave flux to take into account the variation of
  768. ! the distance between the sun and the earth during the year (Oberhuber 1988)
  769. zdist = zxday * 2. * rpi / REAL(nyear_len(1), wp)
  770. zdaycor = 1.0 + 0.0013 * SIN( zdist ) + 0.0342 * COS( zdist )
  771. !CDIR NOVERRCHK
  772. DO jj = 1, jpj
  773. !CDIR NOVERRCHK
  774. DO ji = 1, jpi
  775. ! product of sine (cosine) of latitude and sine (cosine) of solar declination
  776. zps(ji,jj) = SIN( gphit(ji,jj) * rad ) * zsdecl
  777. zpc(ji,jj) = COS( gphit(ji,jj) * rad ) * zcdecl
  778. ! computation of the both local time of sunrise and sunset
  779. zlsrise(ji,jj) = ACOS( - SIGN( 1.e0, zps(ji,jj) ) &
  780. & * MIN( 1.e0, SIGN( 1.e0, zps(ji,jj) ) * ( zps(ji,jj) / zpc(ji,jj) ) ) )
  781. zlsset (ji,jj) = - zlsrise(ji,jj)
  782. ! dividing the solar day into jp24 segments of length zdlha
  783. zdlha (ji,jj) = ( zlsrise(ji,jj) - zlsset(ji,jj) ) / REAL( jp24, wp )
  784. END DO
  785. END DO
  786. !---------------------------------------------!
  787. ! shortwave radiation absorbed by the ocean !
  788. !---------------------------------------------!
  789. pqsr_oce(:,:) = 0.e0 ! set ocean qsr to zero
  790. ! compute and sum ocean qsr over the daylight (i.e. between sunrise and sunset)
  791. !CDIR NOVERRCHK
  792. DO jt = 1, jp24
  793. zcoef = FLOAT( jt ) - 0.5
  794. !CDIR NOVERRCHK
  795. !CDIR COLLAPSE
  796. DO jj = 1, jpj
  797. !CDIR NOVERRCHK
  798. DO ji = 1, jpi
  799. zlha = COS( zlsrise(ji,jj) - zcoef * zdlha(ji,jj) ) ! local hour angle
  800. zcmue = MAX( 0.e0 , zps(ji,jj) + zpc(ji,jj) * zlha ) ! cos of local solar altitude
  801. zcmue2 = 1368.0 * zcmue * zcmue
  802. ! ocean albedo depending on the cloud cover (Payne, 1972)
  803. za_oce = ( 1.0 - sf(jp_ccov)%fnow(ji,jj,1) ) * 0.05 / ( 1.1 * zcmue**1.4 + 0.15 ) & ! clear sky
  804. & + sf(jp_ccov)%fnow(ji,jj,1) * 0.06 ! overcast
  805. ! solar heat flux absorbed by the ocean (Zillman, 1972)
  806. pqsr_oce(ji,jj) = pqsr_oce(ji,jj) &
  807. & + ( 1.0 - za_oce ) * zdlha(ji,jj) * zcmue2 &
  808. & / ( ( zcmue + 2.7 ) * zev(ji,jj) + 1.085 * zcmue + 0.10 )
  809. END DO
  810. END DO
  811. END DO
  812. ! Taking into account the ellipsity of the earth orbit, the clouds AND masked if sea-ice cover > 0%
  813. zcoef1 = srgamma * zdaycor / ( 2. * rpi )
  814. !CDIR COLLAPSE
  815. DO jj = 1, jpj
  816. DO ji = 1, jpi
  817. zlmunoon = ASIN( zps(ji,jj) + zpc(ji,jj) ) / rad ! local noon solar altitude
  818. zcldcor = MIN( 1.e0, ( 1.e0 - 0.62 * sf(jp_ccov)%fnow(ji,jj,1) & ! cloud correction (Reed 1977)
  819. & + 0.0019 * zlmunoon ) )
  820. pqsr_oce(ji,jj) = zcoef1 * zcldcor * pqsr_oce(ji,jj) * tmask(ji,jj,1) ! and zcoef1: ellipsity
  821. END DO
  822. END DO
  823. CALL wrk_dealloc( jpi,jpj, zev, zdlha, zlsrise, zlsset, zps, zpc )
  824. !
  825. IF( nn_timing == 1 ) CALL timing_stop('blk_clio_qsr_oce')
  826. !
  827. END SUBROUTINE blk_clio_qsr_oce
  828. SUBROUTINE blk_clio_qsr_ice( pa_ice_cs, pa_ice_os, pqsr_ice )
  829. !!---------------------------------------------------------------------------
  830. !! *** ROUTINE blk_clio_qsr_ice ***
  831. !!
  832. !! ** Purpose : Computation of the shortwave radiation at the ocean and the
  833. !! snow/ice surfaces.
  834. !!
  835. !! ** Method : - computed qsr from the cloud cover for both ice and ocean
  836. !! - also initialise sbudyko and stauc once for all
  837. !!----------------------------------------------------------------------
  838. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: pa_ice_cs ! albedo of ice under clear sky
  839. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: pa_ice_os ! albedo of ice under overcast sky
  840. REAL(wp), INTENT( out), DIMENSION(:,:,:) :: pqsr_ice ! shortwave radiation over the ice/snow
  841. !!
  842. INTEGER, PARAMETER :: jp24 = 24 ! sampling of the daylight period (sunrise to sunset) into 24 equal parts
  843. !!
  844. INTEGER :: ji, jj, jl, jt ! dummy loop indices
  845. INTEGER :: ijpl ! number of ice categories (3rd dim of pqsr_ice)
  846. INTEGER :: indaet ! = -1, 0, 1 for odd, normal and leap years resp.
  847. INTEGER :: iday ! integer part of day
  848. !!
  849. REAL(wp) :: zcmue, zcmue2, ztamr ! temporary scalars
  850. REAL(wp) :: zmt1, zmt2, zmt3 ! - -
  851. REAL(wp) :: zdecl, zsdecl, zcdecl ! - -
  852. REAL(wp) :: zlha, zdaycor, zes ! - -
  853. REAL(wp) :: zxday, zdist, zcoef, zcoef1 ! - -
  854. REAL(wp) :: zqsr_ice_cs, zqsr_ice_os ! - -
  855. REAL(wp), DIMENSION(:,:), POINTER :: zev ! vapour pressure
  856. REAL(wp), DIMENSION(:,:), POINTER :: zdlha, zlsrise, zlsset ! 2D workspace
  857. REAL(wp), DIMENSION(:,:), POINTER :: zps, zpc ! sine (cosine) of latitude per sine (cosine) of solar declination
  858. !!---------------------------------------------------------------------
  859. !
  860. IF( nn_timing == 1 ) CALL timing_start('blk_clio_qsr_ice')
  861. !
  862. CALL wrk_alloc( jpi,jpj, zev, zdlha, zlsrise, zlsset, zps, zpc )
  863. ijpl = SIZE(pqsr_ice, 3 ) ! number of ice categories
  864. ! Saturated water vapour and vapour pressure
  865. ! ------------------------------------------
  866. !CDIR NOVERRCHK
  867. !CDIR COLLAPSE
  868. DO jj = 1, jpj
  869. !CDIR NOVERRCHK
  870. DO ji = 1, jpi
  871. ztamr = sf(jp_tair)%fnow(ji,jj,1) - rtt
  872. zmt1 = SIGN( 17.269, ztamr )
  873. zmt2 = SIGN( 21.875, ztamr )
  874. zmt3 = SIGN( 28.200, -ztamr )
  875. zes = 611.0 * EXP( ABS( ztamr ) * MIN ( zmt1, zmt2 ) & ! Saturation water vapour
  876. & / ( sf(jp_tair)%fnow(ji,jj,1) - 35.86 + MAX( 0.e0, zmt3 ) ) )
  877. zev(ji,jj) = sf(jp_humi)%fnow(ji,jj,1) * zes * 1.0e-05 ! vapour pressure
  878. END DO
  879. END DO
  880. !-----------------------------------!
  881. ! Computation of solar irradiance !
  882. !-----------------------------------!
  883. !!gm : hard coded leap year ???
  884. indaet = 1 ! = -1, 0, 1 for odd, normal and leap years resp.
  885. zxday = nday_year + rdtbs2 / rday ! day of the year at which the fluxes are calculated
  886. iday = INT( zxday ) ! (centred at the middle of the ice time step)
  887. CALL flx_blk_declin( indaet, iday, zdecl ) ! solar declination of the current day
  888. zsdecl = SIN( zdecl * rad ) ! its sine
  889. zcdecl = COS( zdecl * rad ) ! its cosine
  890. ! correction factor added for computation of shortwave flux to take into account the variation of
  891. ! the distance between the sun and the earth during the year (Oberhuber 1988)
  892. zdist = zxday * 2. * rpi / REAL(nyear_len(1), wp)
  893. zdaycor = 1.0 + 0.0013 * SIN( zdist ) + 0.0342 * COS( zdist )
  894. !CDIR NOVERRCHK
  895. DO jj = 1, jpj
  896. !CDIR NOVERRCHK
  897. DO ji = 1, jpi
  898. ! product of sine (cosine) of latitude and sine (cosine) of solar declination
  899. zps(ji,jj) = SIN( gphit(ji,jj) * rad ) * zsdecl
  900. zpc(ji,jj) = COS( gphit(ji,jj) * rad ) * zcdecl
  901. ! computation of the both local time of sunrise and sunset
  902. zlsrise(ji,jj) = ACOS( - SIGN( 1.e0, zps(ji,jj) ) &
  903. & * MIN( 1.e0, SIGN( 1.e0, zps(ji,jj) ) * ( zps(ji,jj) / zpc(ji,jj) ) ) )
  904. zlsset (ji,jj) = - zlsrise(ji,jj)
  905. ! dividing the solar day into jp24 segments of length zdlha
  906. zdlha (ji,jj) = ( zlsrise(ji,jj) - zlsset(ji,jj) ) / REAL( jp24, wp )
  907. END DO
  908. END DO
  909. !---------------------------------------------!
  910. ! shortwave radiation absorbed by the ice !
  911. !---------------------------------------------!
  912. ! compute and sum ice qsr over the daylight for each ice categories
  913. pqsr_ice(:,:,:) = 0.e0
  914. zcoef1 = zdaycor / ( 2. * rpi ) ! Correction for the ellipsity of the earth orbit
  915. ! !----------------------------!
  916. DO jl = 1, ijpl ! loop over ice categories !
  917. ! !----------------------------!
  918. !CDIR NOVERRCHK
  919. DO jt = 1, jp24
  920. zcoef = FLOAT( jt ) - 0.5
  921. !CDIR NOVERRCHK
  922. !CDIR COLLAPSE
  923. DO jj = 1, jpj
  924. !CDIR NOVERRCHK
  925. DO ji = 1, jpi
  926. zlha = COS( zlsrise(ji,jj) - zcoef * zdlha(ji,jj) ) ! local hour angle
  927. zcmue = MAX( 0.e0 , zps(ji,jj) + zpc(ji,jj) * zlha ) ! cos of local solar altitude
  928. zcmue2 = 1368.0 * zcmue * zcmue
  929. ! solar heat flux absorbed by the ice/snow system (Shine and Crane 1984 adapted to high albedo)
  930. zqsr_ice_cs = ( 1.0 - pa_ice_cs(ji,jj,jl) ) * zdlha(ji,jj) * zcmue2 & ! clear sky
  931. & / ( ( 1.0 + zcmue ) * zev(ji,jj) + 1.2 * zcmue + 0.0455 )
  932. zqsr_ice_os = zdlha(ji,jj) * SQRT( zcmue ) & ! overcast sky
  933. & * ( 53.5 + 1274.5 * zcmue ) * ( 1.0 - 0.996 * pa_ice_os(ji,jj,jl) ) &
  934. & / ( 1.0 + 0.139 * stauc(ji,jj) * ( 1.0 - 0.9435 * pa_ice_os(ji,jj,jl) ) )
  935. pqsr_ice(ji,jj,jl) = pqsr_ice(ji,jj,jl) + ( ( 1.0 - sf(jp_ccov)%fnow(ji,jj,1) ) * zqsr_ice_cs &
  936. & + sf(jp_ccov)%fnow(ji,jj,1) * zqsr_ice_os )
  937. END DO
  938. END DO
  939. END DO
  940. !
  941. ! Correction : Taking into account the ellipsity of the earth orbit
  942. pqsr_ice(:,:,jl) = pqsr_ice(:,:,jl) * zcoef1 * tmask(:,:,1)
  943. !
  944. ! !--------------------------------!
  945. END DO ! end loop over ice categories !
  946. ! !--------------------------------!
  947. !!gm : this should be suppress as input data have been passed through lbc_lnk
  948. DO jl = 1, ijpl
  949. CALL lbc_lnk( pqsr_ice(:,:,jl) , 'T', 1. )
  950. END DO
  951. !
  952. CALL wrk_dealloc( jpi,jpj, zev, zdlha, zlsrise, zlsset, zps, zpc )
  953. !
  954. IF( nn_timing == 1 ) CALL timing_stop('blk_clio_qsr_ice')
  955. !
  956. END SUBROUTINE blk_clio_qsr_ice
  957. SUBROUTINE flx_blk_declin( ky, kday, pdecl )
  958. !!---------------------------------------------------------------------------
  959. !! *** ROUTINE flx_blk_declin ***
  960. !!
  961. !! ** Purpose : Computation of the solar declination for the day
  962. !!
  963. !! ** Method : ???
  964. !!---------------------------------------------------------------------
  965. INTEGER , INTENT(in ) :: ky ! = -1, 0, 1 for odd, normal and leap years resp.
  966. INTEGER , INTENT(in ) :: kday ! day of the year ( kday = 1 on january 1)
  967. REAL(wp), INTENT( out) :: pdecl ! solar declination
  968. !!
  969. REAL(wp) :: a0 = 0.39507671 ! coefficients for solar declinaison computation
  970. REAL(wp) :: a1 = 22.85684301 ! " "" "
  971. REAL(wp) :: a2 = -0.38637317 ! " "" "
  972. REAL(wp) :: a3 = 0.15096535 ! " "" "
  973. REAL(wp) :: a4 = -0.00961411 ! " "" "
  974. REAL(wp) :: b1 = -4.29692073 ! " "" "
  975. REAL(wp) :: b2 = 0.05702074 ! " "" "
  976. REAL(wp) :: b3 = -0.09028607 ! " "" "
  977. REAL(wp) :: b4 = 0.00592797
  978. !!
  979. REAL(wp) :: zday ! corresponding day of type year (cf. ky)
  980. REAL(wp) :: zp ! temporary scalars
  981. !!---------------------------------------------------------------------
  982. IF ( ky == 1 ) THEN ; zday = REAL( kday, wp ) - 0.5
  983. ELSEIF( ky == 3 ) THEN ; zday = REAL( kday, wp ) - 1.
  984. ELSE ; zday = REAL( kday, wp )
  985. ENDIF
  986. zp = rpi * ( 2.0 * zday - 367.0 ) / REAL(nyear_len(1), wp)
  987. pdecl = a0 &
  988. & + a1 * COS( zp ) + a2 * COS( 2. * zp ) + a3 * COS( 3. * zp ) + a4 * COS( 4. * zp ) &
  989. & + b1 * SIN( zp ) + b2 * SIN( 2. * zp ) + b3 * SIN( 3. * zp ) + b4 * SIN( 4. * zp )
  990. !
  991. END SUBROUTINE flx_blk_declin
  992. !!======================================================================
  993. END MODULE sbcblk_clio