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

Last change on this file since 4470 was 4464, checked in by Giersch, 5 years ago

Reset last change (r4463)

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