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

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

Bugfix in output of time-averaged plant-canopy quanities; Output of plant-canopy data only where tall canopy is defined; land-surface model: fix wrong location strings; tests: update urban test case; all source code files: copyright update

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