source: palm/trunk/SOURCE/flow_statistics.f90 @ 4526

Last change on this file since 4526 was 4521, checked in by schwenkel, 5 years ago

add error number

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