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

Last change on this file since 4509 was 4502, checked in by schwenkel, 4 years ago

Implementation of ice microphysics

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