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

Last change on this file since 4659 was 4646, checked in by raasch, 4 years ago

files re-formatted to follow the PALM coding standard

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