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