[888] | 1 | MODULE sbcblk_core |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE sbcblk_core *** |
---|
| 4 | !! Ocean forcing: momentum, heat and freshwater flux formulation |
---|
| 5 | !!===================================================================== |
---|
[1482] | 6 | !! History : 1.0 ! 2004-08 (U. Schweckendiek) Original code |
---|
[4990] | 7 | !! 2.0 ! 2005-04 (L. Brodeau, A.M. Treguier) additions: |
---|
[1482] | 8 | !! - new bulk routine for efficiency |
---|
| 9 | !! - WINDS ARE NOW ASSUMED TO BE AT T POINTS in input files !!!! |
---|
[4990] | 10 | !! - file names and file characteristics in namelist |
---|
| 11 | !! - Implement reading of 6-hourly fields |
---|
| 12 | !! 3.0 ! 2006-06 (G. Madec) sbc rewritting |
---|
| 13 | !! - ! 2006-12 (L. Brodeau) Original code for turb_core_2z |
---|
[1482] | 14 | !! 3.2 ! 2009-04 (B. Lemaire) Introduce iom_put |
---|
[2528] | 15 | !! 3.3 ! 2010-10 (S. Masson) add diurnal cycle |
---|
[3294] | 16 | !! 3.4 ! 2011-11 (C. Harris) Fill arrays required by CICE |
---|
[4990] | 17 | !! 3.7 ! 2014-06 (L. Brodeau) simplification and optimization of CORE bulk |
---|
[888] | 18 | !!---------------------------------------------------------------------- |
---|
| 19 | |
---|
| 20 | !!---------------------------------------------------------------------- |
---|
[4990] | 21 | !! sbc_blk_core : bulk formulation as ocean surface boundary condition (forced mode, CORE bulk formulea) |
---|
| 22 | !! blk_oce_core : computes momentum, heat and freshwater fluxes over ocean |
---|
| 23 | !! blk_ice_core : computes momentum, heat and freshwater fluxes over ice |
---|
| 24 | !! blk_bio_meanqsr : compute daily mean short wave radiation over the ocean |
---|
| 25 | !! blk_ice_meanqsr : compute daily mean short wave radiation over the ice |
---|
| 26 | !! turb_core_2z : Computes turbulent transfert coefficients |
---|
| 27 | !! cd_neutral_10m : Estimate of the neutral drag coefficient at 10m |
---|
| 28 | !! psi_m : universal profile stability function for momentum |
---|
| 29 | !! psi_h : universal profile stability function for temperature and humidity |
---|
[888] | 30 | !!---------------------------------------------------------------------- |
---|
| 31 | USE oce ! ocean dynamics and tracers |
---|
| 32 | USE dom_oce ! ocean space and time domain |
---|
| 33 | USE phycst ! physical constants |
---|
| 34 | USE fldread ! read input fields |
---|
| 35 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
[3680] | 36 | USE cyclone ! Cyclone 10m wind form trac of cyclone centres |
---|
[2528] | 37 | USE sbcdcy ! surface boundary condition: diurnal cycle |
---|
[888] | 38 | USE iom ! I/O manager library |
---|
| 39 | USE in_out_manager ! I/O manager |
---|
| 40 | USE lib_mpp ! distribued memory computing library |
---|
[3294] | 41 | USE wrk_nemo ! work arrays |
---|
| 42 | USE timing ! Timing |
---|
[888] | 43 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 44 | USE prtctl ! Print control |
---|
[4990] | 45 | USE sbcwave, ONLY : cdn_wave ! wave module |
---|
[1465] | 46 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
[4161] | 47 | USE lib_fortran ! to use key_nosignedzero |
---|
[888] | 48 | |
---|
| 49 | IMPLICIT NONE |
---|
| 50 | PRIVATE |
---|
| 51 | |
---|
[2715] | 52 | PUBLIC sbc_blk_core ! routine called in sbcmod module |
---|
| 53 | PUBLIC blk_ice_core ! routine called in sbc_ice_lim module |
---|
[4306] | 54 | PUBLIC blk_ice_meanqsr ! routine called in sbc_ice_lim module |
---|
[3294] | 55 | PUBLIC turb_core_2z ! routine calles in sbcblk_mfs module |
---|
[2715] | 56 | |
---|
[4990] | 57 | INTEGER , PARAMETER :: jpfld = 9 ! maximum number of files to read |
---|
[888] | 58 | INTEGER , PARAMETER :: jp_wndi = 1 ! index of 10m wind velocity (i-component) (m/s) at T-point |
---|
| 59 | INTEGER , PARAMETER :: jp_wndj = 2 ! index of 10m wind velocity (j-component) (m/s) at T-point |
---|
[3625] | 60 | INTEGER , PARAMETER :: jp_humi = 3 ! index of specific humidity ( % ) |
---|
[888] | 61 | INTEGER , PARAMETER :: jp_qsr = 4 ! index of solar heat (W/m2) |
---|
| 62 | INTEGER , PARAMETER :: jp_qlw = 5 ! index of Long wave (W/m2) |
---|
| 63 | INTEGER , PARAMETER :: jp_tair = 6 ! index of 10m air temperature (Kelvin) |
---|
| 64 | INTEGER , PARAMETER :: jp_prec = 7 ! index of total precipitation (rain+snow) (Kg/m2/s) |
---|
| 65 | INTEGER , PARAMETER :: jp_snow = 8 ! index of snow (solid prcipitation) (kg/m2/s) |
---|
[1705] | 66 | INTEGER , PARAMETER :: jp_tdif = 9 ! index of tau diff associated to HF tau (N/m2) at T-point |
---|
[4990] | 67 | |
---|
[888] | 68 | TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf ! structure of input fields (file informations, fields read) |
---|
[4990] | 69 | |
---|
[1601] | 70 | ! !!! CORE bulk parameters |
---|
[888] | 71 | REAL(wp), PARAMETER :: rhoa = 1.22 ! air density |
---|
| 72 | REAL(wp), PARAMETER :: cpa = 1000.5 ! specific heat of air |
---|
| 73 | REAL(wp), PARAMETER :: Lv = 2.5e6 ! latent heat of vaporization |
---|
| 74 | REAL(wp), PARAMETER :: Ls = 2.839e6 ! latent heat of sublimation |
---|
| 75 | REAL(wp), PARAMETER :: Stef = 5.67e-8 ! Stefan Boltzmann constant |
---|
[4161] | 76 | REAL(wp), PARAMETER :: Cice = 1.4e-3 ! iovi 1.63e-3 ! transfer coefficient over ice |
---|
[3625] | 77 | REAL(wp), PARAMETER :: albo = 0.066 ! ocean albedo assumed to be constant |
---|
[888] | 78 | |
---|
[2528] | 79 | ! !!* Namelist namsbc_core : CORE bulk parameters |
---|
[4147] | 80 | LOGICAL :: ln_taudif ! logical flag to use the "mean of stress module - module of mean stress" data |
---|
| 81 | REAL(wp) :: rn_pfac ! multiplication factor for precipitation |
---|
[4161] | 82 | REAL(wp) :: rn_efac ! multiplication factor for evaporation (clem) |
---|
| 83 | REAL(wp) :: rn_vfac ! multiplication factor for ice/ocean velocity in the calculation of wind stress (clem) |
---|
[4245] | 84 | REAL(wp) :: rn_zqt ! z(q,t) : height of humidity and temperature measurements |
---|
| 85 | REAL(wp) :: rn_zu ! z(u) : height of wind measurements |
---|
[888] | 86 | |
---|
| 87 | !! * Substitutions |
---|
| 88 | # include "domzgr_substitute.h90" |
---|
| 89 | # include "vectopt_loop_substitute.h90" |
---|
| 90 | !!---------------------------------------------------------------------- |
---|
[4990] | 91 | !! NEMO/OPA 3.7 , NEMO-consortium (2014) |
---|
[1156] | 92 | !! $Id$ |
---|
[2528] | 93 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[888] | 94 | !!---------------------------------------------------------------------- |
---|
| 95 | CONTAINS |
---|
| 96 | |
---|
| 97 | SUBROUTINE sbc_blk_core( kt ) |
---|
| 98 | !!--------------------------------------------------------------------- |
---|
| 99 | !! *** ROUTINE sbc_blk_core *** |
---|
[4990] | 100 | !! |
---|
[888] | 101 | !! ** Purpose : provide at each time step the surface ocean fluxes |
---|
[4990] | 102 | !! (momentum, heat, freshwater and runoff) |
---|
[888] | 103 | !! |
---|
[1695] | 104 | !! ** Method : (1) READ each fluxes in NetCDF files: |
---|
| 105 | !! the 10m wind velocity (i-component) (m/s) at T-point |
---|
| 106 | !! the 10m wind velocity (j-component) (m/s) at T-point |
---|
[3625] | 107 | !! the 10m or 2m specific humidity ( % ) |
---|
[1695] | 108 | !! the solar heat (W/m2) |
---|
| 109 | !! the Long wave (W/m2) |
---|
[3625] | 110 | !! the 10m or 2m air temperature (Kelvin) |
---|
[1695] | 111 | !! the total precipitation (rain+snow) (Kg/m2/s) |
---|
| 112 | !! the snow (solid prcipitation) (kg/m2/s) |
---|
[3625] | 113 | !! the tau diff associated to HF tau (N/m2) at T-point (ln_taudif=T) |
---|
[1695] | 114 | !! (2) CALL blk_oce_core |
---|
[888] | 115 | !! |
---|
| 116 | !! C A U T I O N : never mask the surface stress fields |
---|
[3625] | 117 | !! the stress is assumed to be in the (i,j) mesh referential |
---|
[888] | 118 | !! |
---|
| 119 | !! ** Action : defined at each time-step at the air-sea interface |
---|
[1695] | 120 | !! - utau, vtau i- and j-component of the wind stress |
---|
[4990] | 121 | !! - taum wind stress module at T-point |
---|
| 122 | !! - wndm wind speed module at T-point over free ocean or leads in presence of sea-ice |
---|
[3625] | 123 | !! - qns, qsr non-solar and solar heat fluxes |
---|
| 124 | !! - emp upward mass flux (evapo. - precip.) |
---|
| 125 | !! - sfx salt flux due to freezing/melting (non-zero only if ice is present) |
---|
| 126 | !! (set in limsbc(_2).F90) |
---|
[4990] | 127 | !! |
---|
| 128 | !! ** References : Large & Yeager, 2004 / Large & Yeager, 2008 |
---|
| 129 | !! Brodeau et al. Ocean Modelling 2010 |
---|
[888] | 130 | !!---------------------------------------------------------------------- |
---|
[2715] | 131 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
[4990] | 132 | ! |
---|
[888] | 133 | INTEGER :: ierror ! return error code |
---|
[1200] | 134 | INTEGER :: ifpr ! dummy loop indice |
---|
[1705] | 135 | INTEGER :: jfld ! dummy loop arguments |
---|
[4147] | 136 | INTEGER :: ios ! Local integer output status for namelist read |
---|
[4990] | 137 | ! |
---|
[888] | 138 | CHARACTER(len=100) :: cn_dir ! Root directory for location of core files |
---|
| 139 | TYPE(FLD_N), DIMENSION(jpfld) :: slf_i ! array of namelist informations on the fields to read |
---|
[4161] | 140 | TYPE(FLD_N) :: sn_wndi, sn_wndj, sn_humi, sn_qsr ! informations about the fields to be read |
---|
| 141 | TYPE(FLD_N) :: sn_qlw , sn_tair, sn_prec, sn_snow ! " " |
---|
| 142 | TYPE(FLD_N) :: sn_tdif ! " " |
---|
[4990] | 143 | NAMELIST/namsbc_core/ cn_dir , ln_taudif, rn_pfac, rn_efac, rn_vfac, & |
---|
[1705] | 144 | & sn_wndi, sn_wndj, sn_humi , sn_qsr , & |
---|
[4245] | 145 | & sn_qlw , sn_tair, sn_prec , sn_snow, & |
---|
[4990] | 146 | & sn_tdif, rn_zqt, rn_zu |
---|
[888] | 147 | !!--------------------------------------------------------------------- |
---|
[4990] | 148 | ! |
---|
[888] | 149 | ! ! ====================== ! |
---|
| 150 | IF( kt == nit000 ) THEN ! First call kt=nit000 ! |
---|
| 151 | ! ! ====================== ! |
---|
[1601] | 152 | ! |
---|
[4147] | 153 | REWIND( numnam_ref ) ! Namelist namsbc_core in reference namelist : CORE bulk parameters |
---|
| 154 | READ ( numnam_ref, namsbc_core, IOSTAT = ios, ERR = 901) |
---|
| 155 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in reference namelist', lwp ) |
---|
[4990] | 156 | ! |
---|
[4147] | 157 | REWIND( numnam_cfg ) ! Namelist namsbc_core in configuration namelist : CORE bulk parameters |
---|
| 158 | READ ( numnam_cfg, namsbc_core, IOSTAT = ios, ERR = 902 ) |
---|
| 159 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in configuration namelist', lwp ) |
---|
| 160 | |
---|
[4990] | 161 | IF(lwm) WRITE( numond, namsbc_core ) |
---|
[2528] | 162 | ! ! check: do we plan to use ln_dm2dc with non-daily forcing? |
---|
[4990] | 163 | IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 ) & |
---|
| 164 | & CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' ) |
---|
[2528] | 165 | IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN |
---|
| 166 | CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr', & |
---|
[4990] | 167 | & ' ==> We force time interpolation = .false. for qsr' ) |
---|
[2528] | 168 | sn_qsr%ln_tint = .false. |
---|
| 169 | ENDIF |
---|
| 170 | ! ! store namelist information in an array |
---|
[888] | 171 | slf_i(jp_wndi) = sn_wndi ; slf_i(jp_wndj) = sn_wndj |
---|
| 172 | slf_i(jp_qsr ) = sn_qsr ; slf_i(jp_qlw ) = sn_qlw |
---|
| 173 | slf_i(jp_tair) = sn_tair ; slf_i(jp_humi) = sn_humi |
---|
| 174 | slf_i(jp_prec) = sn_prec ; slf_i(jp_snow) = sn_snow |
---|
[1705] | 175 | slf_i(jp_tdif) = sn_tdif |
---|
[4990] | 176 | ! |
---|
[2528] | 177 | lhftau = ln_taudif ! do we use HF tau information? |
---|
[1705] | 178 | jfld = jpfld - COUNT( (/.NOT. lhftau/) ) |
---|
| 179 | ! |
---|
[2528] | 180 | ALLOCATE( sf(jfld), STAT=ierror ) ! set sf structure |
---|
[2715] | 181 | IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' ) |
---|
[1705] | 182 | DO ifpr= 1, jfld |
---|
[2528] | 183 | ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) ) |
---|
[2715] | 184 | IF( slf_i(ifpr)%ln_tint ) ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) ) |
---|
[1200] | 185 | END DO |
---|
[2528] | 186 | ! ! fill sf with slf_i and control print |
---|
| 187 | CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' ) |
---|
[1601] | 188 | ! |
---|
[3625] | 189 | sfx(:,:) = 0._wp ! salt flux; zero unless ice is present (computed in limsbc(_2).F90) |
---|
| 190 | ! |
---|
[888] | 191 | ENDIF |
---|
| 192 | |
---|
[3625] | 193 | CALL fld_read( kt, nn_fsbc, sf ) ! input fields provided at the current time-step |
---|
[888] | 194 | |
---|
[3625] | 195 | ! ! compute the surface ocean fluxes using CORE bulk formulea |
---|
[3680] | 196 | IF( MOD( kt - 1, nn_fsbc ) == 0 ) CALL blk_oce_core( kt, sf, sst_m, ssu_m, ssv_m ) |
---|
[3294] | 197 | |
---|
[4990] | 198 | ! If diurnal cycle is activated, compute a daily mean short waves flux for biogeochemistery |
---|
[4306] | 199 | IF( ltrcdm2dc ) CALL blk_bio_meanqsr |
---|
| 200 | |
---|
[3294] | 201 | #if defined key_cice |
---|
| 202 | IF( MOD( kt - 1, nn_fsbc ) == 0 ) THEN |
---|
| 203 | qlw_ice(:,:,1) = sf(jp_qlw)%fnow(:,:,1) |
---|
| 204 | qsr_ice(:,:,1) = sf(jp_qsr)%fnow(:,:,1) |
---|
| 205 | tatm_ice(:,:) = sf(jp_tair)%fnow(:,:,1) |
---|
| 206 | qatm_ice(:,:) = sf(jp_humi)%fnow(:,:,1) |
---|
| 207 | tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac |
---|
| 208 | sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac |
---|
| 209 | wndi_ice(:,:) = sf(jp_wndi)%fnow(:,:,1) |
---|
| 210 | wndj_ice(:,:) = sf(jp_wndj)%fnow(:,:,1) |
---|
| 211 | ENDIF |
---|
| 212 | #endif |
---|
[2528] | 213 | ! |
---|
[888] | 214 | END SUBROUTINE sbc_blk_core |
---|
| 215 | |
---|
| 216 | |
---|
[3680] | 217 | SUBROUTINE blk_oce_core( kt, sf, pst, pu, pv ) |
---|
[888] | 218 | !!--------------------------------------------------------------------- |
---|
| 219 | !! *** ROUTINE blk_core *** |
---|
| 220 | !! |
---|
| 221 | !! ** Purpose : provide the momentum, heat and freshwater fluxes at |
---|
| 222 | !! the ocean surface at each time step |
---|
| 223 | !! |
---|
| 224 | !! ** Method : CORE bulk formulea for the ocean using atmospheric |
---|
| 225 | !! fields read in sbc_read |
---|
| 226 | !! |
---|
| 227 | !! ** Outputs : - utau : i-component of the stress at U-point (N/m2) |
---|
| 228 | !! - vtau : j-component of the stress at V-point (N/m2) |
---|
[1695] | 229 | !! - taum : Wind stress module at T-point (N/m2) |
---|
| 230 | !! - wndm : Wind speed module at T-point (m/s) |
---|
[888] | 231 | !! - qsr : Solar heat flux over the ocean (W/m2) |
---|
| 232 | !! - qns : Non Solar heat flux over the ocean (W/m2) |
---|
[3625] | 233 | !! - emp : evaporation minus precipitation (kg/m2/s) |
---|
[1242] | 234 | !! |
---|
| 235 | !! ** Nota : sf has to be a dummy argument for AGRIF on NEC |
---|
[888] | 236 | !!--------------------------------------------------------------------- |
---|
[3680] | 237 | INTEGER , INTENT(in ) :: kt ! time step index |
---|
| 238 | TYPE(fld), INTENT(inout), DIMENSION(:) :: sf ! input data |
---|
| 239 | REAL(wp) , INTENT(in) , DIMENSION(:,:) :: pst ! surface temperature [Celcius] |
---|
| 240 | REAL(wp) , INTENT(in) , DIMENSION(:,:) :: pu ! surface current at U-point (i-component) [m/s] |
---|
| 241 | REAL(wp) , INTENT(in) , DIMENSION(:,:) :: pv ! surface current at V-point (j-component) [m/s] |
---|
[2715] | 242 | ! |
---|
| 243 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 244 | REAL(wp) :: zcoef_qsatw, zztmp ! local variable |
---|
[3294] | 245 | REAL(wp), DIMENSION(:,:), POINTER :: zwnd_i, zwnd_j ! wind speed components at T-point |
---|
| 246 | REAL(wp), DIMENSION(:,:), POINTER :: zqsatw ! specific humidity at pst |
---|
| 247 | REAL(wp), DIMENSION(:,:), POINTER :: zqlw, zqsb ! long wave and sensible heat fluxes |
---|
| 248 | REAL(wp), DIMENSION(:,:), POINTER :: zqla, zevap ! latent heat fluxes and evaporation |
---|
| 249 | REAL(wp), DIMENSION(:,:), POINTER :: Cd ! transfer coefficient for momentum (tau) |
---|
| 250 | REAL(wp), DIMENSION(:,:), POINTER :: Ch ! transfer coefficient for sensible heat (Q_sens) |
---|
| 251 | REAL(wp), DIMENSION(:,:), POINTER :: Ce ! tansfert coefficient for evaporation (Q_lat) |
---|
| 252 | REAL(wp), DIMENSION(:,:), POINTER :: zst ! surface temperature in Kelvin |
---|
| 253 | REAL(wp), DIMENSION(:,:), POINTER :: zt_zu ! air temperature at wind speed height |
---|
| 254 | REAL(wp), DIMENSION(:,:), POINTER :: zq_zu ! air spec. hum. at wind speed height |
---|
[888] | 255 | !!--------------------------------------------------------------------- |
---|
[2715] | 256 | ! |
---|
[3294] | 257 | IF( nn_timing == 1 ) CALL timing_start('blk_oce_core') |
---|
| 258 | ! |
---|
| 259 | CALL wrk_alloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap ) |
---|
| 260 | CALL wrk_alloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu ) |
---|
| 261 | ! |
---|
[888] | 262 | ! local scalars ( place there for vector optimisation purposes) |
---|
| 263 | zcoef_qsatw = 0.98 * 640380. / rhoa |
---|
| 264 | |
---|
[3625] | 265 | zst(:,:) = pst(:,:) + rt0 ! convert SST from Celcius to Kelvin (and set minimum value far above 0 K) |
---|
[888] | 266 | |
---|
| 267 | ! ----------------------------------------------------------------------------- ! |
---|
| 268 | ! 0 Wind components and module at T-point relative to the moving ocean ! |
---|
| 269 | ! ----------------------------------------------------------------------------- ! |
---|
[1000] | 270 | |
---|
[888] | 271 | ! ... components ( U10m - U_oce ) at T-point (unmasked) |
---|
| 272 | zwnd_i(:,:) = 0.e0 |
---|
| 273 | zwnd_j(:,:) = 0.e0 |
---|
[3680] | 274 | #if defined key_cyclone |
---|
[4990] | 275 | CALL wnd_cyc( kt, zwnd_i, zwnd_j ) ! add analytical tropical cyclone (Vincent et al. JGR 2012) |
---|
[3680] | 276 | DO jj = 2, jpjm1 |
---|
| 277 | DO ji = fs_2, fs_jpim1 ! vect. opt. |
---|
| 278 | sf(jp_wndi)%fnow(ji,jj,1) = sf(jp_wndi)%fnow(ji,jj,1) + zwnd_i(ji,jj) |
---|
| 279 | sf(jp_wndj)%fnow(ji,jj,1) = sf(jp_wndj)%fnow(ji,jj,1) + zwnd_j(ji,jj) |
---|
| 280 | END DO |
---|
| 281 | END DO |
---|
| 282 | #endif |
---|
[888] | 283 | DO jj = 2, jpjm1 |
---|
| 284 | DO ji = fs_2, fs_jpim1 ! vect. opt. |
---|
[4161] | 285 | zwnd_i(ji,jj) = ( sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pu(ji-1,jj ) + pu(ji,jj) ) ) |
---|
| 286 | zwnd_j(ji,jj) = ( sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pv(ji ,jj-1) + pv(ji,jj) ) ) |
---|
[888] | 287 | END DO |
---|
| 288 | END DO |
---|
| 289 | CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. ) |
---|
| 290 | CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. ) |
---|
| 291 | ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked) |
---|
[1025] | 292 | wndm(:,:) = SQRT( zwnd_i(:,:) * zwnd_i(:,:) & |
---|
| 293 | & + zwnd_j(:,:) * zwnd_j(:,:) ) * tmask(:,:,1) |
---|
[888] | 294 | |
---|
| 295 | ! ----------------------------------------------------------------------------- ! |
---|
| 296 | ! I Radiative FLUXES ! |
---|
| 297 | ! ----------------------------------------------------------------------------- ! |
---|
[4990] | 298 | |
---|
[2528] | 299 | ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle ! Short Wave |
---|
| 300 | zztmp = 1. - albo |
---|
| 301 | IF( ln_dm2dc ) THEN ; qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1) |
---|
| 302 | ELSE ; qsr(:,:) = zztmp * sf(jp_qsr)%fnow(:,:,1) * tmask(:,:,1) |
---|
| 303 | ENDIF |
---|
| 304 | zqlw(:,:) = ( sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:) ) * tmask(:,:,1) ! Long Wave |
---|
[888] | 305 | ! ----------------------------------------------------------------------------- ! |
---|
| 306 | ! II Turbulent FLUXES ! |
---|
| 307 | ! ----------------------------------------------------------------------------- ! |
---|
| 308 | |
---|
| 309 | ! ... specific humidity at SST and IST |
---|
[4990] | 310 | zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) |
---|
[888] | 311 | |
---|
| 312 | ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point : |
---|
[4990] | 313 | CALL turb_core_2z( rn_zqt, rn_zu, zst, sf(jp_tair)%fnow, zqsatw, sf(jp_humi)%fnow, wndm, & |
---|
| 314 | & Cd, Ch, Ce, zt_zu, zq_zu ) |
---|
| 315 | |
---|
[1695] | 316 | ! ... tau module, i and j component |
---|
| 317 | DO jj = 1, jpj |
---|
| 318 | DO ji = 1, jpi |
---|
| 319 | zztmp = rhoa * wndm(ji,jj) * Cd(ji,jj) |
---|
| 320 | taum (ji,jj) = zztmp * wndm (ji,jj) |
---|
| 321 | zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj) |
---|
| 322 | zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj) |
---|
| 323 | END DO |
---|
| 324 | END DO |
---|
[1705] | 325 | |
---|
| 326 | ! ... add the HF tau contribution to the wind stress module? |
---|
[4990] | 327 | IF( lhftau ) THEN |
---|
[2528] | 328 | taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1) |
---|
[1705] | 329 | ENDIF |
---|
| 330 | CALL iom_put( "taum_oce", taum ) ! output wind stress module |
---|
| 331 | |
---|
[888] | 332 | ! ... utau, vtau at U- and V_points, resp. |
---|
| 333 | ! Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines |
---|
[4990] | 334 | ! Note the use of MAX(tmask(i,j),tmask(i+1,j) is to mask tau over ice shelves |
---|
[888] | 335 | DO jj = 1, jpjm1 |
---|
| 336 | DO ji = 1, fs_jpim1 |
---|
[4990] | 337 | utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj ) ) & |
---|
| 338 | & * MAX(tmask(ji,jj,1),tmask(ji+1,jj,1)) |
---|
| 339 | vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji ,jj+1) ) & |
---|
| 340 | & * MAX(tmask(ji,jj,1),tmask(ji,jj+1,1)) |
---|
[888] | 341 | END DO |
---|
| 342 | END DO |
---|
| 343 | CALL lbc_lnk( utau(:,:), 'U', -1. ) |
---|
| 344 | CALL lbc_lnk( vtau(:,:), 'V', -1. ) |
---|
| 345 | |
---|
[4990] | 346 | |
---|
[888] | 347 | ! Turbulent fluxes over ocean |
---|
| 348 | ! ----------------------------- |
---|
[4990] | 349 | IF( ABS( rn_zu - rn_zqt) < 0.01_wp ) THEN |
---|
| 350 | !! q_air and t_air are (or "are almost") given at 10m (wind reference height) |
---|
| 351 | zevap(:,:) = rn_efac*MAX( 0._wp, rhoa*Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:,1) )*wndm(:,:) ) ! Evaporation |
---|
| 352 | zqsb (:,:) = cpa*rhoa*Ch(:,:)*( zst (:,:) - sf(jp_tair)%fnow(:,:,1) )*wndm(:,:) ! Sensible Heat |
---|
[888] | 353 | ELSE |
---|
[4990] | 354 | !! q_air and t_air are not given at 10m (wind reference height) |
---|
| 355 | ! Values of temp. and hum. adjusted to height of wind during bulk algorithm iteration must be used!!! |
---|
| 356 | zevap(:,:) = rn_efac*MAX( 0._wp, rhoa*Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) )*wndm(:,:) ) ! Evaporation |
---|
| 357 | zqsb (:,:) = cpa*rhoa*Ch(:,:)*( zst (:,:) - zt_zu(:,:) )*wndm(:,:) ! Sensible Heat |
---|
[888] | 358 | ENDIF |
---|
| 359 | zqla (:,:) = Lv * zevap(:,:) ! Latent Heat |
---|
| 360 | |
---|
| 361 | IF(ln_ctl) THEN |
---|
[1025] | 362 | CALL prt_ctl( tab2d_1=zqla , clinfo1=' blk_oce_core: zqla : ', tab2d_2=Ce , clinfo2=' Ce : ' ) |
---|
| 363 | CALL prt_ctl( tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb : ', tab2d_2=Ch , clinfo2=' Ch : ' ) |
---|
| 364 | CALL prt_ctl( tab2d_1=zqlw , clinfo1=' blk_oce_core: zqlw : ', tab2d_2=qsr, clinfo2=' qsr : ' ) |
---|
| 365 | CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' ) |
---|
| 366 | CALL prt_ctl( tab2d_1=utau , clinfo1=' blk_oce_core: utau : ', mask1=umask, & |
---|
| 367 | & tab2d_2=vtau , clinfo2= ' vtau : ' , mask2=vmask ) |
---|
| 368 | CALL prt_ctl( tab2d_1=wndm , clinfo1=' blk_oce_core: wndm : ') |
---|
| 369 | CALL prt_ctl( tab2d_1=zst , clinfo1=' blk_oce_core: zst : ') |
---|
[888] | 370 | ENDIF |
---|
| 371 | |
---|
| 372 | ! ----------------------------------------------------------------------------- ! |
---|
| 373 | ! III Total FLUXES ! |
---|
| 374 | ! ----------------------------------------------------------------------------- ! |
---|
[4990] | 375 | ! |
---|
[3625] | 376 | emp (:,:) = ( zevap(:,:) & ! mass flux (evap. - precip.) |
---|
| 377 | & - sf(jp_prec)%fnow(:,:,1) * rn_pfac ) * tmask(:,:,1) |
---|
[3772] | 378 | qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:) & ! Downward Non Solar flux |
---|
| 379 | & - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus & ! remove latent melting heat for solid precip |
---|
| 380 | & - zevap(:,:) * pst(:,:) * rcp & ! remove evap heat content at SST |
---|
| 381 | & + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac & ! add liquid precip heat content at Tair |
---|
[4990] | 382 | & * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp & |
---|
[3772] | 383 | & + sf(jp_snow)%fnow(:,:,1) * rn_pfac & ! add solid precip heat content at min(Tair,Tsnow) |
---|
[4990] | 384 | & * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) |
---|
[888] | 385 | ! |
---|
[1482] | 386 | CALL iom_put( "qlw_oce", zqlw ) ! output downward longwave heat over the ocean |
---|
| 387 | CALL iom_put( "qsb_oce", - zqsb ) ! output downward sensible heat over the ocean |
---|
| 388 | CALL iom_put( "qla_oce", - zqla ) ! output downward latent heat over the ocean |
---|
[3625] | 389 | CALL iom_put( "qhc_oce", qns-zqlw+zqsb+zqla ) ! output downward heat content of E-P over the ocean |
---|
[1482] | 390 | CALL iom_put( "qns_oce", qns ) ! output downward non solar heat over the ocean |
---|
| 391 | ! |
---|
| 392 | IF(ln_ctl) THEN |
---|
| 393 | CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb : ', tab2d_2=zqlw , clinfo2=' zqlw : ') |
---|
| 394 | CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_core: zqla : ', tab2d_2=qsr , clinfo2=' qsr : ') |
---|
| 395 | CALL prt_ctl(tab2d_1=pst , clinfo1=' blk_oce_core: pst : ', tab2d_2=emp , clinfo2=' emp : ') |
---|
| 396 | CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_core: utau : ', mask1=umask, & |
---|
| 397 | & tab2d_2=vtau , clinfo2= ' vtau : ' , mask2=vmask ) |
---|
| 398 | ENDIF |
---|
| 399 | ! |
---|
[3294] | 400 | CALL wrk_dealloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap ) |
---|
| 401 | CALL wrk_dealloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu ) |
---|
[2715] | 402 | ! |
---|
[3294] | 403 | IF( nn_timing == 1 ) CALL timing_stop('blk_oce_core') |
---|
| 404 | ! |
---|
[888] | 405 | END SUBROUTINE blk_oce_core |
---|
[4306] | 406 | |
---|
[888] | 407 | |
---|
| 408 | SUBROUTINE blk_ice_core( pst , pui , pvi , palb , & |
---|
| 409 | & p_taui, p_tauj, p_qns , p_qsr, & |
---|
| 410 | & p_qla , p_dqns, p_dqla, & |
---|
| 411 | & p_tpr , p_spr , & |
---|
[1270] | 412 | & p_fr1 , p_fr2 , cd_grid, pdim ) |
---|
[888] | 413 | !!--------------------------------------------------------------------- |
---|
| 414 | !! *** ROUTINE blk_ice_core *** |
---|
| 415 | !! |
---|
| 416 | !! ** Purpose : provide the surface boundary condition over sea-ice |
---|
| 417 | !! |
---|
| 418 | !! ** Method : compute momentum, heat and freshwater exchanged |
---|
| 419 | !! between atmosphere and sea-ice using CORE bulk |
---|
| 420 | !! formulea, ice variables and read atmmospheric fields. |
---|
| 421 | !! NB: ice drag coefficient is assumed to be a constant |
---|
| 422 | !! |
---|
| 423 | !! caution : the net upward water flux has with mm/day unit |
---|
| 424 | !!--------------------------------------------------------------------- |
---|
[2715] | 425 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: pst ! ice surface temperature (>0, =rt0 over land) [Kelvin] |
---|
| 426 | REAL(wp), DIMENSION(:,:) , INTENT(in ) :: pui ! ice surface velocity (i- and i- components [m/s] |
---|
| 427 | REAL(wp), DIMENSION(:,:) , INTENT(in ) :: pvi ! at I-point (B-grid) or U & V-point (C-grid) |
---|
[4990] | 428 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb ! ice albedo (all skies) [%] |
---|
[2715] | 429 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_taui ! i- & j-components of surface ice stress [N/m2] |
---|
| 430 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_tauj ! at I-point (B-grid) or U & V-point (C-grid) |
---|
| 431 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_qns ! non solar heat flux over ice (T-point) [W/m2] |
---|
| 432 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_qsr ! solar heat flux over ice (T-point) [W/m2] |
---|
| 433 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_qla ! latent heat flux over ice (T-point) [W/m2] |
---|
| 434 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_dqns ! non solar heat sensistivity (T-point) [W/m2] |
---|
| 435 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_dqla ! latent heat sensistivity (T-point) [W/m2] |
---|
| 436 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_tpr ! total precipitation (T-point) [Kg/m2/s] |
---|
| 437 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_spr ! solid precipitation (T-point) [Kg/m2/s] |
---|
| 438 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_fr1 ! 1sr fraction of qsr penetration in ice (T-point) [%] |
---|
| 439 | REAL(wp), DIMENSION(:,:) , INTENT( out) :: p_fr2 ! 2nd fraction of qsr penetration in ice (T-point) [%] |
---|
| 440 | CHARACTER(len=1) , INTENT(in ) :: cd_grid ! ice grid ( C or B-grid) |
---|
| 441 | INTEGER , INTENT(in ) :: pdim ! number of ice categories |
---|
| 442 | !! |
---|
[888] | 443 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
| 444 | INTEGER :: ijpl ! number of ice categories (size of 3rd dim of input arrays) |
---|
| 445 | REAL(wp) :: zst2, zst3 |
---|
| 446 | REAL(wp) :: zcoef_wnorm, zcoef_wnorm2, zcoef_dqlw, zcoef_dqla, zcoef_dqsb |
---|
[2528] | 447 | REAL(wp) :: zztmp ! temporary variable |
---|
[888] | 448 | REAL(wp) :: zwnorm_f, zwndi_f , zwndj_f ! relative wind module and components at F-point |
---|
| 449 | REAL(wp) :: zwndi_t , zwndj_t ! relative wind components at T-point |
---|
[2715] | 450 | !! |
---|
[3294] | 451 | REAL(wp), DIMENSION(:,:) , POINTER :: z_wnds_t ! wind speed ( = | U10m - U_ice | ) at T-point |
---|
[2715] | 452 | REAL(wp), DIMENSION(:,:,:), POINTER :: z_qlw ! long wave heat flux over ice |
---|
| 453 | REAL(wp), DIMENSION(:,:,:), POINTER :: z_qsb ! sensible heat flux over ice |
---|
| 454 | REAL(wp), DIMENSION(:,:,:), POINTER :: z_dqlw ! long wave heat sensitivity over ice |
---|
| 455 | REAL(wp), DIMENSION(:,:,:), POINTER :: z_dqsb ! sensible heat sensitivity over ice |
---|
[888] | 456 | !!--------------------------------------------------------------------- |
---|
[3294] | 457 | ! |
---|
| 458 | IF( nn_timing == 1 ) CALL timing_start('blk_ice_core') |
---|
| 459 | ! |
---|
| 460 | CALL wrk_alloc( jpi,jpj, z_wnds_t ) |
---|
| 461 | CALL wrk_alloc( jpi,jpj,pdim, z_qlw, z_qsb, z_dqlw, z_dqsb ) |
---|
[888] | 462 | |
---|
[1270] | 463 | ijpl = pdim ! number of ice categories |
---|
[888] | 464 | |
---|
| 465 | ! local scalars ( place there for vector optimisation purposes) |
---|
[2528] | 466 | zcoef_wnorm = rhoa * Cice |
---|
[888] | 467 | zcoef_wnorm2 = rhoa * Cice * 0.5 |
---|
[2528] | 468 | zcoef_dqlw = 4.0 * 0.95 * Stef |
---|
| 469 | zcoef_dqla = -Ls * Cice * 11637800. * (-5897.8) |
---|
| 470 | zcoef_dqsb = rhoa * cpa * Cice |
---|
[888] | 471 | |
---|
| 472 | !!gm brutal.... |
---|
| 473 | z_wnds_t(:,:) = 0.e0 |
---|
| 474 | p_taui (:,:) = 0.e0 |
---|
| 475 | p_tauj (:,:) = 0.e0 |
---|
| 476 | !!gm end |
---|
| 477 | |
---|
[2777] | 478 | #if defined key_lim3 |
---|
| 479 | tatm_ice(:,:) = sf(jp_tair)%fnow(:,:,1) ! LIM3: make Tair available in sea-ice. WARNING allocated after call to ice_init |
---|
| 480 | #endif |
---|
[888] | 481 | ! ----------------------------------------------------------------------------- ! |
---|
| 482 | ! Wind components and module relative to the moving ocean ( U10m - U_ice ) ! |
---|
| 483 | ! ----------------------------------------------------------------------------- ! |
---|
| 484 | SELECT CASE( cd_grid ) |
---|
[2528] | 485 | CASE( 'I' ) ! B-grid ice dynamics : I-point (i.e. F-point with sea-ice indexation) |
---|
[888] | 486 | ! and scalar wind at T-point ( = | U10m - U_ice | ) (masked) |
---|
| 487 | DO jj = 2, jpjm1 |
---|
[2528] | 488 | DO ji = 2, jpim1 ! B grid : NO vector opt |
---|
[888] | 489 | ! ... scalar wind at I-point (fld being at T-point) |
---|
[2528] | 490 | zwndi_f = 0.25 * ( sf(jp_wndi)%fnow(ji-1,jj ,1) + sf(jp_wndi)%fnow(ji ,jj ,1) & |
---|
[4161] | 491 | & + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji ,jj-1,1) ) - rn_vfac * pui(ji,jj) |
---|
[2528] | 492 | zwndj_f = 0.25 * ( sf(jp_wndj)%fnow(ji-1,jj ,1) + sf(jp_wndj)%fnow(ji ,jj ,1) & |
---|
[4161] | 493 | & + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji ,jj-1,1) ) - rn_vfac * pvi(ji,jj) |
---|
[888] | 494 | zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f ) |
---|
| 495 | ! ... ice stress at I-point |
---|
| 496 | p_taui(ji,jj) = zwnorm_f * zwndi_f |
---|
| 497 | p_tauj(ji,jj) = zwnorm_f * zwndj_f |
---|
| 498 | ! ... scalar wind at T-point (fld being at T-point) |
---|
[4161] | 499 | zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.25 * ( pui(ji,jj+1) + pui(ji+1,jj+1) & |
---|
| 500 | & + pui(ji,jj ) + pui(ji+1,jj ) ) |
---|
| 501 | zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.25 * ( pvi(ji,jj+1) + pvi(ji+1,jj+1) & |
---|
| 502 | & + pvi(ji,jj ) + pvi(ji+1,jj ) ) |
---|
[888] | 503 | z_wnds_t(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1) |
---|
| 504 | END DO |
---|
| 505 | END DO |
---|
| 506 | CALL lbc_lnk( p_taui , 'I', -1. ) |
---|
| 507 | CALL lbc_lnk( p_tauj , 'I', -1. ) |
---|
| 508 | CALL lbc_lnk( z_wnds_t, 'T', 1. ) |
---|
| 509 | ! |
---|
| 510 | CASE( 'C' ) ! C-grid ice dynamics : U & V-points (same as ocean) |
---|
| 511 | DO jj = 2, jpj |
---|
| 512 | DO ji = fs_2, jpi ! vect. opt. |
---|
[4161] | 513 | zwndi_t = ( sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pui(ji-1,jj ) + pui(ji,jj) ) ) |
---|
| 514 | zwndj_t = ( sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pvi(ji ,jj-1) + pvi(ji,jj) ) ) |
---|
[888] | 515 | z_wnds_t(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1) |
---|
| 516 | END DO |
---|
| 517 | END DO |
---|
| 518 | DO jj = 2, jpjm1 |
---|
| 519 | DO ji = fs_2, fs_jpim1 ! vect. opt. |
---|
[2528] | 520 | p_taui(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji+1,jj ) + z_wnds_t(ji,jj) ) & |
---|
[4161] | 521 | & * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - rn_vfac * pui(ji,jj) ) |
---|
[2528] | 522 | p_tauj(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji,jj+1 ) + z_wnds_t(ji,jj) ) & |
---|
[4161] | 523 | & * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - rn_vfac * pvi(ji,jj) ) |
---|
[888] | 524 | END DO |
---|
| 525 | END DO |
---|
| 526 | CALL lbc_lnk( p_taui , 'U', -1. ) |
---|
| 527 | CALL lbc_lnk( p_tauj , 'V', -1. ) |
---|
| 528 | CALL lbc_lnk( z_wnds_t, 'T', 1. ) |
---|
| 529 | ! |
---|
| 530 | END SELECT |
---|
| 531 | |
---|
[2528] | 532 | zztmp = 1. / ( 1. - albo ) |
---|
[888] | 533 | ! ! ========================== ! |
---|
| 534 | DO jl = 1, ijpl ! Loop over ice categories ! |
---|
| 535 | ! ! ========================== ! |
---|
| 536 | DO jj = 1 , jpj |
---|
| 537 | DO ji = 1, jpi |
---|
| 538 | ! ----------------------------! |
---|
| 539 | ! I Radiative FLUXES ! |
---|
| 540 | ! ----------------------------! |
---|
| 541 | zst2 = pst(ji,jj,jl) * pst(ji,jj,jl) |
---|
| 542 | zst3 = pst(ji,jj,jl) * zst2 |
---|
| 543 | ! Short Wave (sw) |
---|
[2528] | 544 | p_qsr(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj) |
---|
[888] | 545 | ! Long Wave (lw) |
---|
[4689] | 546 | z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * pst(ji,jj,jl) * zst3 ) * tmask(ji,jj,1) |
---|
[888] | 547 | ! lw sensitivity |
---|
| 548 | z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3 |
---|
| 549 | |
---|
| 550 | ! ----------------------------! |
---|
| 551 | ! II Turbulent FLUXES ! |
---|
| 552 | ! ----------------------------! |
---|
| 553 | |
---|
| 554 | ! ... turbulent heat fluxes |
---|
| 555 | ! Sensible Heat |
---|
[2528] | 556 | z_qsb(ji,jj,jl) = rhoa * cpa * Cice * z_wnds_t(ji,jj) * ( pst(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) ) |
---|
[888] | 557 | ! Latent Heat |
---|
[4161] | 558 | p_qla(ji,jj,jl) = rn_efac * MAX( 0.e0, rhoa * Ls * Cice * z_wnds_t(ji,jj) & |
---|
| 559 | & * ( 11637800. * EXP( -5897.8 / pst(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1) ) ) |
---|
[888] | 560 | ! Latent heat sensitivity for ice (Dqla/Dt) |
---|
[4990] | 561 | IF( p_qla(ji,jj,jl) > 0._wp ) THEN |
---|
[4689] | 562 | p_dqla(ji,jj,jl) = rn_efac * zcoef_dqla * z_wnds_t(ji,jj) / ( zst2 ) * EXP( -5897.8 / pst(ji,jj,jl) ) |
---|
| 563 | ELSE |
---|
[4990] | 564 | p_dqla(ji,jj,jl) = 0._wp |
---|
[4689] | 565 | ENDIF |
---|
[4990] | 566 | |
---|
[888] | 567 | ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice) |
---|
| 568 | z_dqsb(ji,jj,jl) = zcoef_dqsb * z_wnds_t(ji,jj) |
---|
| 569 | |
---|
| 570 | ! ----------------------------! |
---|
| 571 | ! III Total FLUXES ! |
---|
| 572 | ! ----------------------------! |
---|
| 573 | ! Downward Non Solar flux |
---|
[4990] | 574 | p_qns (ji,jj,jl) = z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - p_qla (ji,jj,jl) |
---|
[888] | 575 | ! Total non solar heat flux sensitivity for ice |
---|
[4990] | 576 | p_dqns(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + p_dqla(ji,jj,jl) ) |
---|
[888] | 577 | END DO |
---|
| 578 | ! |
---|
| 579 | END DO |
---|
| 580 | ! |
---|
| 581 | END DO |
---|
| 582 | ! |
---|
| 583 | !-------------------------------------------------------------------- |
---|
| 584 | ! FRACTIONs of net shortwave radiation which is not absorbed in the |
---|
| 585 | ! thin surface layer and penetrates inside the ice cover |
---|
| 586 | ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 ) |
---|
[4990] | 587 | ! |
---|
| 588 | p_fr1(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice ) |
---|
| 589 | p_fr2(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice ) |
---|
| 590 | ! |
---|
[2528] | 591 | p_tpr(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac ! total precipitation [kg/m2/s] |
---|
| 592 | p_spr(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac ! solid precipitation [kg/m2/s] |
---|
[4990] | 593 | CALL iom_put( 'snowpre', p_spr * 86400. ) ! Snow precipitation |
---|
| 594 | CALL iom_put( 'precip' , p_tpr * 86400. ) ! Total precipitation |
---|
[888] | 595 | ! |
---|
| 596 | IF(ln_ctl) THEN |
---|
| 597 | CALL prt_ctl(tab3d_1=p_qla , clinfo1=' blk_ice_core: p_qla : ', tab3d_2=z_qsb , clinfo2=' z_qsb : ', kdim=ijpl) |
---|
| 598 | CALL prt_ctl(tab3d_1=z_qlw , clinfo1=' blk_ice_core: z_qlw : ', tab3d_2=p_dqla , clinfo2=' p_dqla : ', kdim=ijpl) |
---|
| 599 | CALL prt_ctl(tab3d_1=z_dqsb , clinfo1=' blk_ice_core: z_dqsb : ', tab3d_2=z_dqlw , clinfo2=' z_dqlw : ', kdim=ijpl) |
---|
| 600 | CALL prt_ctl(tab3d_1=p_dqns , clinfo1=' blk_ice_core: p_dqns : ', tab3d_2=p_qsr , clinfo2=' p_qsr : ', kdim=ijpl) |
---|
| 601 | CALL prt_ctl(tab3d_1=pst , clinfo1=' blk_ice_core: pst : ', tab3d_2=p_qns , clinfo2=' p_qns : ', kdim=ijpl) |
---|
| 602 | CALL prt_ctl(tab2d_1=p_tpr , clinfo1=' blk_ice_core: p_tpr : ', tab2d_2=p_spr , clinfo2=' p_spr : ') |
---|
| 603 | CALL prt_ctl(tab2d_1=p_taui , clinfo1=' blk_ice_core: p_taui : ', tab2d_2=p_tauj , clinfo2=' p_tauj : ') |
---|
| 604 | CALL prt_ctl(tab2d_1=z_wnds_t, clinfo1=' blk_ice_core: z_wnds_t : ') |
---|
| 605 | ENDIF |
---|
| 606 | |
---|
[4990] | 607 | CALL wrk_dealloc( jpi,jpj, z_wnds_t ) |
---|
| 608 | CALL wrk_dealloc( jpi,jpj, pdim, z_qlw, z_qsb, z_dqlw, z_dqsb ) |
---|
[2715] | 609 | ! |
---|
[3294] | 610 | IF( nn_timing == 1 ) CALL timing_stop('blk_ice_core') |
---|
| 611 | ! |
---|
[888] | 612 | END SUBROUTINE blk_ice_core |
---|
| 613 | |
---|
[3294] | 614 | |
---|
[4990] | 615 | SUBROUTINE blk_bio_meanqsr |
---|
| 616 | !!--------------------------------------------------------------------- |
---|
| 617 | !! *** ROUTINE blk_bio_meanqsr |
---|
| 618 | !! |
---|
| 619 | !! ** Purpose : provide daily qsr_mean for PISCES when |
---|
| 620 | !! analytic diurnal cycle is applied in physic |
---|
| 621 | !! |
---|
| 622 | !! ** Method : add part where there is no ice |
---|
| 623 | !! |
---|
| 624 | !!--------------------------------------------------------------------- |
---|
| 625 | IF( nn_timing == 1 ) CALL timing_start('blk_bio_meanqsr') |
---|
[3294] | 626 | ! |
---|
[4990] | 627 | qsr_mean(:,:) = (1. - albo ) * sf(jp_qsr)%fnow(:,:,1) |
---|
[3294] | 628 | ! |
---|
[4990] | 629 | IF( nn_timing == 1 ) CALL timing_stop('blk_bio_meanqsr') |
---|
| 630 | ! |
---|
| 631 | END SUBROUTINE blk_bio_meanqsr |
---|
| 632 | |
---|
| 633 | |
---|
| 634 | SUBROUTINE blk_ice_meanqsr( palb, p_qsr_mean, pdim ) |
---|
| 635 | !!--------------------------------------------------------------------- |
---|
[888] | 636 | !! |
---|
[4990] | 637 | !! ** Purpose : provide the daily qsr_mean over sea_ice for PISCES when |
---|
| 638 | !! analytic diurnal cycle is applied in physic |
---|
[888] | 639 | !! |
---|
[4990] | 640 | !! ** Method : compute qsr |
---|
| 641 | !! |
---|
| 642 | !!--------------------------------------------------------------------- |
---|
| 643 | REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb ! ice albedo (clear sky) (alb_ice_cs) [%] |
---|
| 644 | REAL(wp), DIMENSION(:,:,:), INTENT( out) :: p_qsr_mean ! solar heat flux over ice (T-point) [W/m2] |
---|
| 645 | INTEGER , INTENT(in ) :: pdim ! number of ice categories |
---|
| 646 | ! |
---|
| 647 | INTEGER :: ijpl ! number of ice categories (size of 3rd dim of input arrays) |
---|
| 648 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
| 649 | REAL(wp) :: zztmp ! temporary variable |
---|
| 650 | !!--------------------------------------------------------------------- |
---|
| 651 | IF( nn_timing == 1 ) CALL timing_start('blk_ice_meanqsr') |
---|
| 652 | ! |
---|
| 653 | ijpl = pdim ! number of ice categories |
---|
| 654 | zztmp = 1. / ( 1. - albo ) |
---|
| 655 | ! ! ========================== ! |
---|
| 656 | DO jl = 1, ijpl ! Loop over ice categories ! |
---|
| 657 | ! ! ========================== ! |
---|
| 658 | DO jj = 1 , jpj |
---|
| 659 | DO ji = 1, jpi |
---|
| 660 | p_qsr_mean(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr_mean(ji,jj) |
---|
| 661 | END DO |
---|
| 662 | END DO |
---|
[888] | 663 | END DO |
---|
[2715] | 664 | ! |
---|
[4990] | 665 | IF( nn_timing == 1 ) CALL timing_stop('blk_ice_meanqsr') |
---|
[3294] | 666 | ! |
---|
[4990] | 667 | END SUBROUTINE blk_ice_meanqsr |
---|
[888] | 668 | |
---|
| 669 | |
---|
[4990] | 670 | SUBROUTINE turb_core_2z( zt, zu, sst, T_zt, q_sat, q_zt, dU, & |
---|
| 671 | & Cd, Ch, Ce , T_zu, q_zu ) |
---|
[888] | 672 | !!---------------------------------------------------------------------- |
---|
| 673 | !! *** ROUTINE turb_core *** |
---|
| 674 | !! |
---|
[4990] | 675 | !! ** Purpose : Computes turbulent transfert coefficients of surface |
---|
| 676 | !! fluxes according to Large & Yeager (2004) and Large & Yeager (2008) |
---|
| 677 | !! If relevant (zt /= zu), adjust temperature and humidity from height zt to zu |
---|
[888] | 678 | !! |
---|
[4990] | 679 | !! ** Method : Monin Obukhov Similarity Theory |
---|
| 680 | !! + Large & Yeager (2004,2008) closure: CD_n10 = f(U_n10) |
---|
[888] | 681 | !! |
---|
[4990] | 682 | !! ** References : Large & Yeager, 2004 / Large & Yeager, 2008 |
---|
| 683 | !! |
---|
| 684 | !! ** Last update: Laurent Brodeau, June 2014: |
---|
| 685 | !! - handles both cases zt=zu and zt/=zu |
---|
| 686 | !! - optimized: less 2D arrays allocated and less operations |
---|
| 687 | !! - better first guess of stability by checking air-sea difference of virtual temperature |
---|
| 688 | !! rather than temperature difference only... |
---|
| 689 | !! - added function "cd_neutral_10m" that uses the improved parametrization of |
---|
| 690 | !! Large & Yeager 2008. Drag-coefficient reduction for Cyclone conditions! |
---|
| 691 | !! - using code-wide physical constants defined into "phycst.mod" rather than redifining them |
---|
| 692 | !! => 'vkarmn' and 'grav' |
---|
[888] | 693 | !!---------------------------------------------------------------------- |
---|
[3294] | 694 | REAL(wp), INTENT(in ) :: zt ! height for T_zt and q_zt [m] |
---|
| 695 | REAL(wp), INTENT(in ) :: zu ! height for dU [m] |
---|
| 696 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: sst ! sea surface temperature [Kelvin] |
---|
| 697 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: T_zt ! potential air temperature [Kelvin] |
---|
| 698 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: q_sat ! sea surface specific humidity [kg/kg] |
---|
| 699 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: q_zt ! specific air humidity [kg/kg] |
---|
[4990] | 700 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: dU ! relative wind module at zu [m/s] |
---|
[3294] | 701 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: Cd ! transfer coefficient for momentum (tau) |
---|
| 702 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: Ch ! transfer coefficient for sensible heat (Q_sens) |
---|
| 703 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: Ce ! transfert coefficient for evaporation (Q_lat) |
---|
| 704 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: T_zu ! air temp. shifted at zu [K] |
---|
| 705 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: q_zu ! spec. hum. shifted at zu [kg/kg] |
---|
[4990] | 706 | ! |
---|
| 707 | INTEGER :: j_itt |
---|
| 708 | INTEGER , PARAMETER :: nb_itt = 5 ! number of itterations |
---|
| 709 | LOGICAL :: l_zt_equal_zu = .FALSE. ! if q and t are given at different height than U |
---|
| 710 | ! |
---|
| 711 | REAL(wp), DIMENSION(:,:), POINTER :: U_zu ! relative wind at zu [m/s] |
---|
[3294] | 712 | REAL(wp), DIMENSION(:,:), POINTER :: Ce_n10 ! 10m neutral latent coefficient |
---|
| 713 | REAL(wp), DIMENSION(:,:), POINTER :: Ch_n10 ! 10m neutral sensible coefficient |
---|
| 714 | REAL(wp), DIMENSION(:,:), POINTER :: sqrt_Cd_n10 ! root square of Cd_n10 |
---|
| 715 | REAL(wp), DIMENSION(:,:), POINTER :: sqrt_Cd ! root square of Cd |
---|
| 716 | REAL(wp), DIMENSION(:,:), POINTER :: zeta_u ! stability parameter at height zu |
---|
| 717 | REAL(wp), DIMENSION(:,:), POINTER :: zeta_t ! stability parameter at height zt |
---|
[4990] | 718 | REAL(wp), DIMENSION(:,:), POINTER :: zpsi_h_u, zpsi_m_u |
---|
| 719 | REAL(wp), DIMENSION(:,:), POINTER :: ztmp0, ztmp1, ztmp2 |
---|
| 720 | REAL(wp), DIMENSION(:,:), POINTER :: stab ! 1st stability test integer |
---|
[2528] | 721 | !!---------------------------------------------------------------------- |
---|
[888] | 722 | |
---|
[4990] | 723 | IF( nn_timing == 1 ) CALL timing_start('turb_core_2z') |
---|
| 724 | |
---|
| 725 | CALL wrk_alloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd ) |
---|
| 726 | CALL wrk_alloc( jpi,jpj, zeta_u, stab ) |
---|
| 727 | CALL wrk_alloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 ) |
---|
[888] | 728 | |
---|
[4990] | 729 | l_zt_equal_zu = .FALSE. |
---|
| 730 | IF( ABS(zu - zt) < 0.01 ) l_zt_equal_zu = .TRUE. ! testing "zu == zt" is risky with double precision |
---|
| 731 | |
---|
| 732 | IF( .NOT. l_zt_equal_zu ) CALL wrk_alloc( jpi,jpj, zeta_t ) |
---|
| 733 | |
---|
| 734 | U_zu = MAX( 0.5 , dU ) ! relative wind speed at zu (normally 10m), we don't want to fall under 0.5 m/s |
---|
| 735 | |
---|
| 736 | !! First guess of stability: |
---|
| 737 | ztmp0 = T_zt*(1. + 0.608*q_zt) - sst*(1. + 0.608*q_sat) ! air-sea difference of virtual pot. temp. at zt |
---|
| 738 | stab = 0.5 + sign(0.5,ztmp0) ! stab = 1 if dTv > 0 => STABLE, 0 if unstable |
---|
| 739 | |
---|
| 740 | !! Neutral coefficients at 10m: |
---|
| 741 | IF( ln_cdgw ) THEN ! wave drag case |
---|
| 742 | cdn_wave(:,:) = cdn_wave(:,:) + rsmall * ( 1._wp - tmask(:,:,1) ) |
---|
| 743 | ztmp0 (:,:) = cdn_wave(:,:) |
---|
[3294] | 744 | ELSE |
---|
[4990] | 745 | ztmp0 = cd_neutral_10m( U_zu ) |
---|
[3294] | 746 | ENDIF |
---|
[4990] | 747 | sqrt_Cd_n10 = SQRT( ztmp0 ) |
---|
[4333] | 748 | Ce_n10 = 1.e-3*( 34.6 * sqrt_Cd_n10 ) |
---|
| 749 | Ch_n10 = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab)) |
---|
[4990] | 750 | |
---|
[888] | 751 | !! Initializing transf. coeff. with their first guess neutral equivalents : |
---|
[4990] | 752 | Cd = ztmp0 ; Ce = Ce_n10 ; Ch = Ch_n10 ; sqrt_Cd = sqrt_Cd_n10 |
---|
[888] | 753 | |
---|
[4990] | 754 | !! Initializing values at z_u with z_t values: |
---|
| 755 | T_zu = T_zt ; q_zu = q_zt |
---|
[888] | 756 | |
---|
| 757 | !! * Now starting iteration loop |
---|
| 758 | DO j_itt=1, nb_itt |
---|
[4990] | 759 | ! |
---|
| 760 | ztmp1 = T_zu - sst ! Updating air/sea differences |
---|
| 761 | ztmp2 = q_zu - q_sat |
---|
| 762 | |
---|
| 763 | ! Updating turbulent scales : (L&Y 2004 eq. (7)) |
---|
| 764 | ztmp1 = Ch/sqrt_Cd*ztmp1 ! theta* |
---|
| 765 | ztmp2 = Ce/sqrt_Cd*ztmp2 ! q* |
---|
| 766 | |
---|
| 767 | ztmp0 = T_zu*(1. + 0.608*q_zu) ! virtual potential temperature at zu |
---|
| 768 | |
---|
| 769 | ! Estimate the inverse of Monin-Obukov length (1/L) at height zu: |
---|
| 770 | ztmp0 = (vkarmn*grav/ztmp0*(ztmp1*(1.+0.608*q_zu) + 0.608*T_zu*ztmp2)) / (Cd*U_zu*U_zu) |
---|
| 771 | ! ( Cd*U_zu*U_zu is U*^2 at zu) |
---|
| 772 | |
---|
[888] | 773 | !! Stability parameters : |
---|
[4990] | 774 | zeta_u = zu*ztmp0 ; zeta_u = sign( min(abs(zeta_u),10.0), zeta_u ) |
---|
| 775 | zpsi_h_u = psi_h( zeta_u ) |
---|
| 776 | zpsi_m_u = psi_m( zeta_u ) |
---|
| 777 | |
---|
| 778 | !! Shifting temperature and humidity at zu (L&Y 2004 eq. (9b-9c)) |
---|
| 779 | IF ( .NOT. l_zt_equal_zu ) THEN |
---|
| 780 | zeta_t = zt*ztmp0 ; zeta_t = sign( min(abs(zeta_t),10.0), zeta_t ) |
---|
| 781 | stab = LOG(zu/zt) - zpsi_h_u + psi_h(zeta_t) ! stab just used as temp array!!! |
---|
| 782 | T_zu = T_zt + ztmp1/vkarmn*stab ! ztmp1 is still theta* |
---|
| 783 | q_zu = q_zt + ztmp2/vkarmn*stab ! ztmp2 is still q* |
---|
| 784 | q_zu = max(0., q_zu) |
---|
| 785 | END IF |
---|
| 786 | |
---|
| 787 | IF( ln_cdgw ) THEN ! surface wave case |
---|
| 788 | sqrt_Cd = vkarmn / ( vkarmn / sqrt_Cd_n10 - zpsi_m_u ) |
---|
| 789 | Cd = sqrt_Cd * sqrt_Cd |
---|
[3294] | 790 | ELSE |
---|
[4990] | 791 | ! Update neutral wind speed at 10m and neutral Cd at 10m (L&Y 2004 eq. 9a)... |
---|
| 792 | ! In very rare low-wind conditions, the old way of estimating the |
---|
| 793 | ! neutral wind speed at 10m leads to a negative value that causes the code |
---|
| 794 | ! to crash. To prevent this a threshold of 0.25m/s is imposed. |
---|
| 795 | ztmp0 = MAX( 0.25 , U_zu/(1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u)) ) ! U_n10 |
---|
| 796 | ztmp0 = cd_neutral_10m(ztmp0) ! Cd_n10 |
---|
| 797 | sqrt_Cd_n10 = sqrt(ztmp0) |
---|
| 798 | |
---|
| 799 | Ce_n10 = 1.e-3 * (34.6 * sqrt_Cd_n10) ! L&Y 2004 eq. (6b) |
---|
| 800 | stab = 0.5 + sign(0.5,zeta_u) ! update stability |
---|
| 801 | Ch_n10 = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab)) ! L&Y 2004 eq. (6c-6d) |
---|
| 802 | |
---|
| 803 | !! Update of transfer coefficients: |
---|
| 804 | ztmp1 = 1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u) ! L&Y 2004 eq. (10a) |
---|
| 805 | Cd = ztmp0 / ( ztmp1*ztmp1 ) |
---|
| 806 | sqrt_Cd = SQRT( Cd ) |
---|
[3294] | 807 | ENDIF |
---|
[4990] | 808 | ! |
---|
| 809 | ztmp0 = (LOG(zu/10.) - zpsi_h_u) / vkarmn / sqrt_Cd_n10 |
---|
| 810 | ztmp2 = sqrt_Cd / sqrt_Cd_n10 |
---|
| 811 | ztmp1 = 1. + Ch_n10*ztmp0 |
---|
| 812 | Ch = Ch_n10*ztmp2 / ztmp1 ! L&Y 2004 eq. (10b) |
---|
| 813 | ! |
---|
| 814 | ztmp1 = 1. + Ce_n10*ztmp0 |
---|
| 815 | Ce = Ce_n10*ztmp2 / ztmp1 ! L&Y 2004 eq. (10c) |
---|
| 816 | ! |
---|
[888] | 817 | END DO |
---|
[4990] | 818 | |
---|
| 819 | CALL wrk_dealloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd ) |
---|
| 820 | CALL wrk_dealloc( jpi,jpj, zeta_u, stab ) |
---|
| 821 | CALL wrk_dealloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 ) |
---|
| 822 | |
---|
| 823 | IF( .NOT. l_zt_equal_zu ) CALL wrk_dealloc( jpi,jpj, zeta_t ) |
---|
| 824 | |
---|
| 825 | IF( nn_timing == 1 ) CALL timing_stop('turb_core_2z') |
---|
| 826 | ! |
---|
| 827 | END SUBROUTINE turb_core_2z |
---|
| 828 | |
---|
| 829 | |
---|
| 830 | FUNCTION cd_neutral_10m( zw10 ) |
---|
| 831 | !!---------------------------------------------------------------------- |
---|
| 832 | !! Estimate of the neutral drag coefficient at 10m as a function |
---|
| 833 | !! of neutral wind speed at 10m |
---|
[888] | 834 | !! |
---|
[4990] | 835 | !! Origin: Large & Yeager 2008 eq.(11a) and eq.(11b) |
---|
| 836 | !! |
---|
| 837 | !! Author: L. Brodeau, june 2014 |
---|
| 838 | !!---------------------------------------------------------------------- |
---|
| 839 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zw10 ! scalar wind speed at 10m (m/s) |
---|
| 840 | REAL(wp), DIMENSION(jpi,jpj) :: cd_neutral_10m |
---|
[2715] | 841 | ! |
---|
[4990] | 842 | REAL(wp), DIMENSION(:,:), POINTER :: rgt33 |
---|
| 843 | !!---------------------------------------------------------------------- |
---|
[3294] | 844 | ! |
---|
[4990] | 845 | CALL wrk_alloc( jpi,jpj, rgt33 ) |
---|
| 846 | ! |
---|
| 847 | !! When wind speed > 33 m/s => Cyclone conditions => special treatment |
---|
| 848 | rgt33 = 0.5_wp + SIGN( 0.5_wp, (zw10 - 33._wp) ) ! If zw10 < 33. => 0, else => 1 |
---|
| 849 | cd_neutral_10m = 1.e-3 * ( & |
---|
[5065] | 850 | & (1._wp - rgt33)*( 2.7_wp/zw10 + 0.142_wp + zw10/13.09_wp - 3.14807E-10*zw10**6) & ! zw10< 33. |
---|
[4990] | 851 | & + rgt33 * 2.34 ) ! zw10 >= 33. |
---|
| 852 | ! |
---|
| 853 | CALL wrk_dealloc( jpi,jpj, rgt33) |
---|
| 854 | ! |
---|
| 855 | END FUNCTION cd_neutral_10m |
---|
[888] | 856 | |
---|
| 857 | |
---|
[4990] | 858 | FUNCTION psi_m(pta) !! Psis, L&Y 2004 eq. (8c), (8d), (8e) |
---|
[2715] | 859 | !------------------------------------------------------------------------------- |
---|
[4990] | 860 | ! universal profile stability function for momentum |
---|
| 861 | !------------------------------------------------------------------------------- |
---|
| 862 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pta |
---|
| 863 | ! |
---|
[888] | 864 | REAL(wp), DIMENSION(jpi,jpj) :: psi_m |
---|
[3294] | 865 | REAL(wp), DIMENSION(:,:), POINTER :: X2, X, stabit |
---|
[2715] | 866 | !------------------------------------------------------------------------------- |
---|
[4990] | 867 | ! |
---|
[3294] | 868 | CALL wrk_alloc( jpi,jpj, X2, X, stabit ) |
---|
[4990] | 869 | ! |
---|
| 870 | X2 = SQRT( ABS( 1. - 16.*pta ) ) ; X2 = MAX( X2 , 1. ) ; X = SQRT( X2 ) |
---|
| 871 | stabit = 0.5 + SIGN( 0.5 , pta ) |
---|
| 872 | psi_m = -5.*pta*stabit & ! Stable |
---|
| 873 | & + (1. - stabit)*(2.*LOG((1. + X)*0.5) + LOG((1. + X2)*0.5) - 2.*ATAN(X) + rpi*0.5) ! Unstable |
---|
| 874 | ! |
---|
[3294] | 875 | CALL wrk_dealloc( jpi,jpj, X2, X, stabit ) |
---|
[2715] | 876 | ! |
---|
[4990] | 877 | END FUNCTION psi_m |
---|
[888] | 878 | |
---|
| 879 | |
---|
[4990] | 880 | FUNCTION psi_h( pta ) !! Psis, L&Y 2004 eq. (8c), (8d), (8e) |
---|
[2715] | 881 | !------------------------------------------------------------------------------- |
---|
[4990] | 882 | ! universal profile stability function for temperature and humidity |
---|
| 883 | !------------------------------------------------------------------------------- |
---|
| 884 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pta |
---|
[2715] | 885 | ! |
---|
| 886 | REAL(wp), DIMENSION(jpi,jpj) :: psi_h |
---|
[4990] | 887 | REAL(wp), DIMENSION(:,:), POINTER :: X2, X, stabit |
---|
[2715] | 888 | !------------------------------------------------------------------------------- |
---|
[4990] | 889 | ! |
---|
[3294] | 890 | CALL wrk_alloc( jpi,jpj, X2, X, stabit ) |
---|
[4990] | 891 | ! |
---|
| 892 | X2 = SQRT( ABS( 1. - 16.*pta ) ) ; X2 = MAX( X2 , 1. ) ; X = SQRT( X2 ) |
---|
| 893 | stabit = 0.5 + SIGN( 0.5 , pta ) |
---|
| 894 | psi_h = -5.*pta*stabit & ! Stable |
---|
| 895 | & + (1. - stabit)*(2.*LOG( (1. + X2)*0.5 )) ! Unstable |
---|
| 896 | ! |
---|
[3294] | 897 | CALL wrk_dealloc( jpi,jpj, X2, X, stabit ) |
---|
[2715] | 898 | ! |
---|
[4990] | 899 | END FUNCTION psi_h |
---|
| 900 | |
---|
[888] | 901 | !!====================================================================== |
---|
| 902 | END MODULE sbcblk_core |
---|