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

Last change on this file since 4893 was 4888, checked in by suehring, 4 years ago

Timeseries output of ghf and r_a as sum of all horizontally upward-facing surfaces in LSM+USM simulations; virtual measurements: Reduce number of additional grid point output in z for timeseries data to a reasonable number

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