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

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

Change order of dimension in surface arrays %frac, %emissivity and %albedo to allow for better vectorization in the radiation interactions; Set back turbulent length scale to 8 x grid spacing in the parametrized mode for the synthetic turbulence generator (was accidentally changed in last commit)

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