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

Last change on this file since 4443 was 4442, checked in by suehring, 4 years ago

last commit documented

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