1 | !> @file flow_statistics.f90 |
---|
2 | !--------------------------------------------------------------------------------------------------! |
---|
3 | ! This file is part of the PALM model system. |
---|
4 | ! |
---|
5 | ! PALM is free software: you can redistribute it and/or modify it under the terms of the GNU General |
---|
6 | ! Public License as published by the Free Software Foundation, either version 3 of the License, or |
---|
7 | ! (at your option) any later version. |
---|
8 | ! |
---|
9 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the |
---|
10 | ! implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General |
---|
11 | ! Public License for more details. |
---|
12 | ! |
---|
13 | ! You should have received a copy of the GNU General Public License along with PALM. If not, see |
---|
14 | ! <http://www.gnu.org/licenses/>. |
---|
15 | ! |
---|
16 | ! Copyright 1997-2020 Leibniz Universitaet Hannover |
---|
17 | !--------------------------------------------------------------------------------------------------! |
---|
18 | ! |
---|
19 | ! Current revisions: |
---|
20 | ! ------------------ |
---|
21 | ! |
---|
22 | ! |
---|
23 | ! Former revisions: |
---|
24 | ! ----------------- |
---|
25 | ! $Id: flow_statistics.f90 4742 2020-10-14 15:11:02Z suehring $ |
---|
26 | ! Implement snow and graupel (bulk microphysics) |
---|
27 | ! |
---|
28 | ! 4703 2020-09-28 09:21:45Z suehring |
---|
29 | ! Revise averaging of land-surface quantities |
---|
30 | ! |
---|
31 | ! 4672 2020-09-09 21:27:32Z pavelkrc |
---|
32 | ! OpenACC bugfix |
---|
33 | ! |
---|
34 | ! 4671 2020-09-09 20:27:58Z pavelkrc |
---|
35 | ! Implementation of downward facing USM and LSM surfaces |
---|
36 | ! |
---|
37 | ! 4646 2020-08-24 16:02:40Z raasch |
---|
38 | ! file re-formatted to follow the PALM coding standard |
---|
39 | ! |
---|
40 | ! 4581 2020-06-29 08:49:58Z suehring |
---|
41 | ! Formatting adjustment |
---|
42 | ! |
---|
43 | ! 4551 2020-06-02 10:22:25Z suehring |
---|
44 | ! Bugfix in summation for statistical regions |
---|
45 | ! |
---|
46 | ! 4521 2020-05-06 11:39:49Z schwenkel |
---|
47 | ! Rename variable |
---|
48 | ! |
---|
49 | ! 4502 2020-04-17 16:14:16Z schwenkel |
---|
50 | ! Implementation of ice microphysics |
---|
51 | ! |
---|
52 | ! 4472 2020-03-24 12:21:00Z Giersch |
---|
53 | ! Calculations of the Kolmogorov lengt scale eta implemented |
---|
54 | ! |
---|
55 | ! 4464 2020-03-17 11:08:46Z Giersch |
---|
56 | ! Reset last change (r4463) |
---|
57 | ! |
---|
58 | ! 4463 2020-03-17 09:27:36Z Giersch |
---|
59 | ! Calculate horizontally averaged profiles of all velocity components at the same place |
---|
60 | ! |
---|
61 | ! 4444 2020-03-05 15:59:50Z raasch |
---|
62 | ! bugfix: cpp-directives for serial mode added |
---|
63 | ! |
---|
64 | ! 4442 2020-03-04 19:21:13Z suehring |
---|
65 | ! Change order of dimension in surface array %frac to allow for better vectorization. |
---|
66 | ! |
---|
67 | ! 4441 2020-03-04 19:20:35Z suehring |
---|
68 | ! Introduction of wall_flags_total_0, which currently sets bits based on static topography |
---|
69 | ! information used in wall_flags_static_0 |
---|
70 | ! |
---|
71 | ! 4329 2019-12-10 15:46:36Z motisi |
---|
72 | ! Renamed wall_flags_0 to wall_flags_static_0 |
---|
73 | ! |
---|
74 | ! 4182 2019-08-22 15:20:23Z scharf |
---|
75 | ! Corrected "Former revisions" section |
---|
76 | ! |
---|
77 | ! 4131 2019-08-02 11:06:18Z monakurppa |
---|
78 | ! Allow profile output for salsa variables. |
---|
79 | ! |
---|
80 | ! 4039 2019-06-18 10:32:41Z suehring |
---|
81 | ! Correct conversion to kinematic scalar fluxes in case of pw-scheme and statistic regions |
---|
82 | ! |
---|
83 | ! 3828 2019-03-27 19:36:23Z raasch |
---|
84 | ! unused variables removed |
---|
85 | ! |
---|
86 | ! 3676 2019-01-16 15:07:05Z knoop |
---|
87 | ! Bugfix, terminate OMP Parallel block |
---|
88 | ! |
---|
89 | ! Revision 1.1 1997/08/11 06:15:17 raasch |
---|
90 | ! Initial revision |
---|
91 | ! |
---|
92 | ! |
---|
93 | ! Description: |
---|
94 | ! ------------ |
---|
95 | !> Compute average profiles and further average flow quantities for the different user-defined |
---|
96 | !> (sub-)regions. The region indexed 0 is the total model domain. |
---|
97 | !> |
---|
98 | !> @note For simplicity, nzb_s_inner and nzb_diff_s_inner are used as a lower vertical index for |
---|
99 | !> k-loops for all variables, although strictly speaking the k-loops would have to be split |
---|
100 | !> up according to the staggered grid. However, this implies no error since staggered velocity |
---|
101 | !> components are zero at the walls and inside buildings. |
---|
102 | !--------------------------------------------------------------------------------------------------! |
---|
103 | SUBROUTINE flow_statistics |
---|
104 | |
---|
105 | |
---|
106 | USE arrays_3d, & |
---|
107 | ONLY: ddzu, ddzw, d_exner, e, heatflux_output_conversion, hyp, km, kh, & |
---|
108 | momentumflux_output_conversion, & |
---|
109 | nc, ni, ng, nr, ns, p, prho, prr, pt, q, qc, qi, qg, ql, qr, qs, & |
---|
110 | rho_air, rho_air_zw, rho_ocean, s, sa, u, ug, v, vg, vpt, w, w_subs, & |
---|
111 | waterflux_output_conversion, zw |
---|
112 | |
---|
113 | USE basic_constants_and_equations_mod, & |
---|
114 | ONLY: g, lv_d_cp |
---|
115 | |
---|
116 | USE bulk_cloud_model_mod, & |
---|
117 | ONLY: bulk_cloud_model, graupel, snow, microphysics_morrison, microphysics_seifert, & |
---|
118 | microphysics_ice_phase |
---|
119 | |
---|
120 | USE chem_modules, & |
---|
121 | ONLY: max_pr_cs |
---|
122 | |
---|
123 | USE control_parameters, & |
---|
124 | ONLY: air_chemistry, average_count_pr, cloud_droplets, do_sum, dt_3d, humidity, & |
---|
125 | initializing_actions, kolmogorov_length_scale, land_surface, large_scale_forcing, & |
---|
126 | large_scale_subsidence, max_pr_salsa, max_pr_user, message_string, neutral, & |
---|
127 | ocean_mode, passive_scalar, salsa, simulated_time, simulated_time_at_begin, & |
---|
128 | use_subsidence_tendencies, use_surface_fluxes, use_top_fluxes, ws_scheme_mom, & |
---|
129 | ws_scheme_sca |
---|
130 | |
---|
131 | USE cpulog, & |
---|
132 | ONLY: cpu_log, log_point |
---|
133 | |
---|
134 | USE grid_variables, & |
---|
135 | ONLY: ddx, ddy |
---|
136 | |
---|
137 | USE indices, & |
---|
138 | ONLY: ngp_2dh, ngp_2dh_s_inner, ngp_3d, ngp_3d_inner, nxl, nxr, nyn, nys, nzb, nzt, & |
---|
139 | topo_min_level, wall_flags_total_0 |
---|
140 | |
---|
141 | #if defined( __parallel ) |
---|
142 | USE indices, & |
---|
143 | ONLY: ngp_sums, ngp_sums_ls |
---|
144 | #endif |
---|
145 | |
---|
146 | USE kinds |
---|
147 | |
---|
148 | USE land_surface_model_mod, & |
---|
149 | ONLY: m_soil_h, nzb_soil, nzt_soil, t_soil_h |
---|
150 | |
---|
151 | USE lsf_nudging_mod, & |
---|
152 | ONLY: td_lsa_lpt, td_lsa_q, td_sub_lpt, td_sub_q, time_vert |
---|
153 | |
---|
154 | USE module_interface, & |
---|
155 | ONLY: module_interface_statistics |
---|
156 | |
---|
157 | USE netcdf_interface, & |
---|
158 | ONLY: dots_rad, dots_soil, dots_max |
---|
159 | |
---|
160 | USE pegrid |
---|
161 | |
---|
162 | USE radiation_model_mod, & |
---|
163 | ONLY: radiation, radiation_scheme, & |
---|
164 | rad_lw_in, rad_lw_out, rad_lw_cs_hr, rad_lw_hr, & |
---|
165 | rad_sw_in, rad_sw_out, rad_sw_cs_hr, rad_sw_hr |
---|
166 | |
---|
167 | USE statistics |
---|
168 | |
---|
169 | USE surface_mod, & |
---|
170 | ONLY : surf_def_h, surf_lsm_h, surf_usm_h |
---|
171 | |
---|
172 | |
---|
173 | IMPLICIT NONE |
---|
174 | |
---|
175 | INTEGER(iwp) :: i !< |
---|
176 | INTEGER(iwp) :: j !< |
---|
177 | INTEGER(iwp) :: k !< |
---|
178 | INTEGER(iwp) :: ki !< |
---|
179 | INTEGER(iwp) :: k_surface_level !< |
---|
180 | INTEGER(iwp) :: m !< loop variable over all horizontal wall elements |
---|
181 | INTEGER(iwp) :: l !< loop variable over surface facing -- up- or downward-facing |
---|
182 | INTEGER(iwp) :: nt !< |
---|
183 | !$ INTEGER(iwp) :: omp_get_thread_num !< |
---|
184 | INTEGER(iwp) :: sr !< |
---|
185 | INTEGER(iwp) :: tn !< |
---|
186 | |
---|
187 | LOGICAL :: first !< |
---|
188 | |
---|
189 | REAL(wp) :: dissipation !< dissipation rate |
---|
190 | REAL(wp) :: dptdz_threshold !< |
---|
191 | REAL(wp) :: du_dx !< Derivative of u fluctuations with respect to x |
---|
192 | REAL(wp) :: du_dy !< Derivative of u fluctuations with respect to y |
---|
193 | REAL(wp) :: du_dz !< Derivative of u fluctuations with respect to z |
---|
194 | REAL(wp) :: dv_dx !< Derivative of v fluctuations with respect to x |
---|
195 | REAL(wp) :: dv_dy !< Derivative of v fluctuations with respect to y |
---|
196 | REAL(wp) :: dv_dz !< Derivative of v fluctuations with respect to z |
---|
197 | REAL(wp) :: dw_dx !< Derivative of w fluctuations with respect to x |
---|
198 | REAL(wp) :: dw_dy !< Derivative of w fluctuations with respect to y |
---|
199 | REAL(wp) :: dw_dz !< Derivative of w fluctuations with respect to z |
---|
200 | REAL(wp) :: eta !< Kolmogorov length scale |
---|
201 | REAL(wp) :: fac !< |
---|
202 | REAL(wp) :: flag !< |
---|
203 | REAL(wp) :: height !< |
---|
204 | REAL(wp) :: pts !< |
---|
205 | REAL(wp) :: s11 !< fluctuating rate-of-strain tensor component 11 |
---|
206 | REAL(wp) :: s21 !< fluctuating rate-of-strain tensor component 21 |
---|
207 | REAL(wp) :: s31 !< fluctuating rate-of-strain tensor component 31 |
---|
208 | REAL(wp) :: s12 !< fluctuating rate-of-strain tensor component 12 |
---|
209 | REAL(wp) :: s22 !< fluctuating rate-of-strain tensor component 22 |
---|
210 | REAL(wp) :: s32 !< fluctuating rate-of-strain tensor component 32 |
---|
211 | REAL(wp) :: s13 !< fluctuating rate-of-strain tensor component 13 |
---|
212 | REAL(wp) :: s23 !< fluctuating rate-of-strain tensor component 23 |
---|
213 | REAL(wp) :: s33 !< fluctuating rate-of-strain tensor component 33 |
---|
214 | REAL(wp) :: sums_l_etot !< |
---|
215 | REAL(wp) :: ust !< |
---|
216 | REAL(wp) :: ust2 !< |
---|
217 | REAL(wp) :: u2 !< |
---|
218 | REAL(wp) :: vst !< |
---|
219 | REAL(wp) :: vst2 !< |
---|
220 | REAL(wp) :: v2 !< |
---|
221 | REAL(wp) :: w2 !< |
---|
222 | |
---|
223 | REAL(wp) :: dptdz(nzb+1:nzt+1) !< |
---|
224 | REAL(wp) :: sums_ll(nzb:nzt+1,2) !< |
---|
225 | |
---|
226 | |
---|
227 | CALL cpu_log( log_point(10), 'flow_statistics', 'start' ) |
---|
228 | |
---|
229 | |
---|
230 | ! |
---|
231 | !-- To be on the safe side, check whether flow_statistics has already been called once after the |
---|
232 | !-- current time step. |
---|
233 | IF ( flow_statistics_called ) THEN |
---|
234 | |
---|
235 | message_string = 'flow_statistics is called two times within one ' // 'timestep' |
---|
236 | CALL message( 'flow_statistics', 'PA0190', 1, 2, 0, 6, 0 ) |
---|
237 | |
---|
238 | ENDIF |
---|
239 | |
---|
240 | ! |
---|
241 | !-- Compute statistics for each (sub-)region |
---|
242 | DO sr = 0, statistic_regions |
---|
243 | |
---|
244 | ! |
---|
245 | !-- Initialize (local) summation array |
---|
246 | sums_l = 0.0_wp |
---|
247 | #ifdef _OPENACC |
---|
248 | !$ACC KERNELS PRESENT(sums_l) |
---|
249 | sums_l = 0.0_wp |
---|
250 | !$ACC END KERNELS |
---|
251 | #endif |
---|
252 | |
---|
253 | ! |
---|
254 | !-- Store sums that have been computed in other subroutines in summation array |
---|
255 | sums_l(:,11,:) = sums_l_l(:,sr,:) ! mixing length from diffusivities |
---|
256 | !-- WARNING: next line still has to be adjusted for OpenMP |
---|
257 | sums_l(:,21,0) = sums_wsts_bc_l(:,sr) * & |
---|
258 | heatflux_output_conversion ! heat flux from advec_s_bc |
---|
259 | sums_l(nzb+9,pr_palm,0) = sums_divold_l(sr) ! old divergence from pres |
---|
260 | sums_l(nzb+10,pr_palm,0) = sums_divnew_l(sr) ! new divergence from pres |
---|
261 | |
---|
262 | ! |
---|
263 | !-- When calcuating horizontally-averaged total (resolved- plus subgrid-scale) vertical fluxes |
---|
264 | !-- and velocity variances by using commonly-applied Reynolds-based methods |
---|
265 | !-- ( e.g. <w'pt'> = (w-<w>)*(pt-<pt>) ) in combination with the 5th order advection scheme, |
---|
266 | !-- pronounced artificial kinks could be observed in the vertical profiles near the surface. |
---|
267 | !-- Please note: these kinks were not related to the model truth, i.e. these kinks are just |
---|
268 | !-- related to an evaluation problem. |
---|
269 | !-- In order avoid these kinks, vertical fluxes and horizontal as well vertical velocity |
---|
270 | !-- variances are calculated directly within the advection routines, according to the numerical |
---|
271 | !-- discretization, to evaluate the statistical quantities as they will appear within the |
---|
272 | !-- prognostic equations. |
---|
273 | !-- Copy the turbulent quantities, evaluated in the advection routines to the local array |
---|
274 | !-- sums_l() for further computations. |
---|
275 | IF ( ws_scheme_mom .AND. sr == 0 ) THEN |
---|
276 | |
---|
277 | ! |
---|
278 | !-- According to the Neumann bc for the horizontal velocity components, the corresponding |
---|
279 | !-- fluxes has to satisfiy the same bc. |
---|
280 | IF ( ocean_mode ) THEN |
---|
281 | sums_us2_ws_l(nzt+1,:) = sums_us2_ws_l(nzt,:) |
---|
282 | sums_vs2_ws_l(nzt+1,:) = sums_vs2_ws_l(nzt,:) |
---|
283 | ENDIF |
---|
284 | |
---|
285 | DO i = 0, threads_per_task-1 |
---|
286 | ! |
---|
287 | !-- Swap the turbulent quantities evaluated in advec_ws. |
---|
288 | sums_l(:,13,i) = sums_wsus_ws_l(:,i) * momentumflux_output_conversion ! w*u* |
---|
289 | sums_l(:,15,i) = sums_wsvs_ws_l(:,i) * momentumflux_output_conversion ! w*v* |
---|
290 | sums_l(:,30,i) = sums_us2_ws_l(:,i) ! u*2 |
---|
291 | sums_l(:,31,i) = sums_vs2_ws_l(:,i) ! v*2 |
---|
292 | sums_l(:,32,i) = sums_ws2_ws_l(:,i) ! w*2 |
---|
293 | sums_l(:,34,i) = sums_l(:,34,i) + 0.5_wp * & |
---|
294 | ( sums_us2_ws_l(:,i) + sums_vs2_ws_l(:,i) + sums_ws2_ws_l(:,i) ) ! e* |
---|
295 | ENDDO |
---|
296 | |
---|
297 | ENDIF |
---|
298 | |
---|
299 | IF ( ws_scheme_sca .AND. sr == 0 ) THEN |
---|
300 | |
---|
301 | DO i = 0, threads_per_task-1 |
---|
302 | sums_l(:,17,i) = sums_wspts_ws_l(:,i) & |
---|
303 | * heatflux_output_conversion ! w*pt* |
---|
304 | IF ( ocean_mode ) sums_l(:,66,i) = sums_wssas_ws_l(:,i) ! w*sa* |
---|
305 | IF ( humidity ) sums_l(:,49,i) = sums_wsqs_ws_l(:,i) & |
---|
306 | * waterflux_output_conversion ! w*q* |
---|
307 | IF ( passive_scalar ) sums_l(:,114,i) = sums_wsss_ws_l(:,i) ! w*s* |
---|
308 | ENDDO |
---|
309 | |
---|
310 | ENDIF |
---|
311 | ! |
---|
312 | !-- Horizontally averaged profiles of horizontal velocities and temperature. |
---|
313 | !-- They must have been computed before, because they are already required for other horizontal |
---|
314 | !-- averages. |
---|
315 | tn = 0 |
---|
316 | !$OMP PARALLEL PRIVATE( i, j, k, tn, flag ) |
---|
317 | !$ tn = omp_get_thread_num() |
---|
318 | !$OMP DO |
---|
319 | !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k, flag) & |
---|
320 | !$ACC PRESENT(wall_flags_total_0, u, v, pt, rmask, sums_l) |
---|
321 | DO i = nxl, nxr |
---|
322 | DO j = nys, nyn |
---|
323 | DO k = nzb, nzt+1 |
---|
324 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
325 | !$ACC ATOMIC |
---|
326 | sums_l(k,1,tn) = sums_l(k,1,tn) + u(k,j,i) * rmask(j,i,sr) * flag |
---|
327 | !$ACC ATOMIC |
---|
328 | sums_l(k,2,tn) = sums_l(k,2,tn) + v(k,j,i) * rmask(j,i,sr) * flag |
---|
329 | !$ACC ATOMIC |
---|
330 | sums_l(k,4,tn) = sums_l(k,4,tn) + pt(k,j,i) * rmask(j,i,sr) * flag |
---|
331 | ENDDO |
---|
332 | ENDDO |
---|
333 | ENDDO |
---|
334 | !$ACC UPDATE HOST(sums_l(:,1,tn), sums_l(:,2,tn), sums_l(:,4,tn)) |
---|
335 | |
---|
336 | ! |
---|
337 | !-- Horizontally averaged profile of salinity |
---|
338 | IF ( ocean_mode ) THEN |
---|
339 | !$OMP DO |
---|
340 | DO i = nxl, nxr |
---|
341 | DO j = nys, nyn |
---|
342 | DO k = nzb, nzt+1 |
---|
343 | sums_l(k,23,tn) = sums_l(k,23,tn) + sa(k,j,i) & |
---|
344 | * rmask(j,i,sr) & |
---|
345 | * MERGE( 1.0_wp, 0.0_wp, & |
---|
346 | BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
347 | ENDDO |
---|
348 | ENDDO |
---|
349 | ENDDO |
---|
350 | ENDIF |
---|
351 | |
---|
352 | ! |
---|
353 | !-- Horizontally averaged profiles of virtual potential temperature, total water content, water |
---|
354 | !-- vapor mixing ratio and liquid water potential temperature |
---|
355 | IF ( humidity ) THEN |
---|
356 | !$OMP DO |
---|
357 | DO i = nxl, nxr |
---|
358 | DO j = nys, nyn |
---|
359 | DO k = nzb, nzt+1 |
---|
360 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
361 | sums_l(k,44,tn) = sums_l(k,44,tn) + vpt(k,j,i) * rmask(j,i,sr) * flag |
---|
362 | sums_l(k,41,tn) = sums_l(k,41,tn) + q(k,j,i) * rmask(j,i,sr) * flag |
---|
363 | ENDDO |
---|
364 | ENDDO |
---|
365 | ENDDO |
---|
366 | IF ( bulk_cloud_model ) THEN |
---|
367 | !$OMP DO |
---|
368 | DO i = nxl, nxr |
---|
369 | DO j = nys, nyn |
---|
370 | DO k = nzb, nzt+1 |
---|
371 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
372 | sums_l(k,42,tn) = sums_l(k,42,tn) + & |
---|
373 | ( q(k,j,i) - ql(k,j,i) ) * rmask(j,i,sr) * flag |
---|
374 | sums_l(k,43,tn) = sums_l(k,43,tn) + ( & |
---|
375 | pt(k,j,i) + lv_d_cp * d_exner(k) * ql(k,j,i) & |
---|
376 | ) * rmask(j,i,sr) * flag |
---|
377 | ENDDO |
---|
378 | ENDDO |
---|
379 | ENDDO |
---|
380 | ENDIF |
---|
381 | ENDIF |
---|
382 | |
---|
383 | ! |
---|
384 | !-- Horizontally averaged profiles of passive scalar |
---|
385 | IF ( passive_scalar ) THEN |
---|
386 | !$OMP DO |
---|
387 | DO i = nxl, nxr |
---|
388 | DO j = nys, nyn |
---|
389 | DO k = nzb, nzt+1 |
---|
390 | sums_l(k,115,tn) = sums_l(k,115,tn) + s(k,j,i) & |
---|
391 | * rmask(j,i,sr) & |
---|
392 | * MERGE( 1.0_wp, 0.0_wp, & |
---|
393 | BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
394 | ENDDO |
---|
395 | ENDDO |
---|
396 | ENDDO |
---|
397 | ENDIF |
---|
398 | !$OMP END PARALLEL |
---|
399 | ! |
---|
400 | !-- Summation of thread sums |
---|
401 | IF ( threads_per_task > 1 ) THEN |
---|
402 | DO i = 1, threads_per_task-1 |
---|
403 | sums_l(:,1,0) = sums_l(:,1,0) + sums_l(:,1,i) |
---|
404 | sums_l(:,2,0) = sums_l(:,2,0) + sums_l(:,2,i) |
---|
405 | sums_l(:,4,0) = sums_l(:,4,0) + sums_l(:,4,i) |
---|
406 | IF ( ocean_mode ) THEN |
---|
407 | sums_l(:,23,0) = sums_l(:,23,0) + sums_l(:,23,i) |
---|
408 | ENDIF |
---|
409 | IF ( humidity ) THEN |
---|
410 | sums_l(:,41,0) = sums_l(:,41,0) + sums_l(:,41,i) |
---|
411 | sums_l(:,44,0) = sums_l(:,44,0) + sums_l(:,44,i) |
---|
412 | IF ( bulk_cloud_model ) THEN |
---|
413 | sums_l(:,42,0) = sums_l(:,42,0) + sums_l(:,42,i) |
---|
414 | sums_l(:,43,0) = sums_l(:,43,0) + sums_l(:,43,i) |
---|
415 | ENDIF |
---|
416 | ENDIF |
---|
417 | IF ( passive_scalar ) THEN |
---|
418 | sums_l(:,115,0) = sums_l(:,115,0) + sums_l(:,115,i) |
---|
419 | ENDIF |
---|
420 | ENDDO |
---|
421 | ENDIF |
---|
422 | |
---|
423 | #if defined( __parallel ) |
---|
424 | ! |
---|
425 | !-- Compute total sum from local sums |
---|
426 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
427 | CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d, & |
---|
428 | ierr ) |
---|
429 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
430 | CALL MPI_ALLREDUCE( sums_l(nzb,2,0), sums(nzb,2), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d, & |
---|
431 | ierr ) |
---|
432 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
433 | CALL MPI_ALLREDUCE( sums_l(nzb,4,0), sums(nzb,4), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d, & |
---|
434 | ierr ) |
---|
435 | IF ( ocean_mode ) THEN |
---|
436 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
437 | CALL MPI_ALLREDUCE( sums_l(nzb,23,0), sums(nzb,23), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d,& |
---|
438 | ierr ) |
---|
439 | ENDIF |
---|
440 | IF ( humidity ) THEN |
---|
441 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
442 | CALL MPI_ALLREDUCE( sums_l(nzb,44,0), sums(nzb,44), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d,& |
---|
443 | ierr ) |
---|
444 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
445 | CALL MPI_ALLREDUCE( sums_l(nzb,41,0), sums(nzb,41), nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d,& |
---|
446 | ierr ) |
---|
447 | IF ( bulk_cloud_model ) THEN |
---|
448 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
449 | CALL MPI_ALLREDUCE( sums_l(nzb,42,0), sums(nzb,42), nzt+2-nzb, MPI_REAL, MPI_SUM, & |
---|
450 | comm2d, ierr ) |
---|
451 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
452 | CALL MPI_ALLREDUCE( sums_l(nzb,43,0), sums(nzb,43), nzt+2-nzb, MPI_REAL, MPI_SUM, & |
---|
453 | comm2d, ierr ) |
---|
454 | ENDIF |
---|
455 | ENDIF |
---|
456 | |
---|
457 | IF ( passive_scalar ) THEN |
---|
458 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
459 | CALL MPI_ALLREDUCE( sums_l(nzb,115,0), sums(nzb,115), nzt+2-nzb, MPI_REAL, MPI_SUM, & |
---|
460 | comm2d, ierr ) |
---|
461 | ENDIF |
---|
462 | #else |
---|
463 | sums(:,1) = sums_l(:,1,0) |
---|
464 | sums(:,2) = sums_l(:,2,0) |
---|
465 | sums(:,4) = sums_l(:,4,0) |
---|
466 | IF ( ocean_mode ) sums(:,23) = sums_l(:,23,0) |
---|
467 | IF ( humidity ) THEN |
---|
468 | sums(:,44) = sums_l(:,44,0) |
---|
469 | sums(:,41) = sums_l(:,41,0) |
---|
470 | IF ( bulk_cloud_model ) THEN |
---|
471 | sums(:,42) = sums_l(:,42,0) |
---|
472 | sums(:,43) = sums_l(:,43,0) |
---|
473 | ENDIF |
---|
474 | ENDIF |
---|
475 | IF ( passive_scalar ) sums(:,115) = sums_l(:,115,0) |
---|
476 | #endif |
---|
477 | |
---|
478 | ! |
---|
479 | !-- Final values are obtained by division by the total number of grid points used for summation. |
---|
480 | !-- After that store profiles. |
---|
481 | sums(:,1) = sums(:,1) / ngp_2dh(sr) |
---|
482 | sums(:,2) = sums(:,2) / ngp_2dh(sr) |
---|
483 | sums(:,4) = sums(:,4) / ngp_2dh_s_inner(:,sr) |
---|
484 | hom(:,1,1,sr) = sums(:,1) ! u |
---|
485 | hom(:,1,2,sr) = sums(:,2) ! v |
---|
486 | hom(:,1,4,sr) = sums(:,4) ! pt |
---|
487 | !$ACC UPDATE DEVICE(hom(:,1,1,sr), hom(:,1,2,sr), hom(:,1,4,sr)) |
---|
488 | |
---|
489 | |
---|
490 | ! |
---|
491 | !-- Salinity |
---|
492 | IF ( ocean_mode ) THEN |
---|
493 | sums(:,23) = sums(:,23) / ngp_2dh_s_inner(:,sr) |
---|
494 | hom(:,1,23,sr) = sums(:,23) ! sa |
---|
495 | ENDIF |
---|
496 | |
---|
497 | ! |
---|
498 | !-- Humidity and cloud parameters |
---|
499 | IF ( humidity ) THEN |
---|
500 | sums(:,44) = sums(:,44) / ngp_2dh_s_inner(:,sr) |
---|
501 | sums(:,41) = sums(:,41) / ngp_2dh_s_inner(:,sr) |
---|
502 | hom(:,1,44,sr) = sums(:,44) ! vpt |
---|
503 | hom(:,1,41,sr) = sums(:,41) ! qv (q) |
---|
504 | IF ( bulk_cloud_model ) THEN |
---|
505 | sums(:,42) = sums(:,42) / ngp_2dh_s_inner(:,sr) |
---|
506 | sums(:,43) = sums(:,43) / ngp_2dh_s_inner(:,sr) |
---|
507 | hom(:,1,42,sr) = sums(:,42) ! qv |
---|
508 | hom(:,1,43,sr) = sums(:,43) ! pt |
---|
509 | ENDIF |
---|
510 | ENDIF |
---|
511 | |
---|
512 | ! |
---|
513 | !-- Passive scalar |
---|
514 | IF ( passive_scalar ) hom(:,1,115,sr) = sums(:,115) / ngp_2dh_s_inner(:,sr) ! s |
---|
515 | |
---|
516 | ! |
---|
517 | !-- Horizontally averaged profiles of the remaining prognostic variables, variances, the total |
---|
518 | !-- and the perturbation energy (single values in last column of sums_l) and some diagnostic |
---|
519 | !-- quantities. |
---|
520 | !-- NOTE: for simplicity, nzb_s_inner is used below, although strictly speaking the following |
---|
521 | !-- ---- k-loop would have to be split up and rearranged according to the staggered grid. |
---|
522 | !-- However, this implies no error since staggered velocity components are zero at the |
---|
523 | !-- walls and inside buildings. |
---|
524 | tn = 0 |
---|
525 | !$OMP PARALLEL PRIVATE( i, j, k, pts, sums_ll, & |
---|
526 | !$OMP sums_l_etot, tn, ust, ust2, u2, vst, vst2, v2, & |
---|
527 | !$OMP w2, flag, m, ki, l ) |
---|
528 | !$ tn = omp_get_thread_num() |
---|
529 | !$OMP DO |
---|
530 | !$ACC PARALLEL LOOP COLLAPSE(2) PRIVATE(i, j, k, m) & |
---|
531 | !$ACC PRIVATE(sums_l_etot, flag, du_dx, du_dy, du_dz) & |
---|
532 | !$ACC PRIVATE(dv_dx, dv_dy, dv_dz, dw_dx, dw_dy, dw_dz) & |
---|
533 | !$ACC PRIVATE(s11, s21, s31, s12, s22, s32, s13, s23, s33) & |
---|
534 | !$ACC PRIVATE(dissipation, eta) & |
---|
535 | !$ACC PRESENT(wall_flags_total_0, rmask, momentumflux_output_conversion) & |
---|
536 | !$ACC PRESENT(hom(:,1,1:2,sr), hom(:,1,4,sr)) & |
---|
537 | !$ACC PRESENT(e, u, v, w, km, kh, p, pt) & |
---|
538 | !$ACC PRESENT(ddx, ddy, ddzu, ddzw) & |
---|
539 | !$ACC PRESENT(surf_def_h(0), surf_lsm_h(0), surf_usm_h(0)) & |
---|
540 | !$ACC PRESENT(sums_l) |
---|
541 | DO i = nxl, nxr |
---|
542 | DO j = nys, nyn |
---|
543 | sums_l_etot = 0.0_wp |
---|
544 | DO k = nzb, nzt+1 |
---|
545 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
546 | ! |
---|
547 | !-- Prognostic and diagnostic variables |
---|
548 | !$ACC ATOMIC |
---|
549 | sums_l(k,3,tn) = sums_l(k,3,tn) + w(k,j,i) * rmask(j,i,sr) * flag |
---|
550 | !$ACC ATOMIC |
---|
551 | sums_l(k,8,tn) = sums_l(k,8,tn) + e(k,j,i) * rmask(j,i,sr) * flag |
---|
552 | !$ACC ATOMIC |
---|
553 | sums_l(k,9,tn) = sums_l(k,9,tn) + km(k,j,i) * rmask(j,i,sr) * flag |
---|
554 | !$ACC ATOMIC |
---|
555 | sums_l(k,10,tn) = sums_l(k,10,tn) + kh(k,j,i) * rmask(j,i,sr) * flag |
---|
556 | !$ACC ATOMIC |
---|
557 | sums_l(k,40,tn) = sums_l(k,40,tn) + ( p(k,j,i) & |
---|
558 | / momentumflux_output_conversion(k) ) * flag |
---|
559 | |
---|
560 | !$ACC ATOMIC |
---|
561 | sums_l(k,33,tn) = sums_l(k,33,tn) + & |
---|
562 | ( pt(k,j,i)-hom(k,1,4,sr) )**2 * rmask(j,i,sr) * flag |
---|
563 | #ifndef _OPENACC |
---|
564 | IF ( humidity ) THEN |
---|
565 | sums_l(k,70,tn) = sums_l(k,70,tn) + & |
---|
566 | ( q(k,j,i)-hom(k,1,41,sr) )**2 * rmask(j,i,sr) * flag |
---|
567 | ENDIF |
---|
568 | IF ( passive_scalar ) THEN |
---|
569 | sums_l(k,116,tn) = sums_l(k,116,tn) + & |
---|
570 | ( s(k,j,i)-hom(k,1,115,sr) )**2 * rmask(j,i,sr) * flag |
---|
571 | ENDIF |
---|
572 | #endif |
---|
573 | ! |
---|
574 | !-- Higher moments |
---|
575 | !-- (Computation of the skewness of w further below) |
---|
576 | !$ACC ATOMIC |
---|
577 | sums_l(k,38,tn) = sums_l(k,38,tn) + w(k,j,i)**3 * rmask(j,i,sr) * flag |
---|
578 | |
---|
579 | sums_l_etot = sums_l_etot + 0.5_wp * ( u(k,j,i)**2 + v(k,j,i)**2 + w(k,j,i)**2 ) & |
---|
580 | * rmask(j,i,sr) * flag |
---|
581 | |
---|
582 | ! |
---|
583 | !-- Computation of the Kolmogorov length scale. Calculation is based on gradients of the |
---|
584 | !-- deviations from the horizontal mean. |
---|
585 | !-- Kolmogorov scale at the boundaries (k=0/z=0m and k=nzt+1) is set to zero. |
---|
586 | IF ( kolmogorov_length_scale .AND. k /= nzb .AND. k /= nzt+1) THEN |
---|
587 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
588 | |
---|
589 | ! |
---|
590 | !-- Calculate components of the fluctuating rate-of-strain tensor |
---|
591 | !-- (0.5*(del u'_i/del x_j + del u'_j/del x_i)) defined in the center of each grid |
---|
592 | !-- box. |
---|
593 | du_dx = ( ( u(k,j,i+1) - hom(k,1,1,sr) ) - & |
---|
594 | ( u(k,j,i) - hom(k,1,1,sr) ) ) * ddx |
---|
595 | du_dy = 0.25_wp * ddy * & |
---|
596 | ( ( u(k,j+1,i) - hom(k,1,1,sr) ) - & |
---|
597 | ( u(k,j-1,i) - hom(k,1,1,sr) ) + & |
---|
598 | ( u(k,j+1,i+1) - hom(k,1,1,sr) ) - & |
---|
599 | ( u(k,j-1,i+1) - hom(k,1,1,sr) ) ) |
---|
600 | du_dz = 0.25_wp * ( ( ( u(k+1,j,i) - hom(k+1,1,1,sr) ) - & |
---|
601 | ( u(k,j,i) - hom(k,1,1,sr) ) ) * & |
---|
602 | ddzu(k+1) + & |
---|
603 | ( ( u(k,j,i) - hom(k,1,1,sr) ) - & |
---|
604 | ( u(k-1,j,i) - hom(k-1,1,1,sr) ) ) * & |
---|
605 | ddzu(k) + & |
---|
606 | ( ( u(k+1,j,i+1) - hom(k+1,1,1,sr) ) - & |
---|
607 | ( u(k,j,i+1) - hom(k,1,1,sr) ) ) * & |
---|
608 | ddzu(k+1) + & |
---|
609 | ( ( u(k,j,i+1) - hom(k,1,1,sr) ) - & |
---|
610 | ( u(k-1,j,i+1) - hom(k-1,1,1,sr) ) ) * & |
---|
611 | ddzu(k) ) |
---|
612 | |
---|
613 | dv_dx = 0.25_wp * ddx * & |
---|
614 | ( ( v(k,j,i+1) - hom(k,1,2,sr) ) - & |
---|
615 | ( v(k,j,i-1) - hom(k,1,2,sr) ) + & |
---|
616 | ( v(k,j+1,i+1) - hom(k,1,2,sr) ) - & |
---|
617 | ( v(k,j+1,i-1) - hom(k,1,2,sr) ) ) |
---|
618 | dv_dy = ( ( v(k,j+1,i) - hom(k,1,2,sr) ) - ( v(k,j,i) - hom(k,1,2,sr) ) ) * ddy |
---|
619 | dv_dz = 0.25_wp * ( ( ( v(k+1,j,i) - hom(k+1,1,2,sr) ) - & |
---|
620 | ( v(k,j,i) - hom(k,1,2,sr) ) ) * & |
---|
621 | ddzu(k+1) + & |
---|
622 | ( ( v(k,j,i) - hom(k,1,2,sr) ) - & |
---|
623 | ( v(k-1,j,i) - hom(k-1,1,2,sr) ) ) * & |
---|
624 | ddzu(k) + & |
---|
625 | ( ( v(k+1,j+1,i) - hom(k+1,1,2,sr) ) - & |
---|
626 | ( v(k,j+1,i) - hom(k,1,2,sr) ) ) * & |
---|
627 | ddzu(k+1) + & |
---|
628 | ( ( v(k,j+1,i) - hom(k,1,2,sr) ) - & |
---|
629 | ( v(k-1,j+1,i) - hom(k-1,1,2,sr) ) ) * & |
---|
630 | ddzu(k) ) |
---|
631 | |
---|
632 | dw_dx = 0.25_wp * ddx * ( w(k,j,i+1) - w(k,j,i-1) + w(k-1,j,i+1) - w(k-1,j,i-1) ) |
---|
633 | dw_dy = 0.25_wp * ddy * ( w(k,j+1,i) - w(k,j-1,i) + w(k-1,j+1,i) - w(k-1,j-1,i) ) |
---|
634 | dw_dz = ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) |
---|
635 | |
---|
636 | s11 = 0.5_wp * ( du_dx + du_dx ) |
---|
637 | s21 = 0.5_wp * ( dv_dx + du_dy ) |
---|
638 | s31 = 0.5_wp * ( dw_dx + du_dz ) |
---|
639 | |
---|
640 | s12 = s21 |
---|
641 | s22 = 0.5 * ( dv_dy + dv_dy ) |
---|
642 | s32 = 0.5 * ( dw_dy + dv_dz ) |
---|
643 | |
---|
644 | s13 = s31 |
---|
645 | s23 = s32 |
---|
646 | s33 = 0.5_wp * ( dw_dz + dw_dz ) |
---|
647 | |
---|
648 | !-- Calculate 3D instantaneous energy dissipation rate following Pope (2000): |
---|
649 | !-- Turbulent flows, p.259. It is defined in the center of each grid volume. |
---|
650 | dissipation = 2.0_wp * km(k,j,i) * & |
---|
651 | ( s11*s11 + s21*s21 + s31*s31 + & |
---|
652 | s12*s12 + s22*s22 + s32*s32 + & |
---|
653 | s13*s13 + s23*s23 + s33*s33 ) |
---|
654 | eta = ( km(k,j,i)**3.0_wp / ( dissipation+1.0E-12 ) )**(1.0_wp/4.0_wp) |
---|
655 | |
---|
656 | !$ACC ATOMIC |
---|
657 | sums_l(k,121,tn) = sums_l(k,121,tn) + eta * rmask(j,i,sr) * flag |
---|
658 | |
---|
659 | |
---|
660 | ENDIF !Kolmogorov length scale |
---|
661 | |
---|
662 | ENDDO !k-loop |
---|
663 | ! |
---|
664 | !-- Total and perturbation energy for the total domain (being collected in the last column |
---|
665 | !-- of sums_l). Summation of these quantities is seperated from the previous loop in order |
---|
666 | !-- to allow vectorization of that loop. |
---|
667 | !$ACC ATOMIC |
---|
668 | sums_l(nzb+4,pr_palm,tn) = sums_l(nzb+4,pr_palm,tn) + sums_l_etot |
---|
669 | ! |
---|
670 | !-- 2D-arrays (being collected in the last column of sums_l) |
---|
671 | IF ( surf_def_h(0)%end_index(j,i) >= surf_def_h(0)%start_index(j,i) ) THEN |
---|
672 | m = surf_def_h(0)%start_index(j,i) |
---|
673 | !$ACC ATOMIC |
---|
674 | sums_l(nzb,pr_palm,tn) = sums_l(nzb,pr_palm,tn) + & |
---|
675 | surf_def_h(0)%us(m) * rmask(j,i,sr) |
---|
676 | !$ACC ATOMIC |
---|
677 | sums_l(nzb+1,pr_palm,tn) = sums_l(nzb+1,pr_palm,tn) + & |
---|
678 | surf_def_h(0)%usws(m) * rmask(j,i,sr) |
---|
679 | !$ACC ATOMIC |
---|
680 | sums_l(nzb+2,pr_palm,tn) = sums_l(nzb+2,pr_palm,tn) + & |
---|
681 | surf_def_h(0)%vsws(m) * rmask(j,i,sr) |
---|
682 | !$ACC ATOMIC |
---|
683 | sums_l(nzb+3,pr_palm,tn) = sums_l(nzb+3,pr_palm,tn) + & |
---|
684 | surf_def_h(0)%ts(m) * rmask(j,i,sr) |
---|
685 | #ifndef _OPENACC |
---|
686 | IF ( humidity ) THEN |
---|
687 | sums_l(nzb+12,pr_palm,tn) = sums_l(nzb+12,pr_palm,tn) + & |
---|
688 | surf_def_h(0)%qs(m) * rmask(j,i,sr) |
---|
689 | ENDIF |
---|
690 | IF ( passive_scalar ) THEN |
---|
691 | sums_l(nzb+13,pr_palm,tn) = sums_l(nzb+13,pr_palm,tn) + & |
---|
692 | surf_def_h(0)%ss(m) * rmask(j,i,sr) |
---|
693 | ENDIF |
---|
694 | #endif |
---|
695 | ! |
---|
696 | !-- Summation of surface temperature. |
---|
697 | !$ACC ATOMIC |
---|
698 | sums_l(nzb+14,pr_palm,tn) = sums_l(nzb+14,pr_palm,tn) + & |
---|
699 | surf_def_h(0)%pt_surface(m) * rmask(j,i,sr) |
---|
700 | ENDIF |
---|
701 | IF ( surf_lsm_h(0)%end_index(j,i) >= surf_lsm_h(0)%start_index(j,i) ) THEN |
---|
702 | m = surf_lsm_h(0)%start_index(j,i) |
---|
703 | !$ACC ATOMIC |
---|
704 | sums_l(nzb,pr_palm,tn) = sums_l(nzb,pr_palm,tn) + & |
---|
705 | surf_lsm_h(0)%us(m) * rmask(j,i,sr) |
---|
706 | !$ACC ATOMIC |
---|
707 | sums_l(nzb+1,pr_palm,tn) = sums_l(nzb+1,pr_palm,tn) + & |
---|
708 | surf_lsm_h(0)%usws(m) * rmask(j,i,sr) |
---|
709 | !$ACC ATOMIC |
---|
710 | sums_l(nzb+2,pr_palm,tn) = sums_l(nzb+2,pr_palm,tn) + & |
---|
711 | surf_lsm_h(0)%vsws(m) * rmask(j,i,sr) |
---|
712 | !$ACC ATOMIC |
---|
713 | sums_l(nzb+3,pr_palm,tn) = sums_l(nzb+3,pr_palm,tn) + & |
---|
714 | surf_lsm_h(0)%ts(m) * rmask(j,i,sr) |
---|
715 | #ifndef _OPENACC |
---|
716 | IF ( humidity ) THEN |
---|
717 | sums_l(nzb+12,pr_palm,tn) = sums_l(nzb+12,pr_palm,tn) + & |
---|
718 | surf_lsm_h(0)%qs(m) * rmask(j,i,sr) |
---|
719 | ENDIF |
---|
720 | IF ( passive_scalar ) THEN |
---|
721 | sums_l(nzb+13,pr_palm,tn) = sums_l(nzb+13,pr_palm,tn) + & |
---|
722 | surf_lsm_h(0)%ss(m) * rmask(j,i,sr) |
---|
723 | ENDIF |
---|
724 | #endif |
---|
725 | ! |
---|
726 | !-- Summation of surface temperature. |
---|
727 | !$ACC ATOMIC |
---|
728 | sums_l(nzb+14,pr_palm,tn) = sums_l(nzb+14,pr_palm,tn) + & |
---|
729 | surf_lsm_h(0)%pt_surface(m) * rmask(j,i,sr) |
---|
730 | ENDIF |
---|
731 | IF ( surf_usm_h(0)%end_index(j,i) >= surf_usm_h(0)%start_index(j,i) ) THEN |
---|
732 | m = surf_usm_h(0)%start_index(j,i) |
---|
733 | !$ACC ATOMIC |
---|
734 | sums_l(nzb,pr_palm,tn) = sums_l(nzb,pr_palm,tn) + & |
---|
735 | surf_usm_h(0)%us(m) * rmask(j,i,sr) |
---|
736 | !$ACC ATOMIC |
---|
737 | sums_l(nzb+1,pr_palm,tn) = sums_l(nzb+1,pr_palm,tn) + & |
---|
738 | surf_usm_h(0)%usws(m) * rmask(j,i,sr) |
---|
739 | !$ACC ATOMIC |
---|
740 | sums_l(nzb+2,pr_palm,tn) = sums_l(nzb+2,pr_palm,tn) + & |
---|
741 | surf_usm_h(0)%vsws(m) * rmask(j,i,sr) |
---|
742 | !$ACC ATOMIC |
---|
743 | sums_l(nzb+3,pr_palm,tn) = sums_l(nzb+3,pr_palm,tn) + & |
---|
744 | surf_usm_h(0)%ts(m) * rmask(j,i,sr) |
---|
745 | #ifndef _OPENACC |
---|
746 | IF ( humidity ) THEN |
---|
747 | sums_l(nzb+12,pr_palm,tn) = sums_l(nzb+12,pr_palm,tn) + & |
---|
748 | surf_usm_h(0)%qs(m) * rmask(j,i,sr) |
---|
749 | ENDIF |
---|
750 | IF ( passive_scalar ) THEN |
---|
751 | sums_l(nzb+13,pr_palm,tn) = sums_l(nzb+13,pr_palm,tn) + & |
---|
752 | surf_usm_h(0)%ss(m) * rmask(j,i,sr) |
---|
753 | ENDIF |
---|
754 | #endif |
---|
755 | ! |
---|
756 | !-- Summation of surface temperature. |
---|
757 | !$ACC ATOMIC |
---|
758 | sums_l(nzb+14,pr_palm,tn) = sums_l(nzb+14,pr_palm,tn) + & |
---|
759 | surf_usm_h(0)%pt_surface(m) * rmask(j,i,sr) |
---|
760 | ENDIF |
---|
761 | ENDDO !j-loop |
---|
762 | ENDDO !i-loop |
---|
763 | !$ACC UPDATE & |
---|
764 | !$ACC HOST(sums_l(:,3,tn), sums_l(:,8,tn), sums_l(:,9,tn)) & |
---|
765 | !$ACC HOST(sums_l(:,10,tn), sums_l(:,40,tn), sums_l(:,33,tn)) & |
---|
766 | !$ACC HOST(sums_l(:,38,tn), sums_l(:,121,tn)) & |
---|
767 | !$ACC HOST(sums_l(nzb:nzb+4,pr_palm,tn), sums_l(nzb+14:nzb+14,pr_palm,tn)) |
---|
768 | |
---|
769 | ! |
---|
770 | !-- Computation of statistics when ws-scheme is not used. Else these |
---|
771 | !-- quantities are evaluated in the advection routines. |
---|
772 | IF ( .NOT. ws_scheme_mom .OR. sr /= 0 .OR. simulated_time == 0.0_wp ) THEN |
---|
773 | !$OMP DO |
---|
774 | DO i = nxl, nxr |
---|
775 | DO j = nys, nyn |
---|
776 | DO k = nzb, nzt+1 |
---|
777 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
778 | |
---|
779 | u2 = u(k,j,i)**2 |
---|
780 | v2 = v(k,j,i)**2 |
---|
781 | w2 = w(k,j,i)**2 |
---|
782 | ust2 = ( u(k,j,i) - hom(k,1,1,sr) )**2 |
---|
783 | vst2 = ( v(k,j,i) - hom(k,1,2,sr) )**2 |
---|
784 | |
---|
785 | sums_l(k,30,tn) = sums_l(k,30,tn) + ust2 * rmask(j,i,sr) * flag |
---|
786 | sums_l(k,31,tn) = sums_l(k,31,tn) + vst2 * rmask(j,i,sr) * flag |
---|
787 | sums_l(k,32,tn) = sums_l(k,32,tn) + w2 * rmask(j,i,sr) * flag |
---|
788 | ! |
---|
789 | !-- Perturbation energy |
---|
790 | sums_l(k,34,tn) = sums_l(k,34,tn) + & |
---|
791 | 0.5_wp * ( ust2 + vst2 + w2 ) * rmask(j,i,sr) * flag |
---|
792 | ENDDO |
---|
793 | ENDDO |
---|
794 | ENDDO |
---|
795 | ENDIF |
---|
796 | ! |
---|
797 | !-- Computaion of domain-averaged perturbation energy. Please note, to prevent that perturbation |
---|
798 | !-- energy is larger (even if only slightly) than the total kinetic energy, calculation is based |
---|
799 | !-- on deviations from the horizontal mean, instead of spatial descretization of the advection |
---|
800 | !-- term. |
---|
801 | !$OMP DO |
---|
802 | !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k, flag, w2, ust2, vst2) & |
---|
803 | !$ACC PRESENT(wall_flags_total_0, u, v, w, rmask, hom(:,1,1:2,sr)) & |
---|
804 | !$ACC PRESENT(sums_l) |
---|
805 | DO i = nxl, nxr |
---|
806 | DO j = nys, nyn |
---|
807 | DO k = nzb, nzt+1 |
---|
808 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
809 | |
---|
810 | w2 = w(k,j,i)**2 |
---|
811 | ust2 = ( u(k,j,i) - hom(k,1,1,sr) )**2 |
---|
812 | vst2 = ( v(k,j,i) - hom(k,1,2,sr) )**2 |
---|
813 | w2 = w(k,j,i)**2 |
---|
814 | |
---|
815 | !$ACC ATOMIC |
---|
816 | sums_l(nzb+5,pr_palm,tn) = sums_l(nzb+5,pr_palm,tn) & |
---|
817 | + 0.5_wp * ( ust2 + vst2 + w2 ) * rmask(j,i,sr) * flag |
---|
818 | |
---|
819 | ENDDO |
---|
820 | ENDDO |
---|
821 | ENDDO |
---|
822 | !$ACC UPDATE HOST(sums_l(nzb+5:nzb+5,pr_palm,tn)) |
---|
823 | |
---|
824 | ! |
---|
825 | !-- Horizontally averaged profiles of the vertical fluxes |
---|
826 | !$OMP DO |
---|
827 | !$ACC PARALLEL LOOP COLLAPSE(2) PRIVATE(i, j, k, l, m) & |
---|
828 | !$ACC PRIVATE(ki, flag, ust, vst, pts) & |
---|
829 | !$ACC PRESENT(kh, km, u, v, w, pt) & |
---|
830 | !$ACC PRESENT(wall_flags_total_0, rmask, ddzu, rho_air_zw, hom(:,1,1:4,sr)) & |
---|
831 | !$ACC PRESENT(heatflux_output_conversion, momentumflux_output_conversion) & |
---|
832 | !$ACC PRESENT(surf_def_h(0:2), surf_lsm_h(0:1), surf_usm_h(0:1)) & |
---|
833 | !$ACC PRESENT(sums_l) |
---|
834 | DO i = nxl, nxr |
---|
835 | DO j = nys, nyn |
---|
836 | ! |
---|
837 | !-- Subgridscale fluxes (without Prandtl layer from k=nzb, oterwise from k=nzb+1) |
---|
838 | !-- NOTE: for simplicity, nzb_diff_s_inner is used below, although strictly speaking the |
---|
839 | !-- ---- following k-loop would have to be split up according to the staggered grid. |
---|
840 | !-- However, this implies no error since staggered velocity components are zero at |
---|
841 | !-- the walls and inside buildings. |
---|
842 | !-- Flag 23 is used to mask surface fluxes as well as model-top fluxes, which are added |
---|
843 | !-- further below. |
---|
844 | DO k = nzb, nzt |
---|
845 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 23 ) ) * & |
---|
846 | MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 9 ) ) |
---|
847 | ! |
---|
848 | !-- Momentum flux w"u" |
---|
849 | !$ACC ATOMIC |
---|
850 | sums_l(k,12,tn) = sums_l(k,12,tn) - 0.25_wp * ( & |
---|
851 | km(k,j,i) + km(k+1,j,i) + km(k,j,i-1) + km(k+1,j,i-1) & |
---|
852 | ) * ( & |
---|
853 | ( u(k+1,j,i) - u(k,j,i) ) * ddzu(k+1) & |
---|
854 | + ( w(k,j,i) - w(k,j,i-1) ) * ddx & |
---|
855 | ) * rmask(j,i,sr) & |
---|
856 | * rho_air_zw(k) & |
---|
857 | * momentumflux_output_conversion(k) & |
---|
858 | * flag |
---|
859 | ! |
---|
860 | !-- Momentum flux w"v" |
---|
861 | !$ACC ATOMIC |
---|
862 | sums_l(k,14,tn) = sums_l(k,14,tn) - 0.25_wp * ( & |
---|
863 | km(k,j,i) + km(k+1,j,i) + km(k,j-1,i) + km(k+1,j-1,i) & |
---|
864 | ) * ( & |
---|
865 | ( v(k+1,j,i) - v(k,j,i) ) * ddzu(k+1) & |
---|
866 | + ( w(k,j,i) - w(k,j-1,i) ) * ddy & |
---|
867 | ) * rmask(j,i,sr) & |
---|
868 | * rho_air_zw(k) & |
---|
869 | * momentumflux_output_conversion(k) & |
---|
870 | * flag |
---|
871 | ! |
---|
872 | !-- Heat flux w"pt" |
---|
873 | !$ACC ATOMIC |
---|
874 | sums_l(k,16,tn) = sums_l(k,16,tn) & |
---|
875 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
876 | * ( pt(k+1,j,i) - pt(k,j,i) ) & |
---|
877 | * rho_air_zw(k) & |
---|
878 | * heatflux_output_conversion(k) & |
---|
879 | * ddzu(k+1) * rmask(j,i,sr) & |
---|
880 | * flag |
---|
881 | |
---|
882 | ! |
---|
883 | !-- Salinity flux w"sa" |
---|
884 | #ifndef _OPENACC |
---|
885 | IF ( ocean_mode ) THEN |
---|
886 | sums_l(k,65,tn) = sums_l(k,65,tn) & |
---|
887 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
888 | * ( sa(k+1,j,i) - sa(k,j,i) ) & |
---|
889 | * ddzu(k+1) * rmask(j,i,sr) & |
---|
890 | * flag |
---|
891 | ENDIF |
---|
892 | |
---|
893 | ! |
---|
894 | !-- Buoyancy flux, water flux (humidity flux) w"q" |
---|
895 | IF ( humidity ) THEN |
---|
896 | sums_l(k,45,tn) = sums_l(k,45,tn) & |
---|
897 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
898 | * ( vpt(k+1,j,i) - vpt(k,j,i) ) & |
---|
899 | * rho_air_zw(k) & |
---|
900 | * heatflux_output_conversion(k) & |
---|
901 | * ddzu(k+1) * rmask(j,i,sr) * flag |
---|
902 | sums_l(k,48,tn) = sums_l(k,48,tn) & |
---|
903 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
904 | * ( q(k+1,j,i) - q(k,j,i) ) & |
---|
905 | * rho_air_zw(k) & |
---|
906 | * waterflux_output_conversion(k) & |
---|
907 | * ddzu(k+1) * rmask(j,i,sr) * flag |
---|
908 | |
---|
909 | IF ( bulk_cloud_model ) THEN |
---|
910 | sums_l(k,51,tn) = sums_l(k,51,tn) & |
---|
911 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
912 | * ( ( q(k+1,j,i) - ql(k+1,j,i) ) & |
---|
913 | - ( q(k,j,i) - ql(k,j,i) ) ) & |
---|
914 | * rho_air_zw(k) & |
---|
915 | * waterflux_output_conversion(k) & |
---|
916 | * ddzu(k+1) * rmask(j,i,sr) * flag |
---|
917 | ENDIF |
---|
918 | ENDIF |
---|
919 | |
---|
920 | ! |
---|
921 | !-- Passive scalar flux |
---|
922 | IF ( passive_scalar ) THEN |
---|
923 | sums_l(k,117,tn) = sums_l(k,117,tn) & |
---|
924 | - 0.5_wp * ( kh(k,j,i) + kh(k+1,j,i) ) & |
---|
925 | * ( s(k+1,j,i) - s(k,j,i) ) & |
---|
926 | * ddzu(k+1) * rmask(j,i,sr) & |
---|
927 | * flag |
---|
928 | ENDIF |
---|
929 | #endif |
---|
930 | |
---|
931 | ENDDO |
---|
932 | |
---|
933 | ! |
---|
934 | !-- Subgridscale fluxes in the Prandtl layer |
---|
935 | IF ( use_surface_fluxes ) THEN |
---|
936 | DO l = 0, 1 |
---|
937 | ! The original code using MERGE doesn't work with the PGI |
---|
938 | ! compiler when running on the GPU. |
---|
939 | ! This is submitted as a compiler Bug in PGI ticket TPR#26718 |
---|
940 | ! ki = MERGE( -1, 0, l == 0 ) |
---|
941 | ki = -1 + l |
---|
942 | IF ( surf_def_h(l)%ns >= 1 ) THEN |
---|
943 | DO m = surf_def_h(l)%start_index(j,i), & |
---|
944 | surf_def_h(l)%end_index(j,i) |
---|
945 | k = surf_def_h(l)%k(m) |
---|
946 | |
---|
947 | !$ACC ATOMIC |
---|
948 | sums_l(k+ki,12,tn) = sums_l(k+ki,12,tn) + & |
---|
949 | momentumflux_output_conversion(k+ki) * & |
---|
950 | surf_def_h(l)%usws(m) * rmask(j,i,sr) ! w"u" |
---|
951 | !$ACC ATOMIC |
---|
952 | sums_l(k+ki,14,tn) = sums_l(k+ki,14,tn) + & |
---|
953 | momentumflux_output_conversion(k+ki) * & |
---|
954 | surf_def_h(l)%vsws(m) * rmask(j,i,sr) ! w"v" |
---|
955 | !$ACC ATOMIC |
---|
956 | sums_l(k+ki,16,tn) = sums_l(k+ki,16,tn) + & |
---|
957 | heatflux_output_conversion(k+ki) * & |
---|
958 | surf_def_h(l)%shf(m) * rmask(j,i,sr) ! w"pt" |
---|
959 | #if 0 |
---|
960 | sums_l(k+ki,58,tn) = sums_l(k+ki,58,tn) + & |
---|
961 | 0.0_wp * rmask(j,i,sr) ! u"pt" |
---|
962 | sums_l(k+ki,61,tn) = sums_l(k+ki,61,tn) + & |
---|
963 | 0.0_wp * rmask(j,i,sr) ! v"pt" |
---|
964 | #endif |
---|
965 | #ifndef _OPENACC |
---|
966 | IF ( ocean_mode ) THEN |
---|
967 | sums_l(k+ki,65,tn) = sums_l(k+ki,65,tn) + & |
---|
968 | surf_def_h(l)%sasws(m) * rmask(j,i,sr) ! w"sa" |
---|
969 | ENDIF |
---|
970 | IF ( humidity ) THEN |
---|
971 | sums_l(k+ki,48,tn) = sums_l(k+ki,48,tn) + & |
---|
972 | waterflux_output_conversion(k+ki) * & |
---|
973 | surf_def_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
974 | sums_l(k+ki,45,tn) = sums_l(k+ki,45,tn) + ( & |
---|
975 | ( 1.0_wp + 0.61_wp * q(k+ki,j,i) ) * & |
---|
976 | surf_def_h(l)%shf(m) + 0.61_wp * pt(k+ki,j,i) * & |
---|
977 | surf_def_h(l)%qsws(m) ) & |
---|
978 | * heatflux_output_conversion(k+ki) |
---|
979 | IF ( cloud_droplets ) THEN |
---|
980 | sums_l(k+ki,45,tn) = sums_l(k+ki,45,tn) + ( & |
---|
981 | ( 1.0_wp + 0.61_wp * q(k+ki,j,i) - & |
---|
982 | ql(k+ki,j,i) ) * surf_def_h(l)%shf(m) + & |
---|
983 | 0.61_wp * pt(k+ki,j,i) & |
---|
984 | * surf_def_h(l)%qsws(m) ) & |
---|
985 | * heatflux_output_conversion(k+ki) |
---|
986 | ENDIF |
---|
987 | IF ( bulk_cloud_model ) THEN |
---|
988 | ! |
---|
989 | !-- Formula does not work if ql(k+ki) /= 0.0 |
---|
990 | sums_l(k+ki,51,tn) = sums_l(k+ki,51,tn) + & |
---|
991 | waterflux_output_conversion(k+ki) * & |
---|
992 | surf_def_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
993 | ENDIF |
---|
994 | ENDIF |
---|
995 | IF ( passive_scalar ) THEN |
---|
996 | sums_l(k+ki,117,tn) = sums_l(k+ki,117,tn) + & |
---|
997 | surf_def_h(l)%ssws(m) * rmask(j,i,sr) ! w"s" |
---|
998 | ENDIF |
---|
999 | #endif |
---|
1000 | |
---|
1001 | ENDDO |
---|
1002 | |
---|
1003 | ENDIF |
---|
1004 | IF ( surf_lsm_h(l)%end_index(j,i) >= surf_lsm_h(l)%start_index(j,i) ) THEN |
---|
1005 | m = surf_lsm_h(l)%start_index(j,i) |
---|
1006 | !$ACC ATOMIC |
---|
1007 | sums_l(nzb,12,tn) = sums_l(nzb,12,tn) + & |
---|
1008 | momentumflux_output_conversion(nzb) * & |
---|
1009 | surf_lsm_h(l)%usws(m) * rmask(j,i,sr) ! w"u" |
---|
1010 | !$ACC ATOMIC |
---|
1011 | sums_l(nzb,14,tn) = sums_l(nzb,14,tn) + & |
---|
1012 | momentumflux_output_conversion(nzb) * & |
---|
1013 | surf_lsm_h(l)%vsws(m) * rmask(j,i,sr) ! w"v" |
---|
1014 | !$ACC ATOMIC |
---|
1015 | sums_l(nzb,16,tn) = sums_l(nzb,16,tn) + & |
---|
1016 | heatflux_output_conversion(nzb) * & |
---|
1017 | surf_lsm_h(l)%shf(m) * rmask(j,i,sr) ! w"pt" |
---|
1018 | #if 0 |
---|
1019 | sums_l(nzb,58,tn) = sums_l(nzb,58,tn) + & |
---|
1020 | 0.0_wp * rmask(j,i,sr) ! u"pt" |
---|
1021 | sums_l(nzb,61,tn) = sums_l(nzb,61,tn) + & |
---|
1022 | 0.0_wp * rmask(j,i,sr) ! v"pt" |
---|
1023 | #endif |
---|
1024 | #ifndef _OPENACC |
---|
1025 | IF ( ocean_mode ) THEN |
---|
1026 | sums_l(nzb,65,tn) = sums_l(nzb,65,tn) + & |
---|
1027 | surf_lsm_h(l)%sasws(m) * rmask(j,i,sr) ! w"sa" |
---|
1028 | ENDIF |
---|
1029 | IF ( humidity ) THEN |
---|
1030 | sums_l(nzb,48,tn) = sums_l(nzb,48,tn) + & |
---|
1031 | waterflux_output_conversion(nzb) * & |
---|
1032 | surf_lsm_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
1033 | sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & |
---|
1034 | ( 1.0_wp + 0.61_wp * q(nzb,j,i) ) * surf_lsm_h(l)%shf(m) + & |
---|
1035 | 0.61_wp * pt(nzb,j,i) * surf_lsm_h(l)%qsws(m) ) & |
---|
1036 | * heatflux_output_conversion(nzb) |
---|
1037 | IF ( cloud_droplets ) THEN |
---|
1038 | sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & |
---|
1039 | ( 1.0_wp + 0.61_wp * q(nzb,j,i) - & |
---|
1040 | ql(nzb,j,i) ) * surf_lsm_h(l)%shf(m) + & |
---|
1041 | 0.61_wp * pt(nzb,j,i) * surf_lsm_h(l)%qsws(m) ) & |
---|
1042 | * heatflux_output_conversion(nzb) |
---|
1043 | ENDIF |
---|
1044 | IF ( bulk_cloud_model ) THEN |
---|
1045 | ! |
---|
1046 | !-- Formula does not work if ql(nzb) /= 0.0 |
---|
1047 | sums_l(nzb,51,tn) = sums_l(nzb,51,tn) + & |
---|
1048 | waterflux_output_conversion(nzb) * & |
---|
1049 | surf_lsm_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
1050 | ENDIF |
---|
1051 | ENDIF |
---|
1052 | IF ( passive_scalar ) THEN |
---|
1053 | sums_l(nzb,117,tn) = sums_l(nzb,117,tn) + & |
---|
1054 | surf_lsm_h(l)%ssws(m) * rmask(j,i,sr) ! w"s" |
---|
1055 | ENDIF |
---|
1056 | #endif |
---|
1057 | |
---|
1058 | ENDIF |
---|
1059 | IF ( surf_usm_h(l)%end_index(j,i) >= surf_usm_h(l)%start_index(j,i) ) THEN |
---|
1060 | m = surf_usm_h(l)%start_index(j,i) |
---|
1061 | !$ACC ATOMIC |
---|
1062 | sums_l(nzb,12,tn) = sums_l(nzb,12,tn) + & |
---|
1063 | momentumflux_output_conversion(nzb) * & |
---|
1064 | surf_usm_h(l)%usws(m) * rmask(j,i,sr) ! w"u" |
---|
1065 | !$ACC ATOMIC |
---|
1066 | sums_l(nzb,14,tn) = sums_l(nzb,14,tn) + & |
---|
1067 | momentumflux_output_conversion(nzb) * & |
---|
1068 | surf_usm_h(l)%vsws(m) * rmask(j,i,sr) ! w"v" |
---|
1069 | !$ACC ATOMIC |
---|
1070 | sums_l(nzb,16,tn) = sums_l(nzb,16,tn) + & |
---|
1071 | heatflux_output_conversion(nzb) * & |
---|
1072 | surf_usm_h(l)%shf(m) * rmask(j,i,sr) ! w"pt" |
---|
1073 | #if 0 |
---|
1074 | sums_l(nzb,58,tn) = sums_l(nzb,58,tn) + 0.0_wp * rmask(j,i,sr) ! u"pt" |
---|
1075 | sums_l(nzb,61,tn) = sums_l(nzb,61,tn) + 0.0_wp * rmask(j,i,sr) ! v"pt" |
---|
1076 | #endif |
---|
1077 | #ifndef _OPENACC |
---|
1078 | IF ( ocean_mode ) THEN |
---|
1079 | sums_l(nzb,65,tn) = sums_l(nzb,65,tn) + & |
---|
1080 | surf_usm_h(l)%sasws(m) * rmask(j,i,sr) ! w"sa" |
---|
1081 | ENDIF |
---|
1082 | IF ( humidity ) THEN |
---|
1083 | sums_l(nzb,48,tn) = sums_l(nzb,48,tn) + & |
---|
1084 | waterflux_output_conversion(nzb) * & |
---|
1085 | surf_usm_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
1086 | sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & |
---|
1087 | ( 1.0_wp + 0.61_wp * q(nzb,j,i) ) * & |
---|
1088 | surf_usm_h(l)%shf(m) + 0.61_wp * pt(nzb,j,i) * & |
---|
1089 | surf_usm_h(l)%qsws(m) ) & |
---|
1090 | * heatflux_output_conversion(nzb) |
---|
1091 | IF ( cloud_droplets ) THEN |
---|
1092 | sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & |
---|
1093 | ( 1.0_wp + 0.61_wp * q(nzb,j,i) - & |
---|
1094 | ql(nzb,j,i) ) * surf_usm_h(l)%shf(m) + & |
---|
1095 | 0.61_wp * pt(nzb,j,i) * surf_usm_h(l)%qsws(m) ) & |
---|
1096 | * heatflux_output_conversion(nzb) |
---|
1097 | ENDIF |
---|
1098 | IF ( bulk_cloud_model ) THEN |
---|
1099 | ! |
---|
1100 | !-- Formula does not work if ql(nzb) /= 0.0 |
---|
1101 | sums_l(nzb,51,tn) = sums_l(nzb,51,tn) + & |
---|
1102 | waterflux_output_conversion(nzb) * & |
---|
1103 | surf_usm_h(l)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
1104 | ENDIF |
---|
1105 | ENDIF |
---|
1106 | IF ( passive_scalar ) THEN |
---|
1107 | sums_l(nzb,117,tn) = sums_l(nzb,117,tn) + & |
---|
1108 | surf_usm_h(l)%ssws(m) * rmask(j,i,sr) ! w"s" |
---|
1109 | ENDIF |
---|
1110 | #endif |
---|
1111 | |
---|
1112 | ENDIF |
---|
1113 | ENDDO |
---|
1114 | ENDIF |
---|
1115 | |
---|
1116 | #ifndef _OPENACC |
---|
1117 | IF ( .NOT. neutral ) THEN |
---|
1118 | IF ( surf_def_h(0)%end_index(j,i) >= surf_def_h(0)%start_index(j,i) ) THEN |
---|
1119 | m = surf_def_h(0)%start_index(j,i) |
---|
1120 | sums_l(nzb,112,tn) = sums_l(nzb,112,tn) + surf_def_h(0)%ol(m) * rmask(j,i,sr) ! L |
---|
1121 | ENDIF |
---|
1122 | IF ( surf_lsm_h(0)%end_index(j,i) >= surf_lsm_h(0)%start_index(j,i) ) THEN |
---|
1123 | m = surf_lsm_h(0)%start_index(j,i) |
---|
1124 | sums_l(nzb,112,tn) = sums_l(nzb,112,tn) + surf_lsm_h(0)%ol(m) * rmask(j,i,sr) ! L |
---|
1125 | ENDIF |
---|
1126 | IF ( surf_usm_h(0)%end_index(j,i) >= surf_usm_h(0)%start_index(j,i) ) THEN |
---|
1127 | m = surf_usm_h(0)%start_index(j,i) |
---|
1128 | sums_l(nzb,112,tn) = sums_l(nzb,112,tn) + surf_usm_h(0)%ol(m) * rmask(j,i,sr) ! L |
---|
1129 | ENDIF |
---|
1130 | ENDIF |
---|
1131 | |
---|
1132 | IF ( radiation ) THEN |
---|
1133 | IF ( surf_def_h(0)%end_index(j,i) >= surf_def_h(0)%start_index(j,i) ) THEN |
---|
1134 | m = surf_def_h(0)%start_index(j,i) |
---|
1135 | sums_l(nzb,99,tn) = sums_l(nzb,99,tn) + & |
---|
1136 | surf_def_h(0)%rad_net(m) * rmask(j,i,sr) |
---|
1137 | sums_l(nzb,100,tn) = sums_l(nzb,100,tn) + & |
---|
1138 | surf_def_h(0)%rad_lw_in(m) * rmask(j,i,sr) |
---|
1139 | sums_l(nzb,101,tn) = sums_l(nzb,101,tn) + & |
---|
1140 | surf_def_h(0)%rad_lw_out(m) * rmask(j,i,sr) |
---|
1141 | sums_l(nzb,102,tn) = sums_l(nzb,102,tn) + & |
---|
1142 | surf_def_h(0)%rad_sw_in(m) * rmask(j,i,sr) |
---|
1143 | sums_l(nzb,103,tn) = sums_l(nzb,103,tn) + & |
---|
1144 | surf_def_h(0)%rad_sw_out(m) * rmask(j,i,sr) |
---|
1145 | ENDIF |
---|
1146 | IF ( surf_lsm_h(0)%end_index(j,i) >= surf_lsm_h(0)%start_index(j,i) ) THEN |
---|
1147 | m = surf_lsm_h(0)%start_index(j,i) |
---|
1148 | sums_l(nzb,99,tn) = sums_l(nzb,99,tn) + & |
---|
1149 | surf_lsm_h(0)%rad_net(m) * rmask(j,i,sr) |
---|
1150 | sums_l(nzb,100,tn) = sums_l(nzb,100,tn) + & |
---|
1151 | surf_lsm_h(0)%rad_lw_in(m) * rmask(j,i,sr) |
---|
1152 | sums_l(nzb,101,tn) = sums_l(nzb,101,tn) + & |
---|
1153 | surf_lsm_h(0)%rad_lw_out(m) * rmask(j,i,sr) |
---|
1154 | sums_l(nzb,102,tn) = sums_l(nzb,102,tn) + & |
---|
1155 | surf_lsm_h(0)%rad_sw_in(m) * rmask(j,i,sr) |
---|
1156 | sums_l(nzb,103,tn) = sums_l(nzb,103,tn) + & |
---|
1157 | surf_lsm_h(0)%rad_sw_out(m) * rmask(j,i,sr) |
---|
1158 | ENDIF |
---|
1159 | IF ( surf_usm_h(0)%end_index(j,i) >= surf_usm_h(0)%start_index(j,i) ) THEN |
---|
1160 | m = surf_usm_h(0)%start_index(j,i) |
---|
1161 | sums_l(nzb,99,tn) = sums_l(nzb,99,tn) + & |
---|
1162 | surf_usm_h(0)%rad_net(m) * rmask(j,i,sr) |
---|
1163 | sums_l(nzb,100,tn) = sums_l(nzb,100,tn) + & |
---|
1164 | surf_usm_h(0)%rad_lw_in(m) * rmask(j,i,sr) |
---|
1165 | sums_l(nzb,101,tn) = sums_l(nzb,101,tn) + & |
---|
1166 | surf_usm_h(0)%rad_lw_out(m) * rmask(j,i,sr) |
---|
1167 | sums_l(nzb,102,tn) = sums_l(nzb,102,tn) + & |
---|
1168 | surf_usm_h(0)%rad_sw_in(m) * rmask(j,i,sr) |
---|
1169 | sums_l(nzb,103,tn) = sums_l(nzb,103,tn) + & |
---|
1170 | surf_usm_h(0)%rad_sw_out(m) * rmask(j,i,sr) |
---|
1171 | ENDIF |
---|
1172 | |
---|
1173 | #if defined ( __rrtmg ) |
---|
1174 | IF ( radiation_scheme == 'rrtmg' ) THEN |
---|
1175 | |
---|
1176 | IF ( surf_def_h(0)%end_index(j,i) >= surf_def_h(0)%start_index(j,i) ) THEN |
---|
1177 | m = surf_def_h(0)%start_index(j,i) |
---|
1178 | sums_l(nzb,108,tn) = sums_l(nzb,108,tn) + & |
---|
1179 | surf_def_h(0)%rrtm_aldif(m,0) * rmask(j,i,sr) |
---|
1180 | sums_l(nzb,109,tn) = sums_l(nzb,109,tn) + & |
---|
1181 | surf_def_h(0)%rrtm_aldir(m,0) * rmask(j,i,sr) |
---|
1182 | sums_l(nzb,110,tn) = sums_l(nzb,110,tn) + & |
---|
1183 | surf_def_h(0)%rrtm_asdif(m,0) * rmask(j,i,sr) |
---|
1184 | sums_l(nzb,111,tn) = sums_l(nzb,111,tn) + & |
---|
1185 | surf_def_h(0)%rrtm_asdir(m,0) * rmask(j,i,sr) |
---|
1186 | ENDIF |
---|
1187 | IF ( surf_lsm_h(0)%end_index(j,i) >= surf_lsm_h(0)%start_index(j,i) ) THEN |
---|
1188 | m = surf_lsm_h(0)%start_index(j,i) |
---|
1189 | sums_l(nzb,108,tn) = sums_l(nzb,108,tn) + & |
---|
1190 | SUM( surf_lsm_h(0)%frac(m,:) * & |
---|
1191 | surf_lsm_h(0)%rrtm_aldif(m,:) ) * rmask(j,i,sr) |
---|
1192 | sums_l(nzb,109,tn) = sums_l(nzb,109,tn) + & |
---|
1193 | SUM( surf_lsm_h(0)%frac(m,:) * & |
---|
1194 | surf_lsm_h(0)%rrtm_aldir(m,:) ) * rmask(j,i,sr) |
---|
1195 | sums_l(nzb,110,tn) = sums_l(nzb,110,tn) + & |
---|
1196 | SUM( surf_lsm_h(0)%frac(m,:) * & |
---|
1197 | surf_lsm_h(0)%rrtm_asdif(m,:) ) * rmask(j,i,sr) |
---|
1198 | sums_l(nzb,111,tn) = sums_l(nzb,111,tn) + & |
---|
1199 | SUM( surf_lsm_h(0)%frac(m,:) * & |
---|
1200 | surf_lsm_h(0)%rrtm_asdir(m,:) ) * rmask(j,i,sr) |
---|
1201 | ENDIF |
---|
1202 | IF ( surf_usm_h(0)%end_index(j,i) >= surf_usm_h(0)%start_index(j,i) ) THEN |
---|
1203 | m = surf_usm_h(0)%start_index(j,i) |
---|
1204 | sums_l(nzb,108,tn) = sums_l(nzb,108,tn) + & |
---|
1205 | SUM( surf_usm_h(0)%frac(m,:) * & |
---|
1206 | surf_usm_h(0)%rrtm_aldif(m,:) ) * rmask(j,i,sr) |
---|
1207 | sums_l(nzb,109,tn) = sums_l(nzb,109,tn) + & |
---|
1208 | SUM( surf_usm_h(0)%frac(m,:) * & |
---|
1209 | surf_usm_h(0)%rrtm_aldir(m,:) ) * rmask(j,i,sr) |
---|
1210 | sums_l(nzb,110,tn) = sums_l(nzb,110,tn) + & |
---|
1211 | SUM( surf_usm_h(0)%frac(m,:) * & |
---|
1212 | surf_usm_h(0)%rrtm_asdif(m,:) ) * rmask(j,i,sr) |
---|
1213 | sums_l(nzb,111,tn) = sums_l(nzb,111,tn) + & |
---|
1214 | SUM( surf_usm_h(0)%frac(m,:) * & |
---|
1215 | surf_usm_h(0)%rrtm_asdir(m,:) ) * rmask(j,i,sr) |
---|
1216 | ENDIF |
---|
1217 | |
---|
1218 | ENDIF |
---|
1219 | #endif |
---|
1220 | ENDIF |
---|
1221 | #endif |
---|
1222 | ! |
---|
1223 | !-- Subgridscale fluxes at the top surface |
---|
1224 | IF ( use_top_fluxes ) THEN |
---|
1225 | m = surf_def_h(2)%start_index(j,i) |
---|
1226 | !$ACC ATOMIC |
---|
1227 | sums_l(nzt,12,tn) = sums_l(nzt,12,tn) + & |
---|
1228 | momentumflux_output_conversion(nzt) * & |
---|
1229 | surf_def_h(2)%usws(m) * rmask(j,i,sr) ! w"u" |
---|
1230 | !$ACC ATOMIC |
---|
1231 | sums_l(nzt+1,12,tn) = sums_l(nzt+1,12,tn) + & |
---|
1232 | momentumflux_output_conversion(nzt+1) * & |
---|
1233 | surf_def_h(2)%usws(m) * rmask(j,i,sr) ! w"u" |
---|
1234 | !$ACC ATOMIC |
---|
1235 | sums_l(nzt,14,tn) = sums_l(nzt,14,tn) + & |
---|
1236 | momentumflux_output_conversion(nzt) * & |
---|
1237 | surf_def_h(2)%vsws(m) * rmask(j,i,sr) ! w"v" |
---|
1238 | !$ACC ATOMIC |
---|
1239 | sums_l(nzt+1,14,tn) = sums_l(nzt+1,14,tn) + & |
---|
1240 | momentumflux_output_conversion(nzt+1) * & |
---|
1241 | surf_def_h(2)%vsws(m) * rmask(j,i,sr) ! w"v" |
---|
1242 | !$ACC ATOMIC |
---|
1243 | sums_l(nzt,16,tn) = sums_l(nzt,16,tn) + & |
---|
1244 | heatflux_output_conversion(nzt) * & |
---|
1245 | surf_def_h(2)%shf(m) * rmask(j,i,sr) ! w"pt" |
---|
1246 | !$ACC ATOMIC |
---|
1247 | sums_l(nzt+1,16,tn) = sums_l(nzt+1,16,tn) + & |
---|
1248 | heatflux_output_conversion(nzt+1) * & |
---|
1249 | surf_def_h(2)%shf(m) * rmask(j,i,sr) ! w"pt" |
---|
1250 | #if 0 |
---|
1251 | sums_l(nzt:nzt+1,58,tn) = sums_l(nzt:nzt+1,58,tn) + & |
---|
1252 | 0.0_wp * rmask(j,i,sr) ! u"pt" |
---|
1253 | sums_l(nzt:nzt+1,61,tn) = sums_l(nzt:nzt+1,61,tn) + & |
---|
1254 | 0.0_wp * rmask(j,i,sr) ! v"pt" |
---|
1255 | #endif |
---|
1256 | #ifndef _OPENACC |
---|
1257 | IF ( ocean_mode ) THEN |
---|
1258 | sums_l(nzt,65,tn) = sums_l(nzt,65,tn) + & |
---|
1259 | surf_def_h(2)%sasws(m) * rmask(j,i,sr) ! w"sa" |
---|
1260 | ENDIF |
---|
1261 | IF ( humidity ) THEN |
---|
1262 | sums_l(nzt,48,tn) = sums_l(nzt,48,tn) + & |
---|
1263 | waterflux_output_conversion(nzt) * & |
---|
1264 | surf_def_h(2)%qsws(m) * rmask(j,i,sr) ! w"q" (w"qv") |
---|
1265 | sums_l(nzt,45,tn) = sums_l(nzt,45,tn) + ( & |
---|
1266 | ( 1.0_wp + 0.61_wp * q(nzt,j,i) ) * & |
---|
1267 | surf_def_h(2)%shf(m) + & |
---|
1268 | 0.61_wp * pt(nzt,j,i) * & |
---|
1269 | surf_def_h(2)%qsws(m) ) & |
---|
1270 | * heatflux_output_conversion(nzt) |
---|
1271 | IF ( cloud_droplets ) THEN |
---|
1272 | sums_l(nzt,45,tn) = sums_l(nzt,45,tn) + ( & |
---|
1273 | ( 1.0_wp + 0.61_wp * q(nzt,j,i) - & |
---|
1274 | ql(nzt,j,i) ) * & |
---|
1275 | surf_def_h(2)%shf(m) + & |
---|
1276 | 0.61_wp * pt(nzt,j,i) * & |
---|
1277 | surf_def_h(2)%qsws(m) ) & |
---|
1278 | * heatflux_output_conversion(nzt) |
---|
1279 | ENDIF |
---|
1280 | IF ( bulk_cloud_model ) THEN |
---|
1281 | ! |
---|
1282 | !-- Formula does not work if ql(nzb) /= 0.0 |
---|
1283 | sums_l(nzt,51,tn) = sums_l(nzt,51,tn) + & ! w"q" (w"qv") |
---|
1284 | waterflux_output_conversion(nzt) * & |
---|
1285 | surf_def_h(2)%qsws(m) * rmask(j,i,sr) |
---|
1286 | ENDIF |
---|
1287 | ENDIF |
---|
1288 | IF ( passive_scalar ) THEN |
---|
1289 | sums_l(nzt,117,tn) = sums_l(nzt,117,tn) + & |
---|
1290 | surf_def_h(2)%ssws(m) * rmask(j,i,sr) ! w"s" |
---|
1291 | ENDIF |
---|
1292 | #endif |
---|
1293 | ENDIF |
---|
1294 | |
---|
1295 | ! |
---|
1296 | !-- Resolved fluxes (can be computed for all horizontal points) |
---|
1297 | !-- NOTE: for simplicity, nzb_s_inner is used below, although strictly speaking the |
---|
1298 | !-- ---- following k-loop would have to be split up and rearranged according to the |
---|
1299 | !-- staggered grid. |
---|
1300 | DO k = nzb, nzt |
---|
1301 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 22 ) ) |
---|
1302 | ust = 0.5_wp * ( u(k,j,i) - hom(k,1,1,sr) + & |
---|
1303 | u(k+1,j,i) - hom(k+1,1,1,sr) ) |
---|
1304 | vst = 0.5_wp * ( v(k,j,i) - hom(k,1,2,sr) + & |
---|
1305 | v(k+1,j,i) - hom(k+1,1,2,sr) ) |
---|
1306 | pts = 0.5_wp * ( pt(k,j,i) - hom(k,1,4,sr) + & |
---|
1307 | pt(k+1,j,i) - hom(k+1,1,4,sr) ) |
---|
1308 | ! |
---|
1309 | !-- Higher moments |
---|
1310 | !$ACC ATOMIC |
---|
1311 | sums_l(k,35,tn) = sums_l(k,35,tn) + pts * w(k,j,i)**2 * rmask(j,i,sr) * flag |
---|
1312 | !$ACC ATOMIC |
---|
1313 | sums_l(k,36,tn) = sums_l(k,36,tn) + pts**2 * w(k,j,i) * rmask(j,i,sr) * flag |
---|
1314 | |
---|
1315 | ! |
---|
1316 | !-- Salinity flux and density (density does not belong to here, but so far there is no |
---|
1317 | !-- other suitable place to calculate) |
---|
1318 | #ifndef _OPENACC |
---|
1319 | IF ( ocean_mode ) THEN |
---|
1320 | IF( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN |
---|
1321 | pts = 0.5_wp * ( sa(k,j,i) - hom(k,1,23,sr) + & |
---|
1322 | sa(k+1,j,i) - hom(k+1,1,23,sr) ) |
---|
1323 | sums_l(k,66,tn) = sums_l(k,66,tn) + pts * w(k,j,i) * & |
---|
1324 | rmask(j,i,sr) * flag |
---|
1325 | ENDIF |
---|
1326 | sums_l(k,64,tn) = sums_l(k,64,tn) + rho_ocean(k,j,i) * rmask(j,i,sr) * flag |
---|
1327 | sums_l(k,71,tn) = sums_l(k,71,tn) + prho(k,j,i) * rmask(j,i,sr) * flag |
---|
1328 | ENDIF |
---|
1329 | |
---|
1330 | ! |
---|
1331 | !-- Buoyancy flux, water flux, humidity flux, liquid water content, rain drop |
---|
1332 | !-- concentration and rain water content |
---|
1333 | IF ( humidity ) THEN |
---|
1334 | IF ( bulk_cloud_model .OR. cloud_droplets ) THEN |
---|
1335 | pts = 0.5_wp * ( vpt(k,j,i) - hom(k,1,44,sr) + & |
---|
1336 | vpt(k+1,j,i) - hom(k+1,1,44,sr) ) |
---|
1337 | sums_l(k,46,tn) = sums_l(k,46,tn) + pts * w(k,j,i) * & |
---|
1338 | rho_air_zw(k) * & |
---|
1339 | heatflux_output_conversion(k) * & |
---|
1340 | rmask(j,i,sr) * flag |
---|
1341 | sums_l(k,54,tn) = sums_l(k,54,tn) + ql(k,j,i) * rmask(j,i,sr) * flag |
---|
1342 | |
---|
1343 | IF ( .NOT. cloud_droplets ) THEN |
---|
1344 | pts = 0.5_wp * & |
---|
1345 | ( ( q(k,j,i) - ql(k,j,i) ) - & |
---|
1346 | hom(k,1,42,sr) + & |
---|
1347 | ( q(k+1,j,i) - ql(k+1,j,i) ) - & |
---|
1348 | hom(k+1,1,42,sr) ) |
---|
1349 | sums_l(k,52,tn) = sums_l(k,52,tn) + pts * w(k,j,i) * & |
---|
1350 | rho_air_zw(k) * & |
---|
1351 | waterflux_output_conversion(k) * & |
---|
1352 | rmask(j,i,sr) * flag |
---|
1353 | sums_l(k,75,tn) = sums_l(k,75,tn) + qc(k,j,i) * rmask(j,i,sr) * flag |
---|
1354 | sums_l(k,76,tn) = sums_l(k,76,tn) + prr(k,j,i) * rmask(j,i,sr) * flag |
---|
1355 | IF ( microphysics_morrison ) THEN |
---|
1356 | sums_l(k,123,tn) = sums_l(k,123,tn) + nc(k,j,i) * rmask(j,i,sr) * flag |
---|
1357 | ENDIF |
---|
1358 | IF ( microphysics_ice_phase ) THEN |
---|
1359 | sums_l(k,124,tn) = sums_l(k,124,tn) + ni(k,j,i) * rmask(j,i,sr) * flag |
---|
1360 | sums_l(k,125,tn) = sums_l(k,125,tn) + qi(k,j,i) * rmask(j,i,sr) * flag |
---|
1361 | IF ( graupel .AND. snow ) THEN |
---|
1362 | sums_l(k,126,tn) = sums_l(k,126,tn) + ng(k,j,i) * rmask(j,i,sr) * flag |
---|
1363 | sums_l(k,127,tn) = sums_l(k,127,tn) + qg(k,j,i) * rmask(j,i,sr) * flag |
---|
1364 | sums_l(k,128,tn) = sums_l(k,128,tn) + ns(k,j,i) * rmask(j,i,sr) * flag |
---|
1365 | sums_l(k,129,tn) = sums_l(k,129,tn) + qs(k,j,i) * rmask(j,i,sr) * flag |
---|
1366 | ENDIF |
---|
1367 | ENDIF |
---|
1368 | |
---|
1369 | IF ( microphysics_seifert ) THEN |
---|
1370 | sums_l(k,73,tn) = sums_l(k,73,tn) + nr(k,j,i) * rmask(j,i,sr) * flag |
---|
1371 | sums_l(k,74,tn) = sums_l(k,74,tn) + qr(k,j,i) * rmask(j,i,sr) * flag |
---|
1372 | ENDIF |
---|
1373 | ENDIF |
---|
1374 | |
---|
1375 | ELSE |
---|
1376 | IF( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN |
---|
1377 | pts = 0.5_wp * ( vpt(k,j,i) - hom(k,1,44,sr) + & |
---|
1378 | vpt(k+1,j,i) - hom(k+1,1,44,sr) ) |
---|
1379 | sums_l(k,46,tn) = sums_l(k,46,tn) + pts * w(k,j,i) * & |
---|
1380 | rho_air_zw(k) * & |
---|
1381 | heatflux_output_conversion(k) * & |
---|
1382 | rmask(j,i,sr) * flag |
---|
1383 | ELSE IF ( ws_scheme_sca .AND. sr == 0 ) THEN |
---|
1384 | sums_l(k,46,tn) = ( ( 1.0_wp + 0.61_wp * & |
---|
1385 | hom(k,1,41,sr) ) * & |
---|
1386 | sums_l(k,17,tn) + & |
---|
1387 | 0.61_wp * hom(k,1,4,sr) * & |
---|
1388 | sums_l(k,49,tn) & |
---|
1389 | ) * heatflux_output_conversion(k) * flag |
---|
1390 | END IF |
---|
1391 | END IF |
---|
1392 | ENDIF |
---|
1393 | ! |
---|
1394 | !-- Passive scalar flux |
---|
1395 | IF ( passive_scalar .AND. ( .NOT. ws_scheme_sca .OR. sr /= 0 ) ) THEN |
---|
1396 | pts = 0.5_wp * ( s(k,j,i) - hom(k,1,115,sr) + & |
---|
1397 | s(k+1,j,i) - hom(k+1,1,115,sr) ) |
---|
1398 | sums_l(k,114,tn) = sums_l(k,114,tn) + pts * w(k,j,i) * rmask(j,i,sr) * flag |
---|
1399 | ENDIF |
---|
1400 | #endif |
---|
1401 | |
---|
1402 | ! |
---|
1403 | !-- Energy flux w*e* |
---|
1404 | !-- has to be adjusted |
---|
1405 | !$ACC ATOMIC |
---|
1406 | sums_l(k,37,tn) = sums_l(k,37,tn) + w(k,j,i) * 0.5_wp * & |
---|
1407 | ( ust**2 + vst**2 + w(k,j,i)**2 ) & |
---|
1408 | * rho_air_zw(k) & |
---|
1409 | * momentumflux_output_conversion(k) & |
---|
1410 | * rmask(j,i,sr) * flag |
---|
1411 | ENDDO |
---|
1412 | ENDDO |
---|
1413 | ENDDO |
---|
1414 | !$OMP END PARALLEL |
---|
1415 | |
---|
1416 | !$ACC UPDATE & |
---|
1417 | !$ACC HOST(sums_l(:,12,tn), sums_l(:,14,tn), sums_l(:,16,tn)) & |
---|
1418 | !$ACC HOST(sums_l(:,35,tn), sums_l(:,36,tn), sums_l(:,37,tn)) |
---|
1419 | ! |
---|
1420 | !-- Treat land-surface quantities according to new wall model structure. |
---|
1421 | IF ( land_surface ) THEN |
---|
1422 | tn = 0 |
---|
1423 | !$OMP PARALLEL PRIVATE( i, j, m, tn ) |
---|
1424 | !$ tn = omp_get_thread_num() |
---|
1425 | !$OMP DO |
---|
1426 | DO m = 1, surf_lsm_h(0)%ns |
---|
1427 | i = surf_lsm_h(0)%i(m) |
---|
1428 | j = surf_lsm_h(0)%j(m) |
---|
1429 | |
---|
1430 | sums_l(nzb,93,tn) = sums_l(nzb,93,tn) + surf_lsm_h(0)%ghf(m) * rmask(j,i,sr) |
---|
1431 | sums_l(nzb,94,tn) = sums_l(nzb,94,tn) + surf_lsm_h(0)%qsws_liq(m) * rmask(j,i,sr) |
---|
1432 | sums_l(nzb,95,tn) = sums_l(nzb,95,tn) + surf_lsm_h(0)%qsws_soil(m) * rmask(j,i,sr) |
---|
1433 | sums_l(nzb,96,tn) = sums_l(nzb,96,tn) + surf_lsm_h(0)%qsws_veg(m) * rmask(j,i,sr) |
---|
1434 | sums_l(nzb,97,tn) = sums_l(nzb,97,tn) + surf_lsm_h(0)%r_a(m) * rmask(j,i,sr) |
---|
1435 | sums_l(nzb,98,tn) = sums_l(nzb,98,tn) + surf_lsm_h(0)%r_s(m) * rmask(j,i,sr) |
---|
1436 | ENDDO |
---|
1437 | !$OMP END PARALLEL |
---|
1438 | |
---|
1439 | tn = 0 |
---|
1440 | !$OMP PARALLEL PRIVATE( i, j, k, m, tn ) |
---|
1441 | !$ tn = omp_get_thread_num() |
---|
1442 | !$OMP DO |
---|
1443 | DO m = 1, surf_lsm_h(0)%ns |
---|
1444 | |
---|
1445 | i = surf_lsm_h(0)%i(m) |
---|
1446 | j = surf_lsm_h(0)%j(m) |
---|
1447 | DO k = nzb_soil, nzt_soil |
---|
1448 | sums_l(k,89,tn) = sums_l(k,89,tn) + t_soil_h(0)%var_2d(k,m) * rmask(j,i,sr) |
---|
1449 | sums_l(k,91,tn) = sums_l(k,91,tn) + m_soil_h(0)%var_2d(k,m) * rmask(j,i,sr) |
---|
1450 | ENDDO |
---|
1451 | ENDDO |
---|
1452 | !$OMP END PARALLEL |
---|
1453 | ENDIF |
---|
1454 | ! |
---|
1455 | !-- For speed optimization fluxes which have been computed in part directly inside the WS |
---|
1456 | !-- advection routines are treated seperatly. |
---|
1457 | !-- Momentum fluxes first: |
---|
1458 | |
---|
1459 | tn = 0 |
---|
1460 | !$OMP PARALLEL PRIVATE( i, j, k, tn, flag ) |
---|
1461 | !$ tn = omp_get_thread_num() |
---|
1462 | IF ( .NOT. ws_scheme_mom .OR. sr /= 0 ) THEN |
---|
1463 | !$OMP DO |
---|
1464 | DO i = nxl, nxr |
---|
1465 | DO j = nys, nyn |
---|
1466 | DO k = nzb, nzt |
---|
1467 | ! |
---|
1468 | !-- Flag 23 is used to mask surface fluxes as well as model-top fluxes, which are |
---|
1469 | !-- added further below. |
---|
1470 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 23 ) ) * & |
---|
1471 | MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 9 ) ) |
---|
1472 | |
---|
1473 | ust = 0.5_wp * ( u(k,j,i) - hom(k,1,1,sr) + & |
---|
1474 | u(k+1,j,i) - hom(k+1,1,1,sr) ) |
---|
1475 | vst = 0.5_wp * ( v(k,j,i) - hom(k,1,2,sr) + & |
---|
1476 | v(k+1,j,i) - hom(k+1,1,2,sr) ) |
---|
1477 | ! |
---|
1478 | !-- Momentum flux w*u* |
---|
1479 | sums_l(k,13,tn) = sums_l(k,13,tn) + 0.5_wp * & |
---|
1480 | ( w(k,j,i-1) + w(k,j,i) ) & |
---|
1481 | * rho_air_zw(k) & |
---|
1482 | * momentumflux_output_conversion(k) & |
---|
1483 | * ust * rmask(j,i,sr) & |
---|
1484 | * flag |
---|
1485 | ! |
---|
1486 | !-- Momentum flux w*v* |
---|
1487 | sums_l(k,15,tn) = sums_l(k,15,tn) + 0.5_wp * ( w(k,j-1,i) + w(k,j,i) ) & |
---|
1488 | * rho_air_zw(k) & |
---|
1489 | * momentumflux_output_conversion(k) & |
---|
1490 | * vst * rmask(j,i,sr) & |
---|
1491 | * flag |
---|
1492 | ENDDO |
---|
1493 | ENDDO |
---|
1494 | ENDDO |
---|
1495 | |
---|
1496 | ENDIF |
---|
1497 | IF ( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN |
---|
1498 | !$OMP DO |
---|
1499 | DO i = nxl, nxr |
---|
1500 | DO j = nys, nyn |
---|
1501 | DO k = nzb, nzt |
---|
1502 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 23 ) ) * & |
---|
1503 | MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 9 ) ) |
---|
1504 | ! |
---|
1505 | !-- Vertical heat flux |
---|
1506 | sums_l(k,17,tn) = sums_l(k,17,tn) + 0.5_wp * & |
---|
1507 | ( pt(k,j,i) - hom(k,1,4,sr) + & |
---|
1508 | pt(k+1,j,i) - hom(k+1,1,4,sr) ) & |
---|
1509 | * rho_air_zw(k) & |
---|
1510 | * heatflux_output_conversion(k) & |
---|
1511 | * w(k,j,i) * rmask(j,i,sr) * flag |
---|
1512 | IF ( humidity ) THEN |
---|
1513 | pts = 0.5_wp * ( q(k,j,i) - hom(k,1,41,sr) + & |
---|
1514 | q(k+1,j,i) - hom(k+1,1,41,sr) ) |
---|
1515 | sums_l(k,49,tn) = sums_l(k,49,tn) + pts * w(k,j,i) * & |
---|
1516 | rho_air_zw(k) * & |
---|
1517 | waterflux_output_conversion(k) * & |
---|
1518 | rmask(j,i,sr) * flag |
---|
1519 | ENDIF |
---|
1520 | IF ( passive_scalar ) THEN |
---|
1521 | pts = 0.5_wp * ( s(k,j,i) - hom(k,1,115,sr) + & |
---|
1522 | s(k+1,j,i) - hom(k+1,1,115,sr) ) |
---|
1523 | sums_l(k,114,tn) = sums_l(k,114,tn) + pts * w(k,j,i) * rmask(j,i,sr) * flag |
---|
1524 | ENDIF |
---|
1525 | ENDDO |
---|
1526 | ENDDO |
---|
1527 | ENDDO |
---|
1528 | |
---|
1529 | ENDIF |
---|
1530 | |
---|
1531 | ! |
---|
1532 | !-- Density at top follows Neumann condition |
---|
1533 | IF ( ocean_mode ) THEN |
---|
1534 | sums_l(nzt+1,64,tn) = sums_l(nzt,64,tn) |
---|
1535 | sums_l(nzt+1,71,tn) = sums_l(nzt,71,tn) |
---|
1536 | ENDIF |
---|
1537 | |
---|
1538 | ! |
---|
1539 | !-- Divergence of vertical flux of resolved scale energy and pressure fluctuations as well as |
---|
1540 | !-- flux of pressure fluctuation itself (68). |
---|
1541 | !-- First calculate the products, then the divergence. |
---|
1542 | !-- Calculation is time consuming. Do it only, if profiles shall be plotted. |
---|
1543 | IF ( hom(nzb+1,2,55,0) /= 0.0_wp .OR. hom(nzb+1,2,68,0) /= 0.0_wp ) THEN |
---|
1544 | sums_ll = 0.0_wp ! local array |
---|
1545 | |
---|
1546 | !$OMP DO |
---|
1547 | DO i = nxl, nxr |
---|
1548 | DO j = nys, nyn |
---|
1549 | DO k = nzb+1, nzt |
---|
1550 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 0 ) ) |
---|
1551 | |
---|
1552 | sums_ll(k,1) = sums_ll(k,1) + 0.5_wp * w(k,j,i) * ( & |
---|
1553 | ( 0.25_wp * ( u(k,j,i)+u(k+1,j,i)+u(k,j,i+1)+u(k+1,j,i+1) ) & |
---|
1554 | - 0.5_wp * ( hom(k,1,1,sr) + hom(k+1,1,1,sr) ) )**2 & |
---|
1555 | + ( 0.25_wp * ( v(k,j,i)+v(k+1,j,i)+v(k,j+1,i)+v(k+1,j+1,i) ) & |
---|
1556 | - 0.5_wp * ( hom(k,1,2,sr) + hom(k+1,1,2,sr) ) )**2 & |
---|
1557 | + w(k,j,i)**2 ) * flag * rmask(j,i,sr) |
---|
1558 | |
---|
1559 | sums_ll(k,2) = sums_ll(k,2) + 0.5_wp * w(k,j,i) & |
---|
1560 | * ( ( p(k,j,i) + p(k+1,j,i) ) & |
---|
1561 | / momentumflux_output_conversion(k) ) & |
---|
1562 | * flag * rmask(j,i,sr) |
---|
1563 | |
---|
1564 | ENDDO |
---|
1565 | ENDDO |
---|
1566 | ENDDO |
---|
1567 | sums_ll(0,1) = 0.0_wp ! because w is zero at the bottom |
---|
1568 | sums_ll(nzt+1,1) = 0.0_wp |
---|
1569 | sums_ll(0,2) = 0.0_wp |
---|
1570 | sums_ll(nzt+1,2) = 0.0_wp |
---|
1571 | |
---|
1572 | DO k = nzb+1, nzt |
---|
1573 | sums_l(k,55,tn) = ( sums_ll(k,1) - sums_ll(k-1,1) ) * ddzw(k) |
---|
1574 | sums_l(k,56,tn) = ( sums_ll(k,2) - sums_ll(k-1,2) ) * ddzw(k) |
---|
1575 | sums_l(k,68,tn) = sums_ll(k,2) |
---|
1576 | ENDDO |
---|
1577 | sums_l(nzb,55,tn) = sums_l(nzb+1,55,tn) |
---|
1578 | sums_l(nzb,56,tn) = sums_l(nzb+1,56,tn) |
---|
1579 | sums_l(nzb,68,tn) = 0.0_wp ! because w* = 0 at nzb |
---|
1580 | |
---|
1581 | ENDIF |
---|
1582 | |
---|
1583 | ! |
---|
1584 | !-- Divergence of vertical flux of SGS TKE and the flux itself (69) |
---|
1585 | IF ( hom(nzb+1,2,57,0) /= 0.0_wp .OR. hom(nzb+1,2,69,0) /= 0.0_wp ) THEN |
---|
1586 | !$OMP DO |
---|
1587 | DO i = nxl, nxr |
---|
1588 | DO j = nys, nyn |
---|
1589 | DO k = nzb+1, nzt |
---|
1590 | |
---|
1591 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 0 ) ) |
---|
1592 | |
---|
1593 | sums_l(k,57,tn) = sums_l(k,57,tn) - 0.5_wp * ( & |
---|
1594 | (km(k,j,i)+km(k+1,j,i)) * (e(k+1,j,i)-e(k,j,i)) * ddzu(k+1) & |
---|
1595 | - (km(k-1,j,i)+km(k,j,i)) * (e(k,j,i)-e(k-1,j,i)) * ddzu(k) & |
---|
1596 | ) * ddzw(k) & |
---|
1597 | * flag * rmask(j,i,sr) |
---|
1598 | |
---|
1599 | sums_l(k,69,tn) = sums_l(k,69,tn) - 0.5_wp * ( & |
---|
1600 | ( km(k,j,i) + km(k+1,j,i) ) * & |
---|
1601 | ( e(k+1,j,i) - e(k,j,i) ) * ddzu(k+1) & |
---|
1602 | ) * flag * rmask(j,i,sr) |
---|
1603 | |
---|
1604 | ENDDO |
---|
1605 | ENDDO |
---|
1606 | ENDDO |
---|
1607 | sums_l(nzb,57,tn) = sums_l(nzb+1,57,tn) |
---|
1608 | sums_l(nzb,69,tn) = sums_l(nzb+1,69,tn) |
---|
1609 | |
---|
1610 | ENDIF |
---|
1611 | |
---|
1612 | ! |
---|
1613 | !-- Horizontal heat fluxes (subgrid, resolved, total). |
---|
1614 | !-- Do it only, if profiles shall be plotted. |
---|
1615 | IF ( hom(nzb+1,2,58,0) /= 0.0_wp ) THEN |
---|
1616 | |
---|
1617 | !$OMP DO |
---|
1618 | DO i = nxl, nxr |
---|
1619 | DO j = nys, nyn |
---|
1620 | DO k = nzb+1, nzt |
---|
1621 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 0 ) ) |
---|
1622 | ! |
---|
1623 | !-- Subgrid horizontal heat fluxes u"pt", v"pt" |
---|
1624 | sums_l(k,58,tn) = sums_l(k,58,tn) - 0.5_wp * & |
---|
1625 | ( kh(k,j,i) + kh(k,j,i-1) ) & |
---|
1626 | * ( pt(k,j,i-1) - pt(k,j,i) ) & |
---|
1627 | * rho_air_zw(k) & |
---|
1628 | * heatflux_output_conversion(k) & |
---|
1629 | * ddx * rmask(j,i,sr) * flag |
---|
1630 | sums_l(k,61,tn) = sums_l(k,61,tn) - 0.5_wp * & |
---|
1631 | ( kh(k,j,i) + kh(k,j-1,i) ) & |
---|
1632 | * ( pt(k,j-1,i) - pt(k,j,i) ) & |
---|
1633 | * rho_air_zw(k) & |
---|
1634 | * heatflux_output_conversion(k) & |
---|
1635 | * ddy * rmask(j,i,sr) * flag |
---|
1636 | ! |
---|
1637 | !-- Resolved horizontal heat fluxes u*pt*, v*pt* |
---|
1638 | sums_l(k,59,tn) = sums_l(k,59,tn) + ( u(k,j,i) - hom(k,1,1,sr) ) & |
---|
1639 | * 0.5_wp * ( pt(k,j,i-1) - hom(k,1,4,sr) + & |
---|
1640 | pt(k,j,i) - hom(k,1,4,sr) ) & |
---|
1641 | * heatflux_output_conversion(k) & |
---|
1642 | * flag |
---|
1643 | pts = 0.5_wp * ( pt(k,j-1,i) - hom(k,1,4,sr) + & |
---|
1644 | pt(k,j,i) - hom(k,1,4,sr) ) |
---|
1645 | sums_l(k,62,tn) = sums_l(k,62,tn) + ( v(k,j,i) - hom(k,1,2,sr) ) & |
---|
1646 | * 0.5_wp * ( pt(k,j-1,i) - hom(k,1,4,sr) + & |
---|
1647 | pt(k,j,i) - hom(k,1,4,sr) ) & |
---|
1648 | * heatflux_output_conversion(k) & |
---|
1649 | * flag |
---|
1650 | ENDDO |
---|
1651 | ENDDO |
---|
1652 | ENDDO |
---|
1653 | ! |
---|
1654 | !-- Fluxes at the surface must be zero (e.g. due to the Prandtl-layer) |
---|
1655 | sums_l(nzb,58,tn) = 0.0_wp |
---|
1656 | sums_l(nzb,59,tn) = 0.0_wp |
---|
1657 | sums_l(nzb,60,tn) = 0.0_wp |
---|
1658 | sums_l(nzb,61,tn) = 0.0_wp |
---|
1659 | sums_l(nzb,62,tn) = 0.0_wp |
---|
1660 | sums_l(nzb,63,tn) = 0.0_wp |
---|
1661 | |
---|
1662 | ENDIF |
---|
1663 | !$OMP END PARALLEL |
---|
1664 | |
---|
1665 | ! |
---|
1666 | !-- Collect current large scale advection and subsidence tendencies for |
---|
1667 | !-- data output |
---|
1668 | IF ( large_scale_forcing .AND. ( simulated_time > 0.0_wp ) ) THEN |
---|
1669 | ! |
---|
1670 | !-- Interpolation in time of LSF_DATA |
---|
1671 | nt = 1 |
---|
1672 | DO WHILE ( simulated_time - dt_3d > time_vert(nt) ) |
---|
1673 | nt = nt + 1 |
---|
1674 | ENDDO |
---|
1675 | IF ( simulated_time - dt_3d /= time_vert(nt) ) THEN |
---|
1676 | nt = nt - 1 |
---|
1677 | ENDIF |
---|
1678 | |
---|
1679 | fac = ( simulated_time - dt_3d - time_vert(nt) ) / ( time_vert(nt+1)-time_vert(nt) ) |
---|
1680 | |
---|
1681 | |
---|
1682 | DO k = nzb, nzt |
---|
1683 | sums_ls_l(k,0) = td_lsa_lpt(k,nt) + fac * ( td_lsa_lpt(k,nt+1) - td_lsa_lpt(k,nt) ) |
---|
1684 | sums_ls_l(k,1) = td_lsa_q(k,nt) + fac * ( td_lsa_q(k,nt+1) - td_lsa_q(k,nt) ) |
---|
1685 | ENDDO |
---|
1686 | |
---|
1687 | sums_ls_l(nzt+1,0) = sums_ls_l(nzt,0) |
---|
1688 | sums_ls_l(nzt+1,1) = sums_ls_l(nzt,1) |
---|
1689 | |
---|
1690 | IF ( large_scale_subsidence .AND. use_subsidence_tendencies ) THEN |
---|
1691 | |
---|
1692 | DO k = nzb, nzt |
---|
1693 | sums_ls_l(k,2) = td_sub_lpt(k,nt) + fac * ( td_sub_lpt(k,nt+1) - td_sub_lpt(k,nt) ) |
---|
1694 | sums_ls_l(k,3) = td_sub_q(k,nt) + fac * ( td_sub_q(k,nt+1) - td_sub_q(k,nt) ) |
---|
1695 | ENDDO |
---|
1696 | |
---|
1697 | sums_ls_l(nzt+1,2) = sums_ls_l(nzt,2) |
---|
1698 | sums_ls_l(nzt+1,3) = sums_ls_l(nzt,3) |
---|
1699 | |
---|
1700 | ENDIF |
---|
1701 | |
---|
1702 | ENDIF |
---|
1703 | |
---|
1704 | tn = 0 |
---|
1705 | !$OMP PARALLEL PRIVATE( i, j, k, tn ) |
---|
1706 | !$ tn = omp_get_thread_num() |
---|
1707 | IF ( radiation .AND. radiation_scheme == 'rrtmg' ) THEN |
---|
1708 | !$OMP DO |
---|
1709 | DO i = nxl, nxr |
---|
1710 | DO j = nys, nyn |
---|
1711 | DO k = nzb+1, nzt+1 |
---|
1712 | flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 0 ) ) |
---|
1713 | |
---|
1714 | sums_l(k,100,tn) = sums_l(k,100,tn) + rad_lw_in(k,j,i) * rmask(j,i,sr) * flag |
---|
1715 | sums_l(k,101,tn) = sums_l(k,101,tn) + rad_lw_out(k,j,i) * rmask(j,i,sr) * flag |
---|
1716 | sums_l(k,102,tn) = sums_l(k,102,tn) + rad_sw_in(k,j,i) * rmask(j,i,sr) * flag |
---|
1717 | sums_l(k,103,tn) = sums_l(k,103,tn) + rad_sw_out(k,j,i) * rmask(j,i,sr) * flag |
---|
1718 | sums_l(k,104,tn) = sums_l(k,104,tn) + rad_lw_cs_hr(k,j,i) * rmask(j,i,sr) * flag |
---|
1719 | sums_l(k,105,tn) = sums_l(k,105,tn) + rad_lw_hr(k,j,i) * rmask(j,i,sr) * flag |
---|
1720 | sums_l(k,106,tn) = sums_l(k,106,tn) + rad_sw_cs_hr(k,j,i) * rmask(j,i,sr) * flag |
---|
1721 | sums_l(k,107,tn) = sums_l(k,107,tn) + rad_sw_hr(k,j,i) * rmask(j,i,sr) * flag |
---|
1722 | ENDDO |
---|
1723 | ENDDO |
---|
1724 | ENDDO |
---|
1725 | ENDIF |
---|
1726 | |
---|
1727 | ! |
---|
1728 | !-- Calculate the profiles for all other modules |
---|
1729 | CALL module_interface_statistics( 'profiles', sr, tn, dots_max ) |
---|
1730 | !$OMP END PARALLEL |
---|
1731 | |
---|
1732 | ! |
---|
1733 | !-- Summation of thread sums |
---|
1734 | IF ( threads_per_task > 1 ) THEN |
---|
1735 | DO i = 1, threads_per_task-1 |
---|
1736 | sums_l(:,3,0) = sums_l(:,3,0) + sums_l(:,3,i) |
---|
1737 | sums_l(:,4:40,0) = sums_l(:,4:40,0) + sums_l(:,4:40,i) |
---|
1738 | sums_l(:,45:pr_palm,0) = sums_l(:,45:pr_palm,0) + & |
---|
1739 | sums_l(:,45:pr_palm,i) |
---|
1740 | IF ( max_pr_user > 0 ) THEN |
---|
1741 | sums_l(:,pr_palm+1:pr_palm+max_pr_user,0) = & |
---|
1742 | sums_l(:,pr_palm+1:pr_palm+max_pr_user,0) + & |
---|
1743 | sums_l(:,pr_palm+1:pr_palm+max_pr_user,i) |
---|
1744 | ENDIF |
---|
1745 | |
---|
1746 | IF ( air_chemistry ) THEN |
---|
1747 | IF ( max_pr_cs > 0 ) THEN |
---|
1748 | sums_l(:,pr_palm+max_pr_user+1:pr_palm + max_pr_user+ max_pr_cs,0) = & |
---|
1749 | sums_l(:,pr_palm+max_pr_user+1:pr_palm + max_pr_user+max_pr_cs,0) + & |
---|
1750 | sums_l(:,pr_palm+max_pr_user+1:pr_palm + max_pr_user+max_pr_cs,i) |
---|
1751 | |
---|
1752 | ENDIF |
---|
1753 | ENDIF |
---|
1754 | IF ( salsa ) THEN |
---|
1755 | IF ( max_pr_cs > 0 ) THEN |
---|
1756 | sums_l(:,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,0) = & |
---|
1757 | sums_l(:,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,0) + & |
---|
1758 | sums_l(:,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,i) |
---|
1759 | |
---|
1760 | ENDIF |
---|
1761 | ENDIF |
---|
1762 | ENDDO |
---|
1763 | ENDIF |
---|
1764 | |
---|
1765 | #if defined( __parallel ) |
---|
1766 | |
---|
1767 | ! |
---|
1768 | !-- Compute total sum from local sums |
---|
1769 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
1770 | CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), ngp_sums, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
1771 | IF ( large_scale_forcing ) THEN |
---|
1772 | CALL MPI_ALLREDUCE( sums_ls_l(nzb,2), sums(nzb,83), ngp_sums_ls, MPI_REAL, MPI_SUM, & |
---|
1773 | comm2d, ierr ) |
---|
1774 | ENDIF |
---|
1775 | |
---|
1776 | IF ( air_chemistry .AND. max_pr_cs > 0 ) THEN |
---|
1777 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
1778 | DO i = 1, max_pr_cs |
---|
1779 | CALL MPI_ALLREDUCE( sums_l(nzb,pr_palm+max_pr_user+i,0), & |
---|
1780 | sums(nzb,pr_palm+max_pr_user+i), & |
---|
1781 | nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
1782 | ENDDO |
---|
1783 | ENDIF |
---|
1784 | |
---|
1785 | IF ( salsa .AND. max_pr_salsa > 0 ) THEN |
---|
1786 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
1787 | DO i = 1, max_pr_salsa |
---|
1788 | CALL MPI_ALLREDUCE( sums_l(nzb,pr_palm+max_pr_user+max_pr_cs+i,0), & |
---|
1789 | sums(nzb,pr_palm+max_pr_user+max_pr_user+i), & |
---|
1790 | nzt+2-nzb, MPI_REAL, MPI_SUM, comm2d, ierr ) |
---|
1791 | ENDDO |
---|
1792 | ENDIF |
---|
1793 | |
---|
1794 | #else |
---|
1795 | sums = sums_l(:,:,0) |
---|
1796 | IF ( large_scale_forcing ) THEN |
---|
1797 | sums(:,81:88) = sums_ls_l |
---|
1798 | ENDIF |
---|
1799 | #endif |
---|
1800 | |
---|
1801 | ! |
---|
1802 | !-- Final values are obtained by division by the total number of grid points used for summation. |
---|
1803 | !-- After that store profiles. |
---|
1804 | !-- Check, if statistical regions do contain at least one grid point at the respective k-level, |
---|
1805 | !-- otherwise division by zero will lead to undefined values, which may cause e.g. problems with |
---|
1806 | !-- NetCDF output. |
---|
1807 | !-- Profiles: |
---|
1808 | DO k = nzb, nzt+1 |
---|
1809 | sums(k,3) = sums(k,3) / ngp_2dh(sr) |
---|
1810 | sums(k,12:22) = sums(k,12:22) / ngp_2dh(sr) |
---|
1811 | sums(k,30:32) = sums(k,30:32) / ngp_2dh(sr) |
---|
1812 | sums(k,35:39) = sums(k,35:39) / ngp_2dh(sr) |
---|
1813 | sums(k,45:53) = sums(k,45:53) / ngp_2dh(sr) |
---|
1814 | sums(k,55:63) = sums(k,55:63) / ngp_2dh(sr) |
---|
1815 | sums(k,81:88) = sums(k,81:88) / ngp_2dh(sr) |
---|
1816 | |
---|
1817 | IF ( land_surface .AND. surf_lsm_h(0)%ns_tot > 0 ) THEN |
---|
1818 | sums(k,89:98) = sums(k,89:98) / surf_lsm_h(0)%ns_tot |
---|
1819 | ENDIF |
---|
1820 | |
---|
1821 | sums(k,99:112) = sums(k,99:112) / ngp_2dh(sr) |
---|
1822 | sums(k,114) = sums(k,114) / ngp_2dh(sr) |
---|
1823 | sums(k,117) = sums(k,117) / ngp_2dh(sr) |
---|
1824 | IF ( ngp_2dh_s_inner(k,sr) /= 0 ) THEN |
---|
1825 | sums(k,8:11) = sums(k,8:11) / ngp_2dh_s_inner(k,sr) |
---|
1826 | sums(k,23:29) = sums(k,23:29) / ngp_2dh_s_inner(k,sr) |
---|
1827 | sums(k,33:34) = sums(k,33:34) / ngp_2dh_s_inner(k,sr) |
---|
1828 | sums(k,40) = sums(k,40) / ngp_2dh_s_inner(k,sr) |
---|
1829 | sums(k,54) = sums(k,54) / ngp_2dh_s_inner(k,sr) |
---|
1830 | sums(k,64) = sums(k,64) / ngp_2dh_s_inner(k,sr) |
---|
1831 | sums(k,70:80) = sums(k,70:80) / ngp_2dh_s_inner(k,sr) |
---|
1832 | sums(k,116) = sums(k,116) / ngp_2dh_s_inner(k,sr) |
---|
1833 | sums(k,118:pr_palm-2) = sums(k,118:pr_palm-2) / ngp_2dh_s_inner(k,sr) |
---|
1834 | sums(k,123:129) = sums(k,123:129) * ngp_2dh_s_inner(k,sr) / ngp_2dh(sr) |
---|
1835 | ENDIF |
---|
1836 | ENDDO |
---|
1837 | |
---|
1838 | !-- u* and so on |
---|
1839 | !-- As sums(nzb:nzb+3,pr_palm) are full 2D arrays (us, usws, vsws, ts) whose size is always |
---|
1840 | !-- ( nx + 1 ) * ( ny + 1 ), defined at the first grid layer above the topography, they are |
---|
1841 | !-- divided by ngp_2dh(sr) |
---|
1842 | sums(nzb:nzb+3,pr_palm) = sums(nzb:nzb+3,pr_palm) / ngp_2dh(sr) |
---|
1843 | sums(nzb+12,pr_palm) = sums(nzb+12,pr_palm) / ngp_2dh(sr) ! qs |
---|
1844 | sums(nzb+13,pr_palm) = sums(nzb+13,pr_palm) / ngp_2dh(sr) ! ss |
---|
1845 | sums(nzb+14,pr_palm) = sums(nzb+14,pr_palm) / ngp_2dh(sr) ! surface temperature |
---|
1846 | |
---|
1847 | !-- eges, e* |
---|
1848 | sums(nzb+4:nzb+5,pr_palm) = sums(nzb+4:nzb+5,pr_palm) / ngp_3d(sr) |
---|
1849 | !-- Old and new divergence |
---|
1850 | sums(nzb+9:nzb+10,pr_palm) = sums(nzb+9:nzb+10,pr_palm) / ngp_3d_inner(sr) |
---|
1851 | |
---|
1852 | !-- User-defined profiles |
---|
1853 | IF ( max_pr_user > 0 ) THEN |
---|
1854 | DO k = nzb, nzt+1 |
---|
1855 | sums(k,pr_palm+1:pr_palm+max_pr_user) = sums(k,pr_palm+1:pr_palm+max_pr_user) / & |
---|
1856 | ngp_2dh_s_inner(k,sr) |
---|
1857 | ENDDO |
---|
1858 | ENDIF |
---|
1859 | |
---|
1860 | IF ( air_chemistry ) THEN |
---|
1861 | IF ( max_pr_cs > 0 ) THEN |
---|
1862 | DO k = nzb, nzt+1 |
---|
1863 | sums(k, pr_palm+1:pr_palm+max_pr_user+max_pr_cs) = & |
---|
1864 | sums(k, pr_palm+1:pr_palm+max_pr_user+max_pr_cs) / & |
---|
1865 | ngp_2dh_s_inner(k,sr) |
---|
1866 | ENDDO |
---|
1867 | ENDIF |
---|
1868 | ENDIF |
---|
1869 | |
---|
1870 | IF ( salsa ) THEN |
---|
1871 | IF ( max_pr_salsa > 0 ) THEN |
---|
1872 | DO k = nzb, nzt+1 |
---|
1873 | sums(k,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa) = & |
---|
1874 | sums(k,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa) & |
---|
1875 | / ngp_2dh_s_inner(k,sr) |
---|
1876 | ENDDO |
---|
1877 | ENDIF |
---|
1878 | ENDIF |
---|
1879 | |
---|
1880 | ! |
---|
1881 | !-- Collect horizontal average in hom. |
---|
1882 | !-- Compute deduced averages (e.g. total heat flux) |
---|
1883 | hom(:,1,3,sr) = sums(:,3) ! w |
---|
1884 | hom(:,1,8,sr) = sums(:,8) ! e profiles 5-7 are initial profiles |
---|
1885 | hom(:,1,9,sr) = sums(:,9) ! km |
---|
1886 | hom(:,1,10,sr) = sums(:,10) ! kh |
---|
1887 | hom(:,1,11,sr) = sums(:,11) ! l |
---|
1888 | hom(:,1,12,sr) = sums(:,12) ! w"u" |
---|
1889 | hom(:,1,13,sr) = sums(:,13) ! w*u* |
---|
1890 | hom(:,1,14,sr) = sums(:,14) ! w"v" |
---|
1891 | hom(:,1,15,sr) = sums(:,15) ! w*v* |
---|
1892 | hom(:,1,16,sr) = sums(:,16) ! w"pt" |
---|
1893 | hom(:,1,17,sr) = sums(:,17) ! w*pt* |
---|
1894 | hom(:,1,18,sr) = sums(:,16) + sums(:,17) ! wpt |
---|
1895 | hom(:,1,19,sr) = sums(:,12) + sums(:,13) ! wu |
---|
1896 | hom(:,1,20,sr) = sums(:,14) + sums(:,15) ! wv |
---|
1897 | hom(:,1,21,sr) = sums(:,21) ! w*pt*BC |
---|
1898 | hom(:,1,22,sr) = sums(:,16) + sums(:,21) ! wptBC |
---|
1899 | ! profile 24 is initial profile (sa) |
---|
1900 | ! profiles 25-29 left empty for initial |
---|
1901 | ! profiles |
---|
1902 | hom(:,1,30,sr) = sums(:,30) ! u*2 |
---|
1903 | hom(:,1,31,sr) = sums(:,31) ! v*2 |
---|
1904 | hom(:,1,32,sr) = sums(:,32) ! w*2 |
---|
1905 | hom(:,1,33,sr) = sums(:,33) ! pt*2 |
---|
1906 | hom(:,1,34,sr) = sums(:,34) ! e* |
---|
1907 | hom(:,1,35,sr) = sums(:,35) ! w*2pt* |
---|
1908 | hom(:,1,36,sr) = sums(:,36) ! w*pt*2 |
---|
1909 | hom(:,1,37,sr) = sums(:,37) ! w*e* |
---|
1910 | hom(:,1,38,sr) = sums(:,38) ! w*3 |
---|
1911 | hom(:,1,39,sr) = sums(:,38) / ( ABS( sums(:,32) ) + 1E-20_wp )**1.5_wp ! Sw |
---|
1912 | hom(:,1,40,sr) = sums(:,40) ! p |
---|
1913 | hom(:,1,45,sr) = sums(:,45) ! w"vpt" |
---|
1914 | hom(:,1,46,sr) = sums(:,46) ! w*vpt* |
---|
1915 | hom(:,1,47,sr) = sums(:,45) + sums(:,46) ! wvpt |
---|
1916 | hom(:,1,48,sr) = sums(:,48) ! w"q" (w"qv") |
---|
1917 | hom(:,1,49,sr) = sums(:,49) ! w*q* (w*qv*) |
---|
1918 | hom(:,1,50,sr) = sums(:,48) + sums(:,49) ! wq (wqv) |
---|
1919 | hom(:,1,51,sr) = sums(:,51) ! w"qv" |
---|
1920 | hom(:,1,52,sr) = sums(:,52) ! w*qv* |
---|
1921 | hom(:,1,53,sr) = sums(:,52) + sums(:,51) ! wq (wqv) |
---|
1922 | hom(:,1,54,sr) = sums(:,54) ! ql |
---|
1923 | hom(:,1,55,sr) = sums(:,55) ! w*u*u*/dz |
---|
1924 | hom(:,1,56,sr) = sums(:,56) ! w*p*/dz |
---|
1925 | hom(:,1,57,sr) = sums(:,57) ! ( w"e + w"p"/rho_ocean )/dz |
---|
1926 | hom(:,1,58,sr) = sums(:,58) ! u"pt" |
---|
1927 | hom(:,1,59,sr) = sums(:,59) ! u*pt* |
---|
1928 | hom(:,1,60,sr) = sums(:,58) + sums(:,59) ! upt_t |
---|
1929 | hom(:,1,61,sr) = sums(:,61) ! v"pt" |
---|
1930 | hom(:,1,62,sr) = sums(:,62) ! v*pt* |
---|
1931 | hom(:,1,63,sr) = sums(:,61) + sums(:,62) ! vpt_t |
---|
1932 | hom(:,1,64,sr) = sums(:,64) ! rho_ocean |
---|
1933 | hom(:,1,65,sr) = sums(:,65) ! w"sa" |
---|
1934 | hom(:,1,66,sr) = sums(:,66) ! w*sa* |
---|
1935 | hom(:,1,67,sr) = sums(:,65) + sums(:,66) ! wsa |
---|
1936 | hom(:,1,68,sr) = sums(:,68) ! w*p* |
---|
1937 | hom(:,1,69,sr) = sums(:,69) ! w"e + w"p"/rho_ocean |
---|
1938 | hom(:,1,70,sr) = sums(:,70) ! q*2 |
---|
1939 | hom(:,1,71,sr) = sums(:,71) ! prho |
---|
1940 | hom(:,1,72,sr) = hyp * 1E-2_wp ! hyp in hPa |
---|
1941 | hom(:,1,123,sr) = sums(:,123) ! nc |
---|
1942 | hom(:,1,124,sr) = sums(:,124) ! ni |
---|
1943 | hom(:,1,125,sr) = sums(:,125) ! qi |
---|
1944 | hom(:,1,126,sr) = sums(:,126) ! ng |
---|
1945 | hom(:,1,127,sr) = sums(:,127) ! qg |
---|
1946 | hom(:,1,128,sr) = sums(:,128) ! ns |
---|
1947 | hom(:,1,129,sr) = sums(:,129) ! qs |
---|
1948 | hom(:,1,73,sr) = sums(:,73) ! nr |
---|
1949 | hom(:,1,74,sr) = sums(:,74) ! qr |
---|
1950 | hom(:,1,75,sr) = sums(:,75) ! qc |
---|
1951 | hom(:,1,76,sr) = sums(:,76) ! prr (precipitation rate) |
---|
1952 | ! 77 is initial density profile |
---|
1953 | hom(:,1,78,sr) = ug ! ug |
---|
1954 | hom(:,1,79,sr) = vg ! vg |
---|
1955 | hom(:,1,80,sr) = w_subs ! w_subs |
---|
1956 | |
---|
1957 | IF ( large_scale_forcing ) THEN |
---|
1958 | hom(:,1,81,sr) = sums_ls_l(:,0) ! td_lsa_lpt |
---|
1959 | hom(:,1,82,sr) = sums_ls_l(:,1) ! td_lsa_q |
---|
1960 | IF ( use_subsidence_tendencies ) THEN |
---|
1961 | hom(:,1,83,sr) = sums_ls_l(:,2) ! td_sub_lpt |
---|
1962 | hom(:,1,84,sr) = sums_ls_l(:,3) ! td_sub_q |
---|
1963 | ELSE |
---|
1964 | hom(:,1,83,sr) = sums(:,83) ! td_sub_lpt |
---|
1965 | hom(:,1,84,sr) = sums(:,84) ! td_sub_q |
---|
1966 | ENDIF |
---|
1967 | hom(:,1,85,sr) = sums(:,85) ! td_nud_lpt |
---|
1968 | hom(:,1,86,sr) = sums(:,86) ! td_nud_q |
---|
1969 | hom(:,1,87,sr) = sums(:,87) ! td_nud_u |
---|
1970 | hom(:,1,88,sr) = sums(:,88) ! td_nud_v |
---|
1971 | ENDIF |
---|
1972 | |
---|
1973 | IF ( land_surface ) THEN |
---|
1974 | hom(:,1,89,sr) = sums(:,89) ! t_soil |
---|
1975 | ! 90 is initial t_soil profile |
---|
1976 | hom(:,1,91,sr) = sums(:,91) ! m_soil |
---|
1977 | ! 92 is initial m_soil profile |
---|
1978 | hom(:,1,93,sr) = sums(:,93) ! ghf |
---|
1979 | hom(:,1,94,sr) = sums(:,94) ! qsws_liq |
---|
1980 | hom(:,1,95,sr) = sums(:,95) ! qsws_soil |
---|
1981 | hom(:,1,96,sr) = sums(:,96) ! qsws_veg |
---|
1982 | hom(:,1,97,sr) = sums(:,97) ! r_a |
---|
1983 | hom(:,1,98,sr) = sums(:,98) ! r_s |
---|
1984 | |
---|
1985 | ENDIF |
---|
1986 | |
---|
1987 | IF ( radiation ) THEN |
---|
1988 | hom(:,1,99,sr) = sums(:,99) ! rad_net |
---|
1989 | hom(:,1,100,sr) = sums(:,100) ! rad_lw_in |
---|
1990 | hom(:,1,101,sr) = sums(:,101) ! rad_lw_out |
---|
1991 | hom(:,1,102,sr) = sums(:,102) ! rad_sw_in |
---|
1992 | hom(:,1,103,sr) = sums(:,103) ! rad_sw_out |
---|
1993 | |
---|
1994 | IF ( radiation_scheme == 'rrtmg' ) THEN |
---|
1995 | hom(:,1,104,sr) = sums(:,104) ! rad_lw_cs_hr |
---|
1996 | hom(:,1,105,sr) = sums(:,105) ! rad_lw_hr |
---|
1997 | hom(:,1,106,sr) = sums(:,106) ! rad_sw_cs_hr |
---|
1998 | hom(:,1,107,sr) = sums(:,107) ! rad_sw_hr |
---|
1999 | |
---|
2000 | hom(:,1,108,sr) = sums(:,108) ! rrtm_aldif |
---|
2001 | hom(:,1,109,sr) = sums(:,109) ! rrtm_aldir |
---|
2002 | hom(:,1,110,sr) = sums(:,110) ! rrtm_asdif |
---|
2003 | hom(:,1,111,sr) = sums(:,111) ! rrtm_asdir |
---|
2004 | ENDIF |
---|
2005 | ENDIF |
---|
2006 | |
---|
2007 | hom(:,1,112,sr) = sums(:,112) !: L |
---|
2008 | |
---|
2009 | IF ( passive_scalar ) THEN |
---|
2010 | hom(:,1,117,sr) = sums(:,117) ! w"s" |
---|
2011 | hom(:,1,114,sr) = sums(:,114) ! w*s* |
---|
2012 | hom(:,1,118,sr) = sums(:,117) + sums(:,114) ! ws |
---|
2013 | hom(:,1,116,sr) = sums(:,116) ! s*2 |
---|
2014 | ENDIF |
---|
2015 | |
---|
2016 | hom(:,1,119,sr) = rho_air ! rho_air in Kg/m^3 |
---|
2017 | hom(:,1,120,sr) = rho_air_zw ! rho_air_zw in Kg/m^3 |
---|
2018 | |
---|
2019 | IF ( kolmogorov_length_scale ) THEN |
---|
2020 | hom(:,1,121,sr) = sums(:,121) * 1E3_wp ! eta in mm |
---|
2021 | ENDIF |
---|
2022 | |
---|
2023 | |
---|
2024 | hom(:,1,pr_palm,sr) = sums(:,pr_palm) |
---|
2025 | ! u*, w'u', w'v', t* (in last profile) |
---|
2026 | |
---|
2027 | IF ( max_pr_user > 0 ) THEN ! user-defined profiles |
---|
2028 | hom(:,1,pr_palm+1:pr_palm+max_pr_user,sr) = & |
---|
2029 | sums(:,pr_palm+1:pr_palm+max_pr_user) |
---|
2030 | ENDIF |
---|
2031 | |
---|
2032 | IF ( air_chemistry ) THEN |
---|
2033 | IF ( max_pr_cs > 0 ) THEN ! chem_spcs profiles |
---|
2034 | hom(:, 1, pr_palm+max_pr_user+1:pr_palm + max_pr_user+max_pr_cs, sr) = & |
---|
2035 | sums(:, pr_palm+max_pr_user+1:pr_palm+max_pr_user+max_pr_cs) |
---|
2036 | ENDIF |
---|
2037 | ENDIF |
---|
2038 | |
---|
2039 | IF ( salsa ) THEN |
---|
2040 | IF ( max_pr_salsa > 0 ) THEN ! salsa profiles |
---|
2041 | hom(:,1,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa, sr) = & |
---|
2042 | sums(:,pr_palm+max_pr_user+max_pr_cs+1:pr_palm+max_pr_user+max_pr_cs+max_pr_salsa) |
---|
2043 | ENDIF |
---|
2044 | ENDIF |
---|
2045 | ! |
---|
2046 | !-- Determine the boundary layer height using two different schemes. |
---|
2047 | !-- First scheme: Starting from the Earth's (Ocean's) surface, look for the first relative |
---|
2048 | !-- minimum (maximum) of the total heat flux. |
---|
2049 | !-- The corresponding height is assumed as the boundary layer height, if it is less than 1.5 |
---|
2050 | !-- times the height where the heat flux becomes negative (positive) for the first time. |
---|
2051 | !-- Attention: the resolved vertical sensible heat flux (hom(:,1,17,sr) = w*pt*) is not known at |
---|
2052 | !-- the beginning because the calculation happens in advec_s_ws which is called after |
---|
2053 | !-- flow_statistics. Therefore z_i is directly taken from restart data at the beginning of |
---|
2054 | !-- restart runs. |
---|
2055 | IF ( TRIM( initializing_actions ) /= 'read_restart_data' .OR. & |
---|
2056 | simulated_time_at_begin /= simulated_time ) THEN |
---|
2057 | |
---|
2058 | z_i(1) = 0.0_wp |
---|
2059 | first = .TRUE. |
---|
2060 | |
---|
2061 | IF ( ocean_mode ) THEN |
---|
2062 | DO k = nzt, nzb+1, -1 |
---|
2063 | IF ( first .AND. hom(k,1,18,sr) < -1.0E-8_wp ) THEN |
---|
2064 | first = .FALSE. |
---|
2065 | height = zw(k) |
---|
2066 | ENDIF |
---|
2067 | IF ( hom(k,1,18,sr) < -1.0E-8_wp .AND. hom(k-1,1,18,sr) > hom(k,1,18,sr) ) THEN |
---|
2068 | IF ( zw(k) < 1.5_wp * height ) THEN |
---|
2069 | z_i(1) = zw(k) |
---|
2070 | ELSE |
---|
2071 | z_i(1) = height |
---|
2072 | ENDIF |
---|
2073 | EXIT |
---|
2074 | ENDIF |
---|
2075 | ENDDO |
---|
2076 | ELSE |
---|
2077 | DO k = nzb, nzt-1 |
---|
2078 | IF ( first .AND. hom(k,1,18,sr) < -1.0E-8_wp ) THEN |
---|
2079 | first = .FALSE. |
---|
2080 | height = zw(k) |
---|
2081 | ENDIF |
---|
2082 | IF ( hom(k,1,18,sr) < -1.0E-8_wp .AND. hom(k+1,1,18,sr) > hom(k,1,18,sr) ) THEN |
---|
2083 | IF ( zw(k) < 1.5_wp * height ) THEN |
---|
2084 | z_i(1) = zw(k) |
---|
2085 | ELSE |
---|
2086 | z_i(1) = height |
---|
2087 | ENDIF |
---|
2088 | EXIT |
---|
2089 | ENDIF |
---|
2090 | ENDDO |
---|
2091 | ENDIF |
---|
2092 | |
---|
2093 | ! |
---|
2094 | !-- Second scheme: Gradient scheme from Sullivan et al. (1998), modified by Uhlenbrock(2006). |
---|
2095 | !-- The boundary layer height is the height with the maximal local temperature gradient: |
---|
2096 | !-- starting from the second (the last but one) vertical gridpoint, the local gradient must be |
---|
2097 | !-- at least 0.2K/100m and greater than the next four gradients. |
---|
2098 | !-- WARNING: The threshold value of 0.2K/100m must be adjusted for the |
---|
2099 | !-- ocean case! |
---|
2100 | z_i(2) = 0.0_wp |
---|
2101 | DO k = nzb+1, nzt+1 |
---|
2102 | dptdz(k) = ( hom(k,1,4,sr) - hom(k-1,1,4,sr) ) * ddzu(k) |
---|
2103 | ENDDO |
---|
2104 | dptdz_threshold = 0.2_wp / 100.0_wp |
---|
2105 | |
---|
2106 | IF ( ocean_mode ) THEN |
---|
2107 | DO k = nzt+1, nzb+5, -1 |
---|
2108 | IF ( dptdz(k) > dptdz_threshold .AND. & |
---|
2109 | dptdz(k) > dptdz(k-1) .AND. dptdz(k) > dptdz(k-2) .AND. & |
---|
2110 | dptdz(k) > dptdz(k-3) .AND. dptdz(k) > dptdz(k-4) ) THEN |
---|
2111 | z_i(2) = zw(k-1) |
---|
2112 | EXIT |
---|
2113 | ENDIF |
---|
2114 | ENDDO |
---|
2115 | ELSE |
---|
2116 | DO k = nzb+1, nzt-3 |
---|
2117 | IF ( dptdz(k) > dptdz_threshold .AND. & |
---|
2118 | dptdz(k) > dptdz(k+1) .AND. dptdz(k) > dptdz(k+2) .AND. & |
---|
2119 | dptdz(k) > dptdz(k+3) .AND. dptdz(k) > dptdz(k+4) ) THEN |
---|
2120 | z_i(2) = zw(k-1) |
---|
2121 | EXIT |
---|
2122 | ENDIF |
---|
2123 | ENDDO |
---|
2124 | ENDIF |
---|
2125 | |
---|
2126 | ENDIF |
---|
2127 | |
---|
2128 | hom(nzb+6,1,pr_palm,sr) = z_i(1) |
---|
2129 | hom(nzb+7,1,pr_palm,sr) = z_i(2) |
---|
2130 | |
---|
2131 | ! |
---|
2132 | !-- Determine vertical index which is nearest to the mean surface level height of the respective |
---|
2133 | !-- statistic region |
---|
2134 | DO k = nzb, nzt |
---|
2135 | IF ( zw(k) >= mean_surface_level_height(sr) ) THEN |
---|
2136 | k_surface_level = k |
---|
2137 | EXIT |
---|
2138 | ENDIF |
---|
2139 | ENDDO |
---|
2140 | |
---|
2141 | ! |
---|
2142 | !-- Computation of both the characteristic vertical velocity and the characteristic convective |
---|
2143 | !-- boundary layer temperature. |
---|
2144 | !-- The inversion height entering into the equation is defined with respect to the mean surface |
---|
2145 | !-- level height of the respective statistic region. |
---|
2146 | !-- The horizontal average at surface level index + 1 is input for the average temperature. |
---|
2147 | IF ( hom(k_surface_level,1,18,sr) > 1.0E-8_wp .AND. z_i(1) /= 0.0_wp ) THEN |
---|
2148 | hom(nzb+8,1,pr_palm,sr) = & |
---|
2149 | ( g / hom(k_surface_level+1,1,4,sr) * & |
---|
2150 | ( hom(k_surface_level,1,18,sr) / & |
---|
2151 | ( heatflux_output_conversion(nzb) * rho_air(nzb) ) ) & |
---|
2152 | * ABS( z_i(1) - mean_surface_level_height(sr) ) )**0.333333333_wp |
---|
2153 | ELSE |
---|
2154 | hom(nzb+8,1,pr_palm,sr) = 0.0_wp |
---|
2155 | ENDIF |
---|
2156 | |
---|
2157 | ! |
---|
2158 | !-- Collect the time series quantities. Please note, timeseries quantities which are collected |
---|
2159 | !-- from horizontally averaged profiles, e.g. wpt or pt(zp), are treated specially. In case of |
---|
2160 | !-- elevated model surfaces, index nzb+1 might be within topography and data will be zero. |
---|
2161 | !-- Therefore, take value for the first atmosphere index, which is topo_min_level+1. |
---|
2162 | ts_value(1,sr) = hom(nzb+4,1,pr_palm,sr) ! E |
---|
2163 | ts_value(2,sr) = hom(nzb+5,1,pr_palm,sr) ! E* |
---|
2164 | ts_value(3,sr) = dt_3d |
---|
2165 | ts_value(4,sr) = hom(nzb,1,pr_palm,sr) ! u* |
---|
2166 | ts_value(5,sr) = hom(nzb+3,1,pr_palm,sr) ! th* |
---|
2167 | ts_value(6,sr) = u_max |
---|
2168 | ts_value(7,sr) = v_max |
---|
2169 | ts_value(8,sr) = w_max |
---|
2170 | ts_value(9,sr) = hom(nzb+10,1,pr_palm,sr) ! new divergence |
---|
2171 | ts_value(10,sr) = hom(nzb+9,1,pr_palm,sr) ! old Divergence |
---|
2172 | ts_value(11,sr) = hom(nzb+6,1,pr_palm,sr) ! z_i(1) |
---|
2173 | ts_value(12,sr) = hom(nzb+7,1,pr_palm,sr) ! z_i(2) |
---|
2174 | ts_value(13,sr) = hom(nzb+8,1,pr_palm,sr) ! w* |
---|
2175 | ts_value(14,sr) = hom(nzb,1,16,sr) ! w'pt' at k=0 |
---|
2176 | ts_value(15,sr) = hom(topo_min_level+1,1,16,sr) ! w'pt' at k=1 |
---|
2177 | ts_value(16,sr) = hom(topo_min_level+1,1,18,sr) ! wpt at k=1 |
---|
2178 | ts_value(17,sr) = hom(nzb+14,1,pr_palm,sr) ! pt(0) |
---|
2179 | ts_value(18,sr) = hom(topo_min_level+1,1,4,sr) ! pt(zp) |
---|
2180 | ts_value(19,sr) = hom(nzb+1,1,pr_palm,sr) ! u'w' at k=0 |
---|
2181 | ts_value(20,sr) = hom(nzb+2,1,pr_palm,sr) ! v'w' at k=0 |
---|
2182 | ts_value(21,sr) = hom(nzb,1,48,sr) ! w"q" at k=0 |
---|
2183 | |
---|
2184 | IF ( .NOT. neutral ) THEN |
---|
2185 | ts_value(22,sr) = hom(nzb,1,112,sr) ! L |
---|
2186 | ELSE |
---|
2187 | ts_value(22,sr) = 1.0E10_wp |
---|
2188 | ENDIF |
---|
2189 | |
---|
2190 | ts_value(23,sr) = hom(nzb+12,1,pr_palm,sr) ! q* |
---|
2191 | |
---|
2192 | IF ( passive_scalar ) THEN |
---|
2193 | ts_value(24,sr) = hom(nzb+13,1,117,sr) ! w"s" ( to do ! ) |
---|
2194 | ts_value(25,sr) = hom(nzb+13,1,pr_palm,sr) ! s* |
---|
2195 | ENDIF |
---|
2196 | |
---|
2197 | ! |
---|
2198 | !-- Collect land surface model timeseries |
---|
2199 | IF ( land_surface ) THEN |
---|
2200 | ts_value(dots_soil ,sr) = hom(nzb,1,93,sr) ! ghf |
---|
2201 | ts_value(dots_soil+1,sr) = hom(nzb,1,94,sr) ! qsws_liq |
---|
2202 | ts_value(dots_soil+2,sr) = hom(nzb,1,95,sr) ! qsws_soil |
---|
2203 | ts_value(dots_soil+3,sr) = hom(nzb,1,96,sr) ! qsws_veg |
---|
2204 | ts_value(dots_soil+4,sr) = hom(nzb,1,97,sr) ! r_a |
---|
2205 | ts_value(dots_soil+5,sr) = hom(nzb,1,98,sr) ! r_s |
---|
2206 | ENDIF |
---|
2207 | ! |
---|
2208 | !-- Collect radiation model timeseries |
---|
2209 | IF ( radiation ) THEN |
---|
2210 | ts_value(dots_rad,sr) = hom(nzb,1,99,sr) ! rad_net |
---|
2211 | ts_value(dots_rad+1,sr) = hom(nzb,1,100,sr) ! rad_lw_in |
---|
2212 | ts_value(dots_rad+2,sr) = hom(nzb,1,101,sr) ! rad_lw_out |
---|
2213 | ts_value(dots_rad+3,sr) = hom(nzb,1,102,sr) ! rad_sw_in |
---|
2214 | ts_value(dots_rad+4,sr) = hom(nzb,1,103,sr) ! rad_sw_out |
---|
2215 | |
---|
2216 | IF ( radiation_scheme == 'rrtmg' ) THEN |
---|
2217 | ts_value(dots_rad+5,sr) = hom(nzb,1,108,sr) ! rrtm_aldif |
---|
2218 | ts_value(dots_rad+6,sr) = hom(nzb,1,109,sr) ! rrtm_aldir |
---|
2219 | ts_value(dots_rad+7,sr) = hom(nzb,1,110,sr) ! rrtm_asdif |
---|
2220 | ts_value(dots_rad+8,sr) = hom(nzb,1,111,sr) ! rrtm_asdir |
---|
2221 | ENDIF |
---|
2222 | |
---|
2223 | ENDIF |
---|
2224 | |
---|
2225 | ! |
---|
2226 | !-- Calculate additional statistics provided by other modules |
---|
2227 | CALL module_interface_statistics( 'time_series', sr, 0, dots_max ) |
---|
2228 | |
---|
2229 | ENDDO ! loop of the subregions |
---|
2230 | |
---|
2231 | ! |
---|
2232 | !-- If required, sum up horizontal averages for subsequent time averaging. |
---|
2233 | !-- Do not sum, if flow statistics is called before the first initial time step. |
---|
2234 | IF ( do_sum .AND. simulated_time /= 0.0_wp ) THEN |
---|
2235 | IF ( average_count_pr == 0 ) hom_sum = 0.0_wp |
---|
2236 | hom_sum = hom_sum + hom(:,1,:,:) |
---|
2237 | average_count_pr = average_count_pr + 1 |
---|
2238 | do_sum = .FALSE. |
---|
2239 | ENDIF |
---|
2240 | |
---|
2241 | ! |
---|
2242 | !-- Set flag for other UPs (e.g. output routines, but also buoyancy). |
---|
2243 | !-- This flag is reset after each time step in time_integration. |
---|
2244 | flow_statistics_called = .TRUE. |
---|
2245 | |
---|
2246 | CALL cpu_log( log_point(10), 'flow_statistics', 'stop' ) |
---|
2247 | |
---|
2248 | |
---|
2249 | END SUBROUTINE flow_statistics |
---|