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

Last change on this file since 4463 was 4463, checked in by Giersch, 5 years ago

Calculation of horizontally averaged w profile moved to the calculations for the other horizontally averaged velocity profiles

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