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

Last change on this file since 4558 was 4551, checked in by suehring, 4 years ago

Bugfix in summation for statistical regions

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