[702] | 1 | MODULE sbcblk_core |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE sbcblk_core *** |
---|
| 4 | !! Ocean forcing: momentum, heat and freshwater flux formulation |
---|
| 5 | !!===================================================================== |
---|
| 6 | !! History : 9.0 ! 04-08 (U. Schweckendiek) Original code |
---|
| 7 | !! ! 05-04 (L. Brodeau, A.M. Treguier) additions: |
---|
| 8 | !! - new bulk routine for efficiency |
---|
| 9 | !! - WINDS ARE NOW ASSUMED TO BE AT T POINTS in input files !!!! |
---|
| 10 | !! - file names and file characteristics in namelist |
---|
| 11 | !! - Implement reading of 6-hourly fields |
---|
| 12 | !! ! 06-06 (G. Madec) sbc rewritting |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
| 14 | |
---|
| 15 | !!---------------------------------------------------------------------- |
---|
| 16 | !! sbc_blk_core : bulk formulation as ocean surface boundary condition |
---|
| 17 | !! (forced mode, CORE bulk formulea) |
---|
| 18 | !! blk_oce_core : ocean: computes momentum, heat and freshwater fluxes |
---|
| 19 | !! blk_ice_core : ice : computes momentum, heat and freshwater fluxes |
---|
| 20 | !! turb_core : computes the CORE turbulent transfer coefficients |
---|
| 21 | !!---------------------------------------------------------------------- |
---|
| 22 | USE oce ! ocean dynamics and tracers |
---|
| 23 | USE dom_oce ! ocean space and time domain |
---|
| 24 | USE phycst ! physical constants |
---|
| 25 | USE daymod ! calendar |
---|
| 26 | USE ocfzpt ! ocean freezing point |
---|
| 27 | USE fldread ! read input fields |
---|
| 28 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
| 29 | USE iom ! I/O manager library |
---|
| 30 | USE in_out_manager ! I/O manager |
---|
| 31 | USE lib_mpp ! distribued memory computing library |
---|
| 32 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 33 | |
---|
| 34 | IMPLICIT NONE |
---|
| 35 | PRIVATE |
---|
| 36 | |
---|
| 37 | PUBLIC sbc_blk_core ! routine called in sbcmod module |
---|
| 38 | PUBLIC blk_ice_core ! routine called in sbc_ice_lim module |
---|
| 39 | |
---|
| 40 | INTEGER , PARAMETER :: jpfld = 8 ! maximum number of files to read |
---|
| 41 | INTEGER , PARAMETER :: jp_wndi = 1 ! index of 10m wind velocity (i-component) (m/s) at T-point |
---|
| 42 | INTEGER , PARAMETER :: jp_wndj = 2 ! index of 10m wind velocity (j-component) (m/s) at T-point |
---|
| 43 | INTEGER , PARAMETER :: jp_humi = 3 ! index of specific humidity ( % ) |
---|
| 44 | INTEGER , PARAMETER :: jp_qsr = 4 ! index of solar heat (W/m2) |
---|
| 45 | INTEGER , PARAMETER :: jp_qlw = 5 ! index of Long wave (W/m2) |
---|
| 46 | INTEGER , PARAMETER :: jp_tair = 6 ! index of 10m air temperature (Kelvin) |
---|
| 47 | INTEGER , PARAMETER :: jp_prec = 7 ! index of total precipitation (rain+snow) (Kg/m2/s) |
---|
| 48 | INTEGER , PARAMETER :: jp_snow = 8 ! index of snow (solid prcipitation) (kg/m2/s) |
---|
| 49 | TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf ! structure of input fields (file informations, fields read) |
---|
| 50 | |
---|
| 51 | !! * CORE bulk parameters |
---|
| 52 | REAL(wp), PARAMETER :: rhoa = 1.22 ! air density |
---|
| 53 | REAL(wp), PARAMETER :: cpa = 1000.5 ! specific heat of air |
---|
| 54 | REAL(wp), PARAMETER :: Lv = 2.5e6 ! latent heat of vaporization |
---|
| 55 | REAL(wp), PARAMETER :: Ls = 2.839e6 ! latent heat of sublimation |
---|
| 56 | REAL(wp), PARAMETER :: Stef = 5.67e-8 ! Stefan Boltzmann constant |
---|
| 57 | REAL(wp), PARAMETER :: Cice = 1.63e-3 ! transfer coefficient over ice |
---|
| 58 | |
---|
[875] | 59 | LOGICAL :: ln_2m = .FALSE. !: logical flag for height of air temp. and hum |
---|
| 60 | REAL(wp) :: alpha_precip=1. !: multiplication factor for precipitation |
---|
| 61 | |
---|
[702] | 62 | !! * Substitutions |
---|
| 63 | # include "domzgr_substitute.h90" |
---|
| 64 | # include "vectopt_loop_substitute.h90" |
---|
| 65 | !!---------------------------------------------------------------------- |
---|
| 66 | !! OPA 9.0 , LOCEAN-IPSL (2006) |
---|
| 67 | !! $Header: $ |
---|
| 68 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
| 69 | !!---------------------------------------------------------------------- |
---|
| 70 | |
---|
| 71 | CONTAINS |
---|
| 72 | |
---|
| 73 | SUBROUTINE sbc_blk_core( kt ) |
---|
| 74 | !!--------------------------------------------------------------------- |
---|
| 75 | !! *** ROUTINE sbc_blk_core *** |
---|
| 76 | !! |
---|
| 77 | !! ** Purpose : provide at each time step the surface ocean fluxes |
---|
| 78 | !! (momentum, heat, freshwater and runoff) |
---|
| 79 | !! |
---|
| 80 | !! ** Method : READ each fluxes in NetCDF files |
---|
[710] | 81 | !! The i-component of the stress utau (N/m2) |
---|
| 82 | !! The j-component of the stress vtau (N/m2) |
---|
[702] | 83 | !! the net downward heat flux qtot (watt/m2) |
---|
| 84 | !! the net downward radiative flux qsr (watt/m2) |
---|
| 85 | !! the net upward water (evapo - precip) emp (kg/m2/s) |
---|
| 86 | !! Assumptions made: |
---|
| 87 | !! - each file content an entire year (read record, not the time axis) |
---|
| 88 | !! - first and last record are part of the previous and next year |
---|
| 89 | !! (useful for time interpolation) |
---|
| 90 | !! - the number of records is 2 + 365*24 / freqh(jf) |
---|
| 91 | !! or 366 in leap year case |
---|
| 92 | !! |
---|
| 93 | !! C A U T I O N : never mask the surface stress fields |
---|
| 94 | !! the stress is assumed to be in the mesh referential |
---|
| 95 | !! i.e. the (i,j) referential |
---|
| 96 | !! |
---|
| 97 | !! ** Action : defined at each time-step at the air-sea interface |
---|
| 98 | !! - utau & vtau : stress components in geographical ref. |
---|
| 99 | !! - qns & qsr : non solar and solar heat fluxes |
---|
| 100 | !! - emp : evap - precip (volume flux) |
---|
| 101 | !! - emps : evap - precip (concentration/dillution) |
---|
| 102 | !!---------------------------------------------------------------------- |
---|
| 103 | INTEGER, INTENT( in ) :: kt ! ocean time step |
---|
| 104 | !! |
---|
| 105 | INTEGER :: jf ! dummy indices |
---|
| 106 | INTEGER :: ierror ! return error code |
---|
| 107 | !! |
---|
| 108 | CHARACTER(len=100) :: cn_dir ! Root directory for location of core files |
---|
| 109 | TYPE(FLD_N), DIMENSION(jpfld) :: slf_i ! array of namelist informations on the fields to read |
---|
| 110 | TYPE(FLD_N) :: sn_wndi, sn_wndj, sn_humi, sn_qsr ! informations about the fields to be read |
---|
| 111 | TYPE(FLD_N) :: sn_qlw , sn_tair, sn_prec, sn_snow ! " " |
---|
[875] | 112 | NAMELIST/namsbc_core/ cn_dir, ln_2m, alpha_precip, sn_wndi, sn_wndj, sn_humi, sn_qsr, & |
---|
| 113 | & sn_qlw , sn_tair, sn_prec, sn_snow |
---|
[702] | 114 | !!--------------------------------------------------------------------- |
---|
| 115 | |
---|
| 116 | ! ! ====================== ! |
---|
| 117 | IF( kt == nit000 ) THEN ! First call kt=nit000 ! |
---|
| 118 | ! ! ====================== ! |
---|
| 119 | ! set file information (default values) |
---|
| 120 | cn_dir = './' ! directory in which the model is executed |
---|
| 121 | |
---|
| 122 | ! (NB: frequency positive => hours, negative => months) |
---|
| 123 | ! ! file ! frequency ! variable ! time intep ! clim ! starting ! |
---|
| 124 | ! ! name ! (hours) ! name ! (T/F) ! (0/1) ! record ! |
---|
| 125 | sn_wndi = FLD_N( 'uwnd10m' , 24. , 'u_10' , .FALSE. , 0 , 0 ) |
---|
| 126 | sn_wndj = FLD_N( 'vwnd10m' , 24. , 'v_10' , .FALSE. , 0 , 0 ) |
---|
| 127 | sn_qsr = FLD_N( 'qsw' , 24. , 'qsw' , .FALSE. , 0 , 0 ) |
---|
| 128 | sn_qlw = FLD_N( 'qlw' , 24. , 'qlw' , .FALSE. , 0 , 0 ) |
---|
| 129 | sn_tair = FLD_N( 'tair10m' , 24. , 't_10' , .FALSE. , 0 , 0 ) |
---|
| 130 | sn_humi = FLD_N( 'humi10m' , 24. , 'q_10' , .FALSE. , 0 , 0 ) |
---|
| 131 | sn_prec = FLD_N( 'precip' , -12. , 'precip' , .TRUE. , 0 , 0 ) |
---|
| 132 | sn_snow = FLD_N( 'snow' , -12. , 'snow' , .TRUE. , 0 , 0 ) |
---|
| 133 | |
---|
| 134 | REWIND( numnam ) ! ... read in namlist namsbc_core |
---|
| 135 | READ ( numnam, namsbc_core ) |
---|
| 136 | |
---|
| 137 | ! store namelist information in an array |
---|
| 138 | slf_i(jp_wndi) = sn_wndi ; slf_i(jp_wndj) = sn_wndj |
---|
| 139 | slf_i(jp_qsr ) = sn_qsr ; slf_i(jp_qlw ) = sn_qlw |
---|
| 140 | slf_i(jp_tair) = sn_tair ; slf_i(jp_humi) = sn_humi |
---|
| 141 | slf_i(jp_prec) = sn_prec ; slf_i(jp_snow) = sn_snow |
---|
| 142 | |
---|
| 143 | ! set sf structure |
---|
| 144 | ALLOCATE( sf(jpfld), STAT=ierror ) |
---|
| 145 | IF( ierror > 0 ) THEN |
---|
| 146 | CALL ctl_stop( 'sbc_blk_core: unable to allocate sf structure' ) ; RETURN |
---|
| 147 | ENDIF |
---|
| 148 | |
---|
| 149 | DO jf = 1, jpfld |
---|
| 150 | WRITE(sf(jf)%clrootname,'(a,a)' ) TRIM( cn_dir ), TRIM( slf_i(jf)%clname ) |
---|
| 151 | sf(jf)%freqh = slf_i(jf)%freqh |
---|
| 152 | sf(jf)%clvar = slf_i(jf)%clvar |
---|
| 153 | sf(jf)%ln_tint = slf_i(jf)%ln_tint |
---|
| 154 | sf(jf)%nclim = slf_i(jf)%nclim |
---|
| 155 | sf(jf)%nstrec = slf_i(jf)%nstrec |
---|
| 156 | END DO |
---|
| 157 | |
---|
| 158 | IF(lwp) THEN ! control print |
---|
| 159 | WRITE(numout,*) |
---|
| 160 | WRITE(numout,*) 'sbc_blk_core : flux formulattion for ocean surface boundary condition' |
---|
| 161 | WRITE(numout,*) '~~~~~~~~~~~~ ' |
---|
| 162 | WRITE(numout,*) ' namsbc_core Namelist' |
---|
[875] | 163 | WRITE(numout,*) ' ln_2m = ', ln_2m |
---|
| 164 | WRITE(numout,*) ' alpha_precip = ', alpha_precip |
---|
[702] | 165 | WRITE(numout,*) ' list of files and frequency (>0: in hours ; <0 in months)' |
---|
| 166 | DO jf = 1, jpfld |
---|
| 167 | WRITE(numout,*) ' file root name: ' , TRIM( sf(jf)%clrootname ), & |
---|
| 168 | & ' variable name: ' , TRIM( sf(jf)%clvar ) |
---|
| 169 | WRITE(numout,*) ' frequency: ' , sf(jf)%freqh , & |
---|
| 170 | & ' time interp: ' , sf(jf)%ln_tint , & |
---|
| 171 | & ' climatology: ' , sf(jf)%nclim , & |
---|
| 172 | & ' starting record: ', sf(jf)%nstrec |
---|
| 173 | END DO |
---|
[875] | 174 | IF( ln_2m ) THEN |
---|
| 175 | WRITE(numout,*) ' Calling TURB_CORE_2Z for bulk transfert coefficients' |
---|
| 176 | ELSE |
---|
| 177 | WRITE(numout,*) ' Calling TURB_CORE_1Z for bulk transfert coefficients' |
---|
| 178 | ENDIF |
---|
| 179 | WRITE(numout,*) |
---|
| 180 | ! |
---|
[702] | 181 | ENDIF |
---|
| 182 | ! |
---|
| 183 | ENDIF |
---|
| 184 | |
---|
| 185 | CALL fld_read( kt, nn_fsbc, sf ) ! Read input fields and provides the |
---|
| 186 | ! ! input fields at the current time-step |
---|
| 187 | |
---|
[758] | 188 | IF( MOD( kt - 1, nn_fsbc ) == 0 ) THEN |
---|
[702] | 189 | |
---|
[758] | 190 | CALL blk_oce_core( sst_m, ssu_m, ssv_m ) ! set the ocean surface fluxes |
---|
[702] | 191 | |
---|
[758] | 192 | ENDIF |
---|
[702] | 193 | ! ! using CORE bulk formulea |
---|
| 194 | END SUBROUTINE sbc_blk_core |
---|
| 195 | |
---|
| 196 | |
---|
| 197 | SUBROUTINE blk_oce_core( pst, pu, pv ) |
---|
| 198 | !!--------------------------------------------------------------------- |
---|
| 199 | !! *** ROUTINE blk_core *** |
---|
| 200 | !! |
---|
| 201 | !! ** Purpose : provide the momentum, heat and freshwater fluxes at |
---|
| 202 | !! the ocean surface at each time step |
---|
| 203 | !! |
---|
| 204 | !! ** Method : CORE bulk formulea for the ocean using atmospheric |
---|
| 205 | !! fields read in sbc_read |
---|
| 206 | !! |
---|
| 207 | !! ** Outputs : - utau : i-component of the stress at U-point (N/m2) |
---|
| 208 | !! - vtau : j-component of the stress at V-point (N/m2) |
---|
| 209 | !! - qsr_oce : Solar heat flux over the ocean (W/m2) |
---|
| 210 | !! - qns_oce : Non Solar heat flux over the ocean (W/m2) |
---|
| 211 | !! - evap : Evaporation over the ocean (kg/m2/s) |
---|
| 212 | !! - tprecip : Total precipitation (Kg/m2/s) |
---|
| 213 | !! - sprecip : Solid precipitation (Kg/m2/s) |
---|
| 214 | !!--------------------------------------------------------------------- |
---|
| 215 | REAL(wp), INTENT(in), DIMENSION(jpi,jpj) :: pst ! surface temperature [Celcius] |
---|
| 216 | REAL(wp), INTENT(in), DIMENSION(jpi,jpj) :: pu ! surface current at U-point (i-component) [m/s] |
---|
| 217 | REAL(wp), INTENT(in), DIMENSION(jpi,jpj) :: pv ! surface current at V-point (j-component) [m/s] |
---|
| 218 | |
---|
| 219 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 220 | REAL(wp) :: zcoef_qsatw |
---|
| 221 | REAL(wp), DIMENSION(jpi,jpj) :: zwnd_i, zwnd_j ! wind speed components at T-point |
---|
| 222 | REAL(wp), DIMENSION(jpi,jpj) :: zwind_speed_t ! wind speed module at T-point ( = | U10m - Uoce | ) |
---|
| 223 | REAL(wp), DIMENSION(jpi,jpj) :: zqsatw ! specific humidity at pst |
---|
| 224 | REAL(wp), DIMENSION(jpi,jpj) :: zqlw, zqsb ! long wave and sensible heat fluxes |
---|
| 225 | REAL(wp), DIMENSION(jpi,jpj) :: zqla, zevap ! latent heat fluxes and evaporation |
---|
| 226 | REAL(wp), DIMENSION(jpi,jpj) :: Cd ! transfer coefficient for momentum (tau) |
---|
| 227 | REAL(wp), DIMENSION(jpi,jpj) :: Ch ! transfer coefficient for sensible heat (Q_sens) |
---|
| 228 | REAL(wp), DIMENSION(jpi,jpj) :: Ce ! tansfert coefficient for evaporation (Q_lat) |
---|
| 229 | REAL(wp), DIMENSION(jpi,jpj) :: zst ! surface temperature in Kelvin |
---|
[875] | 230 | REAL(wp), DIMENSION(jpi,jpj) :: zt_zu ! air temperature at wind speed height |
---|
| 231 | REAL(wp), DIMENSION(jpi,jpj) :: zq_zu ! air spec. hum. at wind speed height |
---|
[702] | 232 | !!--------------------------------------------------------------------- |
---|
| 233 | |
---|
| 234 | ! local scalars ( place there for vector optimisation purposes) |
---|
| 235 | zcoef_qsatw = 0.98 * 640380. / rhoa |
---|
| 236 | |
---|
| 237 | zst(:,:) = pst(:,:) + rt0 ! converte Celcius to Kelvin (and set minimum value far above 0 K) |
---|
| 238 | |
---|
| 239 | ! ----------------------------------------------------------------------------- ! |
---|
| 240 | ! 0 Wind components and module at T-point relative to the moving ocean ! |
---|
| 241 | ! ----------------------------------------------------------------------------- ! |
---|
| 242 | ! ... components ( U10m - U_oce ) at T-point (unmasked) |
---|
| 243 | zwnd_i(:,:) = 0.e0 |
---|
| 244 | zwnd_j(:,:) = 0.e0 |
---|
| 245 | #if defined key_vectopt_loop |
---|
| 246 | !CDIR COLLAPSE |
---|
| 247 | #endif |
---|
| 248 | DO jj = 2, jpjm1 |
---|
| 249 | DO ji = fs_2, fs_jpim1 ! vect. opt. |
---|
| 250 | zwnd_i(ji,jj) = ( sf(jp_wndi)%fnow(ji,jj) - 0.5 * ( pu(ji-1,jj ) + pu(ji,jj) ) ) |
---|
| 251 | zwnd_j(ji,jj) = ( sf(jp_wndj)%fnow(ji,jj) - 0.5 * ( pv(ji ,jj-1) + pv(ji,jj) ) ) |
---|
| 252 | END DO |
---|
| 253 | END DO |
---|
| 254 | CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. ) |
---|
| 255 | CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. ) |
---|
| 256 | ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked) |
---|
| 257 | !CDIR NOVERRCHK |
---|
| 258 | !CDIR COLLAPSE |
---|
| 259 | zwind_speed_t(:,:) = SQRT( zwnd_i(:,:) * zwnd_i(:,:) & |
---|
| 260 | & + zwnd_j(:,:) * zwnd_j(:,:) ) * tmask(:,:,1) |
---|
| 261 | |
---|
| 262 | ! ----------------------------------------------------------------------------- ! |
---|
| 263 | ! I Radiative FLUXES ! |
---|
| 264 | ! ----------------------------------------------------------------------------- ! |
---|
| 265 | |
---|
| 266 | ! ocean albedo assumed to be 0.066 |
---|
| 267 | !CDIR COLLAPSE |
---|
| 268 | qsr (:,:) = ( 1. - 0.066 ) * sf(jp_qsr)%fnow(:,:) * tmask(:,:,1) ! Short Wave |
---|
| 269 | !CDIR COLLAPSE |
---|
| 270 | zqlw(:,:) = ( sf(jp_qlw)%fnow(:,:) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:) ) * tmask(:,:,1) ! Long Wave |
---|
| 271 | |
---|
| 272 | ! ----------------------------------------------------------------------------- ! |
---|
| 273 | ! II Turbulent FLUXES ! |
---|
| 274 | ! ----------------------------------------------------------------------------- ! |
---|
| 275 | |
---|
| 276 | ! ... specific humidity at SST and IST |
---|
| 277 | !CDIR NOVERRCHK |
---|
| 278 | !CDIR COLLAPSE |
---|
| 279 | zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) |
---|
| 280 | |
---|
| 281 | ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point : |
---|
[875] | 282 | IF( ln_2m ) THEN |
---|
| 283 | !! If air temp. and spec. hum. are given at different height (2m) than wind (10m) : |
---|
| 284 | CALL TURB_CORE_2Z(2.,10., zst , sf(jp_tair)%fnow, & |
---|
| 285 | & zqsatw, sf(jp_humi)%fnow, zwind_speed_t, & |
---|
| 286 | & Cd, Ch, Ce, zt_zu, zq_zu ) |
---|
| 287 | ELSE |
---|
| 288 | !! If air temp. and spec. hum. are given at same height than wind (10m) : |
---|
| 289 | !gm bug? at the compiling phase, add a copy in temporary arrays... ==> check perf |
---|
| 290 | ! CALL TURB_CORE_1Z( 10., zst (:,:), sf(jp_tair)%fnow(:,:), & |
---|
| 291 | ! & zqsatw(:,:), sf(jp_humi)%fnow(:,:), zwind_speed_t(:,:), & |
---|
| 292 | ! & Cd(:,:), Ch(:,:), Ce(:,:) ) |
---|
| 293 | !gm bug |
---|
| 294 | CALL TURB_CORE_1Z( 10., zst , sf(jp_tair)%fnow, & |
---|
| 295 | & zqsatw, sf(jp_humi)%fnow, zwind_speed_t, & |
---|
| 296 | & Cd, Ch, Ce ) |
---|
| 297 | ENDIF |
---|
| 298 | |
---|
[874] | 299 | ! ... utau, vtau at U- and V_points, resp. |
---|
| 300 | ! Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines |
---|
[702] | 301 | zwnd_i(:,:) = rhoa * zwind_speed_t(:,:) * Cd(:,:) * zwnd_i(:,:) |
---|
| 302 | zwnd_j(:,:) = rhoa * zwind_speed_t(:,:) * Cd(:,:) * zwnd_j(:,:) |
---|
| 303 | DO jj = 1, jpjm1 |
---|
| 304 | DO ji = 1, fs_jpim1 |
---|
[874] | 305 | utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj ) ) |
---|
| 306 | vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji ,jj+1) ) |
---|
[702] | 307 | END DO |
---|
| 308 | END DO |
---|
| 309 | CALL lbc_lnk( utau(:,:), 'U', -1. ) |
---|
| 310 | CALL lbc_lnk( vtau(:,:), 'V', -1. ) |
---|
| 311 | |
---|
| 312 | ! Turbulent fluxes over ocean |
---|
| 313 | ! ----------------------------- |
---|
[875] | 314 | IF( ln_2m ) THEN |
---|
| 315 | ! Values of temp. and hum. adjusted to 10m must be used instead of 2m values |
---|
| 316 | zevap(:,:) = MAX( 0.e0, rhoa *Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) ) * zwind_speed_t(:,:) ) ! Evaporation |
---|
| 317 | zqsb (:,:) = rhoa*cpa*Ch(:,:)*( zst (:,:) - zt_zu(:,:) ) * zwind_speed_t(:,:) ! Sensible Heat |
---|
| 318 | ELSE |
---|
[702] | 319 | !CDIR COLLAPSE |
---|
[875] | 320 | zevap(:,:) = MAX( 0.e0, rhoa *Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:) ) * zwind_speed_t(:,:) ) ! Evaporation |
---|
[702] | 321 | !CDIR COLLAPSE |
---|
[875] | 322 | zqsb (:,:) = rhoa*cpa*Ch(:,:)*( zst (:,:) - sf(jp_tair)%fnow(:,:) ) * zwind_speed_t(:,:) ! Sensible Heat |
---|
| 323 | ENDIF |
---|
[702] | 324 | !CDIR COLLAPSE |
---|
[874] | 325 | zqla (:,:) = Lv * zevap(:,:) ! Latent Heat |
---|
[702] | 326 | |
---|
| 327 | ! ----------------------------------------------------------------------------- ! |
---|
| 328 | ! III Total FLUXES ! |
---|
| 329 | ! ----------------------------------------------------------------------------- ! |
---|
| 330 | |
---|
| 331 | !CDIR COLLAPSE |
---|
| 332 | qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:) ! Downward Non Solar flux |
---|
| 333 | !CDIR COLLAPSE |
---|
[875] | 334 | emp (:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:) * alpha_precip * tmask(:,:,1) |
---|
[702] | 335 | !CDIR COLLAPSE |
---|
[875] | 336 | emps(:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:) * alpha_precip * tmask(:,:,1) |
---|
[702] | 337 | ! |
---|
| 338 | END SUBROUTINE blk_oce_core |
---|
| 339 | |
---|
| 340 | |
---|
| 341 | SUBROUTINE blk_ice_core( pst , pui , pvi , palb , & |
---|
| 342 | & p_taui, p_tauj, p_qns , p_qsr, & |
---|
| 343 | & p_qla , p_dqns, p_dqla, & |
---|
| 344 | & p_tpr , p_spr , & |
---|
| 345 | & p_fr1 , p_fr2 ) |
---|
| 346 | !!--------------------------------------------------------------------- |
---|
| 347 | !! *** ROUTINE blk_ice_core *** |
---|
| 348 | !! |
---|
| 349 | !! ** Purpose : provide the surface boundary condition over sea-ice |
---|
| 350 | !! |
---|
| 351 | !! ** Method : compute momentum, heat and freshwater exchanged |
---|
| 352 | !! between atmosphere and sea-ice using CORE bulk |
---|
| 353 | !! formulea, ice variables and read atmmospheric fields. |
---|
| 354 | !! NB: ice drag coefficient is assumed to be a constant |
---|
| 355 | !! |
---|
| 356 | !! caution : the net upward water flux has with mm/day unit |
---|
| 357 | !!--------------------------------------------------------------------- |
---|
| 358 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: pst ! ice surface temperature (>0, =rt0 over land) [Kelvin] |
---|
| 359 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: pui ! ice surface velocity (i-component, I-point) [m/s] |
---|
| 360 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: pvi ! ice surface velocity (j-component, I-point) [m/s] |
---|
| 361 | REAL(wp), INTENT(in ), DIMENSION(jpi,jpj) :: palb ! ice albedo (clear sky) (alb_ice_cs) [%] |
---|
| 362 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_taui ! surface ice stress at I-point (i-component) [N/m2] |
---|
| 363 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_tauj ! surface ice stress at I-point (j-component) [N/m2] |
---|
| 364 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_qns ! non solar heat flux over ice (T-point) [W/m2] |
---|
| 365 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_qsr ! solar heat flux over ice (T-point) [W/m2] |
---|
| 366 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_qla ! latent heat flux over ice (T-point) [W/m2] |
---|
| 367 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_dqns ! non solar heat sensistivity (T-point) [W/m2] |
---|
| 368 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_dqla ! latent heat sensistivity (T-point) [W/m2] |
---|
| 369 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_tpr ! total precipitation (T-point) [Kg/m2/s] |
---|
| 370 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_spr ! solid precipitation (T-point) [Kg/m2/s] |
---|
| 371 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_fr1 ! 1sr fraction of qsr penetration in ice [%] |
---|
| 372 | REAL(wp), INTENT( out), DIMENSION(jpi,jpj) :: p_fr2 ! 2nd fraction of qsr penetration in ice [%] |
---|
| 373 | !! |
---|
| 374 | INTEGER :: ji, jj ! dummy loop indices |
---|
| 375 | REAL(wp) :: zst3 |
---|
[810] | 376 | REAL(wp) :: zcoef_wnorm, zcoef_dqlw, zcoef_dqla, zcoef_dqsb |
---|
| 377 | REAL(wp) :: zcoef_frca ! fractional cloud amount |
---|
[702] | 378 | REAL(wp) :: zwnorm_f, zwndi_f , zwndj_f ! relative wind module and components at F-point |
---|
| 379 | REAL(wp) :: zwndi_t , zwndj_t ! relative wind components at T-point |
---|
| 380 | REAL(wp), DIMENSION(jpi,jpj) :: z_wnds_t ! wind speed ( = | U10m - U_ice | ) at T-point |
---|
| 381 | REAL(wp), DIMENSION(jpi,jpj) :: z_qlw ! long wave heat flux over ice |
---|
| 382 | REAL(wp), DIMENSION(jpi,jpj) :: z_qsb ! sensible heat flux over ice |
---|
| 383 | REAL(wp), DIMENSION(jpi,jpj) :: z_dqlw ! sensible heat flux over ice |
---|
| 384 | REAL(wp), DIMENSION(jpi,jpj) :: z_dqsb ! sensible heat flux over ice |
---|
| 385 | !!--------------------------------------------------------------------- |
---|
| 386 | |
---|
| 387 | ! local scalars ( place there for vector optimisation purposes) |
---|
| 388 | zcoef_wnorm = rhoa * Cice |
---|
| 389 | zcoef_dqlw = 4.0 * 0.95 * Stef |
---|
[810] | 390 | zcoef_dqla = -Ls * Cice * 11637800. * (-5897.8) |
---|
[702] | 391 | zcoef_dqsb = rhoa * cpa * Cice |
---|
[810] | 392 | zcoef_frca = 1.0 - 0.3 |
---|
[702] | 393 | |
---|
| 394 | !!gm brutal.... |
---|
| 395 | z_wnds_t(:,:) = 0.e0 |
---|
| 396 | p_taui (:,:) = 0.e0 |
---|
| 397 | p_tauj (:,:) = 0.e0 |
---|
| 398 | !!gm end |
---|
| 399 | |
---|
| 400 | ! ----------------------------------------------------------------------------- ! |
---|
| 401 | ! Wind components and module relative to the moving ocean at I and T-point ! |
---|
| 402 | ! ----------------------------------------------------------------------------- ! |
---|
| 403 | ! ... components ( U10m - U_oce ) at I-point (F-point with sea-ice indexation) (unmasked) |
---|
| 404 | ! and scalar wind at T-point ( = | U10m - U_ice | ) (masked) |
---|
| 405 | #if defined key_vectopt_loop |
---|
| 406 | !CDIR COLLAPSE |
---|
| 407 | #endif |
---|
| 408 | !CDIR NOVERRCHK |
---|
| 409 | DO jj = 2, jpjm1 |
---|
| 410 | DO ji = fs_2, fs_jpim1 |
---|
| 411 | ! ... scalar wind at I-point (fld being at T-point) |
---|
| 412 | zwndi_f = 0.25 * ( sf(jp_wndi)%fnow(ji-1,jj ) + sf(jp_wndi)%fnow(ji ,jj ) & |
---|
| 413 | & + sf(jp_wndi)%fnow(ji-1,jj-1) + sf(jp_wndi)%fnow(ji ,jj-1) ) - pui(ji,jj) |
---|
| 414 | zwndj_f = 0.25 * ( sf(jp_wndj)%fnow(ji-1,jj ) + sf(jp_wndj)%fnow(ji ,jj ) & |
---|
| 415 | & + sf(jp_wndj)%fnow(ji-1,jj-1) + sf(jp_wndj)%fnow(ji ,jj-1) ) - pvi(ji,jj) |
---|
| 416 | zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f ) |
---|
| 417 | ! ... ice stress at I-point |
---|
| 418 | p_taui(ji,jj) = zwnorm_f * zwndi_f |
---|
| 419 | p_tauj(ji,jj) = zwnorm_f * zwndj_f |
---|
| 420 | ! ... scalar wind at T-point (fld being at T-point) |
---|
| 421 | zwndi_t = sf(jp_wndi)%fnow(ji,jj) - 0.25 * ( pui(ji,jj+1) + pui(ji+1,jj+1) & |
---|
| 422 | & + pui(ji,jj ) + pui(ji+1,jj ) ) |
---|
| 423 | zwndj_t = sf(jp_wndj)%fnow(ji,jj) - 0.25 * ( pvi(ji,jj+1) + pvi(ji+1,jj+1) & |
---|
| 424 | & + pvi(ji,jj ) + pvi(ji+1,jj ) ) |
---|
| 425 | z_wnds_t(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1) |
---|
| 426 | END DO |
---|
| 427 | END DO |
---|
| 428 | CALL lbc_lnk( p_taui , 'I', -1. ) |
---|
| 429 | CALL lbc_lnk( p_tauj , 'I', -1. ) |
---|
| 430 | CALL lbc_lnk( z_wnds_t, 'T', 1. ) |
---|
| 431 | |
---|
| 432 | ! ----------------------------------------------------------------------------- ! |
---|
| 433 | ! I Radiative FLUXES ! |
---|
| 434 | ! ----------------------------------------------------------------------------- ! |
---|
| 435 | !CDIR COLLAPSE |
---|
| 436 | DO jj = 1, jpj |
---|
| 437 | DO ji = 1, jpi |
---|
| 438 | zst3 = pst(ji,jj) * pst(ji,jj) * pst(ji,jj) |
---|
| 439 | p_qsr(ji,jj) = ( 1. - palb(ji,jj) ) * sf(jp_qsr)%fnow(ji,jj) * tmask(ji,jj,1) ! Short Wave (sw) |
---|
| 440 | z_qlw(ji,jj) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj) & ! Long Wave (lw) |
---|
| 441 | & - Stef * pst(ji,jj) * zst3 ) * tmask(ji,jj,1) |
---|
| 442 | z_dqlw(ji,jj) = zcoef_dqlw * zst3 ! lw sensitivity |
---|
| 443 | END DO |
---|
| 444 | END DO |
---|
| 445 | |
---|
| 446 | ! ----------------------------------------------------------------------------- ! |
---|
| 447 | ! II Turbulent FLUXES ! |
---|
| 448 | ! ----------------------------------------------------------------------------- ! |
---|
| 449 | |
---|
| 450 | ! ... turbulent heat fluxes |
---|
| 451 | !CDIR COLLAPSE |
---|
| 452 | z_qsb(:,:) = rhoa * cpa * Cice * z_wnds_t(:,:) * ( pst(:,:) - sf(jp_tair)%fnow(:,:) ) ! Sensible Heat |
---|
| 453 | !CDIR NOVERRCHK |
---|
| 454 | !CDIR COLLAPSE |
---|
[874] | 455 | p_qla(:,:) = MAX( 0.e0, rhoa * Ls * Cice * z_wnds_t(:,:) & ! Latent Heat |
---|
| 456 | & * ( 11637800. * EXP( -5897.8 / pst(:,:) ) / rhoa - sf(jp_humi)%fnow(:,:) ) ) |
---|
[702] | 457 | |
---|
| 458 | ! Latent heat sensitivity for ice (Dqla/Dt) |
---|
| 459 | !CDIR NOVERRCHK |
---|
| 460 | !CDIR COLLAPSE |
---|
| 461 | p_dqla(:,:) = zcoef_dqla * z_wnds_t(:,:) / ( pst(:,:) * pst(:,:) ) * EXP( -5897.8 / pst(:,:) ) |
---|
| 462 | |
---|
| 463 | ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice) |
---|
| 464 | !CDIR COLLAPSE |
---|
| 465 | z_dqsb(:,:) = zcoef_dqsb * z_wnds_t(:,:) |
---|
| 466 | |
---|
| 467 | ! ----------------------------------------------------------------------------- ! |
---|
| 468 | ! III Total FLUXES ! |
---|
| 469 | ! ----------------------------------------------------------------------------- ! |
---|
| 470 | |
---|
| 471 | !CDIR COLLAPSE |
---|
| 472 | p_qns (:,:) = z_qlw (:,:) - z_qsb (:,:) - p_qla (:,:) ! Downward Non Solar flux |
---|
| 473 | !CDIR COLLAPSE |
---|
| 474 | p_dqns(:,:) = - ( z_dqlw(:,:) + z_dqsb(:,:) + p_dqla(:,:) ) ! Total non solar heat flux sensitivity for ice |
---|
| 475 | |
---|
| 476 | |
---|
| 477 | !-------------------------------------------------------------------- |
---|
| 478 | ! FRACTIONs of net shortwave radiation which is not absorbed in the |
---|
| 479 | ! thin surface layer and penetrates inside the ice cover |
---|
[810] | 480 | ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 ) |
---|
[702] | 481 | |
---|
| 482 | !CDIR COLLAPSE |
---|
[810] | 483 | p_fr1(:,:) = ( 0.18 * ( 1.0 - zcoef_frca ) + 0.35 * zcoef_frca ) |
---|
[702] | 484 | !CDIR COLLAPSE |
---|
[810] | 485 | p_fr2(:,:) = ( 0.82 * ( 1.0 - zcoef_frca ) + 0.65 * zcoef_frca ) |
---|
[702] | 486 | |
---|
| 487 | !CDIR COLLAPSE |
---|
[875] | 488 | p_tpr(:,:) = sf(jp_prec)%fnow(:,:) * alpha_precip ! total precipitation [kg/m2/s] |
---|
[702] | 489 | !CDIR COLLAPSE |
---|
[875] | 490 | p_spr(:,:) = sf(jp_snow)%fnow(:,:) * alpha_precip ! solid precipitation [kg/m2/s] |
---|
[702] | 491 | ! |
---|
| 492 | END SUBROUTINE blk_ice_core |
---|
| 493 | |
---|
[875] | 494 | |
---|
[758] | 495 | SUBROUTINE TURB_CORE_1Z(zu, sst, T_a, q_sat, q_a, & |
---|
| 496 | & dU, Cd, Ch, Ce ) |
---|
| 497 | !!---------------------------------------------------------------------- |
---|
| 498 | !! *** ROUTINE turb_core *** |
---|
[702] | 499 | !! |
---|
[758] | 500 | !! ** Purpose : Computes turbulent transfert coefficients of surface |
---|
| 501 | !! fluxes according to Large & Yeager (2004) |
---|
[702] | 502 | !! |
---|
[758] | 503 | !! ** Method : I N E R T I A L D I S S I P A T I O N M E T H O D |
---|
| 504 | !! Momentum, Latent and sensible heat exchange coefficients |
---|
| 505 | !! Caution: this procedure should only be used in cases when air |
---|
| 506 | !! temperature (T_air), air specific humidity (q_air) and wind (dU) |
---|
| 507 | !! are provided at the same height 'zzu'! |
---|
[702] | 508 | !! |
---|
[758] | 509 | !! References : |
---|
| 510 | !! Large & Yeager, 2004 : ??? |
---|
| 511 | !! History : |
---|
| 512 | !! ! XX-XX (??? ) Original code |
---|
| 513 | !! 9.0 ! 05-08 (L. Brodeau) Rewriting and optimization |
---|
| 514 | !!---------------------------------------------------------------------- |
---|
| 515 | !! * Arguments |
---|
[702] | 516 | |
---|
[758] | 517 | REAL(wp), INTENT(in) :: zu ! altitude of wind measurement [m] |
---|
| 518 | REAL(wp), INTENT(in), DIMENSION(jpi,jpj) :: & |
---|
| 519 | sst, & ! sea surface temperature [Kelvin] |
---|
| 520 | T_a, & ! potential air temperature [Kelvin] |
---|
| 521 | q_sat, & ! sea surface specific humidity [kg/kg] |
---|
| 522 | q_a, & ! specific air humidity [kg/kg] |
---|
| 523 | dU ! wind module |U(zu)-U(0)| [m/s] |
---|
| 524 | REAL(wp), intent(out), DIMENSION(jpi,jpj) :: & |
---|
| 525 | Cd, & ! transfert coefficient for momentum (tau) |
---|
| 526 | Ch, & ! transfert coefficient for temperature (Q_sens) |
---|
| 527 | Ce ! transfert coefficient for evaporation (Q_lat) |
---|
[702] | 528 | |
---|
[758] | 529 | !! * Local declarations |
---|
| 530 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 531 | dU10, & ! dU [m/s] |
---|
| 532 | dT, & ! air/sea temperature differeence [K] |
---|
| 533 | dq, & ! air/sea humidity difference [K] |
---|
| 534 | Cd_n10, & ! 10m neutral drag coefficient |
---|
| 535 | Ce_n10, & ! 10m neutral latent coefficient |
---|
| 536 | Ch_n10, & ! 10m neutral sensible coefficient |
---|
| 537 | sqrt_Cd_n10, & ! root square of Cd_n10 |
---|
| 538 | sqrt_Cd, & ! root square of Cd |
---|
| 539 | T_vpot, & ! virtual potential temperature [K] |
---|
| 540 | T_star, & ! turbulent scale of tem. fluct. |
---|
| 541 | q_star, & ! turbulent humidity of temp. fluct. |
---|
| 542 | U_star, & ! turb. scale of velocity fluct. |
---|
| 543 | L, & ! Monin-Obukov length [m] |
---|
| 544 | zeta, & ! stability parameter at height zu |
---|
| 545 | U_n10, & ! neutral wind velocity at 10m [m] |
---|
| 546 | xlogt, xct, zpsi_h, zpsi_m |
---|
| 547 | !! |
---|
| 548 | INTEGER :: j_itt |
---|
| 549 | INTEGER, PARAMETER :: nb_itt = 3 |
---|
| 550 | INTEGER, DIMENSION(jpi,jpj) :: & |
---|
| 551 | stab ! 1st guess stability test integer |
---|
[702] | 552 | |
---|
[758] | 553 | REAL(wp), PARAMETER :: & |
---|
| 554 | grav = 9.8, & ! gravity |
---|
| 555 | kappa = 0.4 ! von Karman s constant |
---|
[702] | 556 | |
---|
[758] | 557 | !! * Start |
---|
| 558 | !! Air/sea differences |
---|
| 559 | dU10 = max(0.5, dU) ! we don't want to fall under 0.5 m/s |
---|
| 560 | dT = T_a - sst ! assuming that T_a is allready the potential temp. at zzu |
---|
| 561 | dq = q_a - q_sat |
---|
| 562 | !! |
---|
| 563 | !! Virtual potential temperature |
---|
| 564 | T_vpot = T_a*(1. + 0.608*q_a) |
---|
| 565 | !! |
---|
| 566 | !! Neutral Drag Coefficient |
---|
| 567 | stab = 0.5 + sign(0.5,dT) ! stable : stab = 1 ; unstable : stab = 0 |
---|
| 568 | Cd_n10 = 1E-3 * ( 2.7/dU10 + 0.142 + dU10/13.09 ) ! L & Y eq. (6a) |
---|
| 569 | sqrt_Cd_n10 = sqrt(Cd_n10) |
---|
| 570 | Ce_n10 = 1E-3 * ( 34.6 * sqrt_Cd_n10 ) ! L & Y eq. (6b) |
---|
| 571 | Ch_n10 = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1-stab)) ! L & Y eq. (6c), (6d) |
---|
| 572 | !! |
---|
| 573 | !! Initializing transfert coefficients with their first guess neutral equivalents : |
---|
| 574 | Cd = Cd_n10 ; Ce = Ce_n10 ; Ch = Ch_n10 ; sqrt_Cd = sqrt(Cd) |
---|
[702] | 575 | |
---|
[758] | 576 | !! * Now starting iteration loop |
---|
| 577 | DO j_itt=1, nb_itt |
---|
[702] | 578 | !! Turbulent scales : |
---|
[758] | 579 | U_star = sqrt_Cd*dU10 ! L & Y eq. (7a) |
---|
| 580 | T_star = Ch/sqrt_Cd*dT ! L & Y eq. (7b) |
---|
| 581 | q_star = Ce/sqrt_Cd*dq ! L & Y eq. (7c) |
---|
[702] | 582 | |
---|
| 583 | !! Estimate the Monin-Obukov length : |
---|
[758] | 584 | L = (U_star**2)/( kappa*grav*(T_star/T_vpot + q_star/(q_a + 1./0.608)) ) |
---|
[702] | 585 | |
---|
| 586 | !! Stability parameters : |
---|
[758] | 587 | zeta = zu/L ; zeta = sign( min(abs(zeta),10.0), zeta ) |
---|
| 588 | zpsi_h = psi_h(zeta) |
---|
| 589 | zpsi_m = psi_m(zeta) |
---|
[702] | 590 | |
---|
[758] | 591 | !! Shifting the wind speed to 10m and neutral stability : |
---|
| 592 | U_n10 = dU10*1./(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) ! L & Y eq. (9a) |
---|
[702] | 593 | |
---|
| 594 | !! Updating the neutral 10m transfer coefficients : |
---|
[758] | 595 | Cd_n10 = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09) ! L & Y eq. (6a) |
---|
| 596 | sqrt_Cd_n10 = sqrt(Cd_n10) |
---|
| 597 | Ce_n10 = 1E-3 * (34.6 * sqrt_Cd_n10) ! L & Y eq. (6b) |
---|
[702] | 598 | stab = 0.5 + sign(0.5,zeta) |
---|
[758] | 599 | Ch_n10 = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c), (6d) |
---|
| 600 | |
---|
| 601 | !! Shifting the neutral 10m transfer coefficients to ( zu , zeta ) : |
---|
| 602 | !! |
---|
| 603 | xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10) - zpsi_m) |
---|
| 604 | Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd) |
---|
| 605 | !! |
---|
| 606 | xlogt = log(zu/10.) - zpsi_h |
---|
| 607 | !! |
---|
| 608 | xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10 |
---|
| 609 | Ch = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct |
---|
| 610 | !! |
---|
| 611 | xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10 |
---|
| 612 | Ce = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct |
---|
| 613 | !! |
---|
[702] | 614 | END DO |
---|
[758] | 615 | !! |
---|
| 616 | END SUBROUTINE TURB_CORE_1Z |
---|
[702] | 617 | |
---|
[875] | 618 | |
---|
| 619 | SUBROUTINE TURB_CORE_2Z(zt, zu, sst, T_zt, q_sat, q_zt, dU, Cd, Ch, Ce, T_zu, q_zu) |
---|
| 620 | !!---------------------------------------------------------------------- |
---|
| 621 | !! *** ROUTINE turb_core *** |
---|
| 622 | !! |
---|
| 623 | !! ** Purpose : Computes turbulent transfert coefficients of surface |
---|
| 624 | !! fluxes according to Large & Yeager (2004). |
---|
| 625 | !! |
---|
| 626 | !! ** Method : I N E R T I A L D I S S I P A T I O N M E T H O D |
---|
| 627 | !! Momentum, Latent and sensible heat exchange coefficients |
---|
| 628 | !! Caution: this procedure should only be used in cases when air |
---|
| 629 | !! temperature (T_air) and air specific humidity (q_air) are at 2m |
---|
| 630 | !! whereas wind (dU) is at 10m. |
---|
| 631 | !! |
---|
| 632 | !! References : |
---|
| 633 | !! Large & Yeager, 2004 : ??? |
---|
| 634 | !! History : |
---|
| 635 | !! 9.0 ! 06-12 (L. Brodeau) Original code for 2Z |
---|
| 636 | !!---------------------------------------------------------------------- |
---|
| 637 | !! * Arguments |
---|
| 638 | REAL(wp), INTENT(in) :: & |
---|
| 639 | zt, & ! height for T_zt and q_zt [m] |
---|
| 640 | zu ! height for dU [m] |
---|
| 641 | REAL(wp), INTENT(in), DIMENSION(jpi,jpj) :: & |
---|
| 642 | sst, & ! sea surface temperature [Kelvin] |
---|
| 643 | T_zt, & ! potential air temperature [Kelvin] |
---|
| 644 | q_sat, & ! sea surface specific humidity [kg/kg] |
---|
| 645 | q_zt, & ! specific air humidity [kg/kg] |
---|
| 646 | dU ! relative wind module |U(zu)-U(0)| [m/s] |
---|
| 647 | REAL(wp), INTENT(out), DIMENSION(jpi,jpj) :: & |
---|
| 648 | Cd, & ! transfer coefficient for momentum (tau) |
---|
| 649 | Ch, & ! transfer coefficient for sensible heat (Q_sens) |
---|
| 650 | Ce, & ! transfert coefficient for evaporation (Q_lat) |
---|
| 651 | T_zu, & ! air temp. shifted at zu [K] |
---|
| 652 | q_zu ! spec. hum. shifted at zu [kg/kg] |
---|
| 653 | |
---|
| 654 | !! * Local declarations |
---|
| 655 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
| 656 | dU10, & ! dU [m/s] |
---|
| 657 | dT, & ! air/sea temperature differeence [K] |
---|
| 658 | dq, & ! air/sea humidity difference [K] |
---|
| 659 | Cd_n10, & ! 10m neutral drag coefficient |
---|
| 660 | Ce_n10, & ! 10m neutral latent coefficient |
---|
| 661 | Ch_n10, & ! 10m neutral sensible coefficient |
---|
| 662 | sqrt_Cd_n10, & ! root square of Cd_n10 |
---|
| 663 | sqrt_Cd, & ! root square of Cd |
---|
| 664 | T_vpot_u, & ! virtual potential temperature [K] |
---|
| 665 | T_star, & ! turbulent scale of tem. fluct. |
---|
| 666 | q_star, & ! turbulent humidity of temp. fluct. |
---|
| 667 | U_star, & ! turb. scale of velocity fluct. |
---|
| 668 | L, & ! Monin-Obukov length [m] |
---|
| 669 | zeta_u, & ! stability parameter at height zu |
---|
| 670 | zeta_t, & ! stability parameter at height zt |
---|
| 671 | U_n10, & ! neutral wind velocity at 10m [m] |
---|
| 672 | xlogt, xct, zpsi_hu, zpsi_ht, zpsi_m |
---|
| 673 | |
---|
| 674 | INTEGER :: j_itt |
---|
| 675 | INTEGER, PARAMETER :: nb_itt = 3 ! number of itterations |
---|
| 676 | INTEGER, DIMENSION(jpi,jpj) :: & |
---|
| 677 | & stab ! 1st stability test integer |
---|
| 678 | REAL(wp), PARAMETER :: & |
---|
| 679 | grav = 9.8, & ! gravity |
---|
| 680 | kappa = 0.4 ! von Karman's constant |
---|
| 681 | |
---|
| 682 | !! * Start |
---|
| 683 | |
---|
| 684 | !! Initial air/sea differences |
---|
| 685 | dU10 = max(0.5, dU) ! we don't want to fall under 0.5 m/s |
---|
| 686 | dT = T_zt - sst |
---|
| 687 | dq = q_zt - q_sat |
---|
| 688 | |
---|
| 689 | !! Neutral Drag Coefficient : |
---|
| 690 | stab = 0.5 + sign(0.5,dT) ! stab = 1 if dT > 0 -> STABLE |
---|
| 691 | Cd_n10 = 1E-3*( 2.7/dU10 + 0.142 + dU10/13.09 ) |
---|
| 692 | sqrt_Cd_n10 = sqrt(Cd_n10) |
---|
| 693 | Ce_n10 = 1E-3*( 34.6 * sqrt_Cd_n10 ) |
---|
| 694 | Ch_n10 = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1 - stab)) |
---|
| 695 | |
---|
| 696 | !! Initializing transf. coeff. with their first guess neutral equivalents : |
---|
| 697 | Cd = Cd_n10 ; Ce = Ce_n10 ; Ch = Ch_n10 ; sqrt_Cd = sqrt(Cd) |
---|
| 698 | |
---|
| 699 | !! Initializing z_u values with z_t values : |
---|
| 700 | T_zu = T_zt ; q_zu = q_zt |
---|
| 701 | |
---|
| 702 | !! * Now starting iteration loop |
---|
| 703 | DO j_itt=1, nb_itt |
---|
| 704 | dT = T_zu - sst ; dq = q_zu - q_sat ! Updating air/sea differences |
---|
| 705 | T_vpot_u = T_zu*(1. + 0.608*q_zu) ! Updating virtual potential temperature at zu |
---|
| 706 | U_star = sqrt_Cd*dU10 ! Updating turbulent scales : (L & Y eq. (7)) |
---|
| 707 | T_star = Ch/sqrt_Cd*dT ! |
---|
| 708 | q_star = Ce/sqrt_Cd*dq ! |
---|
| 709 | !! |
---|
| 710 | L = (U_star*U_star) & ! Estimate the Monin-Obukov length at height zu |
---|
| 711 | & / (kappa*grav/T_vpot_u*(T_star*(1.+0.608*q_zu) + 0.608*T_zu*q_star)) |
---|
| 712 | !! Stability parameters : |
---|
| 713 | zeta_u = zu/L ; zeta_u = sign( min(abs(zeta_u),10.0), zeta_u ) |
---|
| 714 | zeta_t = zt/L ; zeta_t = sign( min(abs(zeta_t),10.0), zeta_t ) |
---|
| 715 | zpsi_hu = psi_h(zeta_u) |
---|
| 716 | zpsi_ht = psi_h(zeta_t) |
---|
| 717 | zpsi_m = psi_m(zeta_u) |
---|
| 718 | !! |
---|
| 719 | !! Shifting the wind speed to 10m and neutral stability : (L & Y eq.(9a)) |
---|
| 720 | ! U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u))) |
---|
| 721 | U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) |
---|
| 722 | !! |
---|
| 723 | !! Shifting temperature and humidity at zu : (L & Y eq. (9b-9c)) |
---|
| 724 | ! T_zu = T_zt - T_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t)) |
---|
| 725 | T_zu = T_zt - T_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht) |
---|
| 726 | ! q_zu = q_zt - q_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t)) |
---|
| 727 | q_zu = q_zt - q_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht) |
---|
| 728 | !! |
---|
| 729 | !! q_zu cannot have a negative value : forcing 0 |
---|
| 730 | stab = 0.5 + sign(0.5,q_zu) ; q_zu = stab*q_zu |
---|
| 731 | !! |
---|
| 732 | !! Updating the neutral 10m transfer coefficients : |
---|
| 733 | Cd_n10 = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09) ! L & Y eq. (6a) |
---|
| 734 | sqrt_Cd_n10 = sqrt(Cd_n10) |
---|
| 735 | Ce_n10 = 1E-3 * (34.6 * sqrt_Cd_n10) ! L & Y eq. (6b) |
---|
| 736 | stab = 0.5 + sign(0.5,zeta_u) |
---|
| 737 | Ch_n10 = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c-6d) |
---|
| 738 | !! |
---|
| 739 | !! |
---|
| 740 | !! Shifting the neutral 10m transfer coefficients to (zu,zeta_u) : |
---|
| 741 | ! xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u)) |
---|
| 742 | xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m) |
---|
| 743 | Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd) |
---|
| 744 | !! |
---|
| 745 | ! xlogt = log(zu/10.) - psi_h(zeta_u) |
---|
| 746 | xlogt = log(zu/10.) - zpsi_hu |
---|
| 747 | !! |
---|
| 748 | xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10 |
---|
| 749 | Ch = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct |
---|
| 750 | !! |
---|
| 751 | xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10 |
---|
| 752 | Ce = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct |
---|
| 753 | !! |
---|
| 754 | !! |
---|
| 755 | END DO |
---|
| 756 | !! |
---|
| 757 | END SUBROUTINE TURB_CORE_2Z |
---|
| 758 | |
---|
| 759 | |
---|
[758] | 760 | FUNCTION psi_m(zta) !! Psis, L & Y eq. (8c), (8d), (8e) |
---|
| 761 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta |
---|
[702] | 762 | |
---|
[758] | 763 | REAL(wp), PARAMETER :: pi = 3.141592653589793_wp |
---|
| 764 | REAL(wp), DIMENSION(jpi,jpj) :: psi_m |
---|
| 765 | REAL(wp), DIMENSION(jpi,jpj) :: X2, X, stabit |
---|
| 766 | X2 = sqrt(abs(1. - 16.*zta)) ; X2 = max(X2 , 1.0) ; X = sqrt(X2) |
---|
| 767 | stabit = 0.5 + sign(0.5,zta) |
---|
| 768 | psi_m = -5.*zta*stabit & ! Stable |
---|
| 769 | & + (1. - stabit)*(2*log((1. + X)/2) + log((1. + X2)/2) - 2*atan(X) + pi/2) ! Unstable |
---|
| 770 | END FUNCTION psi_m |
---|
| 771 | |
---|
| 772 | FUNCTION psi_h(zta) !! Psis, L & Y eq. (8c), (8d), (8e) |
---|
| 773 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta |
---|
| 774 | |
---|
| 775 | REAL(wp), DIMENSION(jpi,jpj) :: psi_h |
---|
| 776 | REAL(wp), DIMENSION(jpi,jpj) :: X2, X, stabit |
---|
| 777 | X2 = sqrt(abs(1. - 16.*zta)) ; X2 = max(X2 , 1.) ; X = sqrt(X2) |
---|
| 778 | stabit = 0.5 + sign(0.5,zta) |
---|
| 779 | psi_h = -5.*zta*stabit & ! Stable |
---|
| 780 | & + (1. - stabit)*(2.*log( (1. + X2)/2. )) ! Unstable |
---|
| 781 | END FUNCTION psi_h |
---|
[702] | 782 | |
---|
[758] | 783 | |
---|
[702] | 784 | !!====================================================================== |
---|
| 785 | END MODULE sbcblk_core |
---|