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

Last change on this file since 4803 was 4757, checked in by schwenkel, 4 years ago

implement relative humidity as diagnostic output quantity

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