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

Last change on this file since 4459 was 4444, checked in by raasch, 5 years ago

bugfix: cpp-directives for serial mode added

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