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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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