SUBROUTINE flow_statistics !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: flow_statistics.f90 1008 2012-09-19 14:49:14Z letzel $ ! ! 1007 2012-09-19 14:30:36Z franke ! Calculation of buoyancy flux for humidity in case of WS-scheme is now using ! turbulent fluxes of WS-scheme ! Bugfix: Calculation of subgridscale buoyancy flux for humidity and cloud ! droplets at nzb and nzt added ! ! 801 2012-01-10 17:30:36Z suehring ! Calculation of turbulent fluxes in advec_ws is now thread-safe. ! ! 743 2011-08-18 16:10:16Z suehring ! Calculation of turbulent fluxes with WS-scheme only for the whole model ! domain, not for user-defined subregions. ! ! 709 2011-03-30 09:31:40Z raasch ! formatting adjustments ! ! 699 2011-03-22 17:52:22Z suehring ! Bugfix in calculation of vertical velocity skewness. The added absolute value ! avoid negative values in the root. Negative values of w'w' can occur at the ! top or bottom of the model domain due to degrading the order of advection ! scheme. Furthermore the calculation will be the same for all advection ! schemes. ! ! 696 2011-03-18 07:03:49Z raasch ! Bugfix: Summation of Wicker-Skamarock scheme fluxes and variances for all ! threads ! ! 678 2011-02-02 14:31:56Z raasch ! Bugfix in calculation of the divergence of vertical flux of resolved scale ! energy, pressure fluctuations, and flux of pressure fluctuation itself ! ! 673 2011-01-18 16:19:48Z suehring ! Top bc for the horizontal velocity variances added for ocean runs. ! Setting the corresponding bottom bc moved to advec_ws. ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! When advection is computed with ws-scheme, turbulent fluxes are already ! computed in the respective advection routines and buffered in arrays ! sums_xx_ws_l(). This is due to a consistent treatment of statistics with the ! numerics and to avoid unphysical kinks near the surface. ! So some if requests has to be done to dicern between fluxes from ws-scheme ! other advection schemes. ! Furthermore the computation of z_i is only done if the heat flux exceeds a ! minimum value. This affects only simulations of a neutral boundary layer and ! is due to reasons of computations in the advection scheme. ! ! 624 2010-12-10 11:46:30Z heinze ! Calculation of q*2 added ! ! 622 2010-12-10 08:08:13Z raasch ! optional barriers included in order to speed up collective operations ! ! 388 2009-09-23 09:40:33Z raasch ! Vertical profiles of potential density and hydrostatic pressure are ! calculated. ! Added missing timeseries calculation of w"q"(0), moved timeseries q* to the ! end. ! Temperature gradient criterion for estimating the boundary layer height ! replaced by the gradient criterion of Sullivan et al. (1998). ! Output of messages replaced by message handling routine. ! ! 197 2008-09-16 15:29:03Z raasch ! Spline timeseries splptx etc. removed, timeseries w'u', w'v', w'q' (k=0) ! added, ! bugfix: divide sums(k,8) (e) and sums(k,34) (e*) by ngp_2dh_s_inner(k,sr) ! (like other scalars) ! ! 133 2007-11-20 10:10:53Z letzel ! Vertical profiles now based on nzb_s_inner; they are divided by ! ngp_2dh_s_inner (scalars, procucts of scalars) and ngp_2dh (staggered ! velocity components and their products, procucts of scalars and velocity ! components), respectively. ! ! 106 2007-08-16 14:30:26Z raasch ! Prescribed momentum fluxes at the top surface are used, ! profiles for w*p* and w"e are calculated ! ! 97 2007-06-21 08:23:15Z raasch ! Statistics for ocean version (salinity, density) added, ! calculation of z_i and Deardorff velocity scale adjusted to be used with ! the ocean version ! ! 87 2007-05-22 15:46:47Z raasch ! Two more arguments added to user_statistics, which is now also called for ! user-defined profiles, ! var_hom and var_sum renamed pr_palm ! ! 82 2007-04-16 15:40:52Z raasch ! Cpp-directive lcmuk changed to intel_openmp_bug ! ! 75 2007-03-22 09:54:05Z raasch ! Collection of time series quantities moved from routine flow_statistics to ! here, routine user_statistics is called for each statistic region, ! moisture renamed humidity ! ! 19 2007-02-23 04:53:48Z raasch ! fluxes at top modified (tswst, qswst) ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.41 2006/08/04 14:37:50 raasch ! Error removed in non-parallel part (sums_l) ! ! Revision 1.1 1997/08/11 06:15:17 raasch ! Initial revision ! ! ! Description: ! ------------ ! Compute average profiles and further average flow quantities for the different ! user-defined (sub-)regions. The region indexed 0 is the total model domain. ! ! NOTE: For simplicity, nzb_s_inner and nzb_diff_s_inner are being used as a ! ---- lower vertical index for k-loops for all variables, although strictly ! speaking the k-loops would have to be split up according to the staggered ! grid. However, this implies no error since staggered velocity components are ! zero at the walls and inside buildings. !------------------------------------------------------------------------------! USE arrays_3d USE cloud_parameters USE control_parameters USE cpulog USE grid_variables USE indices USE interfaces USE pegrid USE statistics IMPLICIT NONE INTEGER :: i, j, k, omp_get_thread_num, sr, tn LOGICAL :: first REAL :: dptdz_threshold, height, pts, sums_l_eper, sums_l_etot, ust, & ust2, u2, vst, vst2, v2, w2, z_i(2) REAL :: dptdz(nzb+1:nzt+1) REAL :: sums_ll(nzb:nzt+1,2) CALL cpu_log( log_point(10), 'flow_statistics', 'start' ) ! !-- To be on the safe side, check whether flow_statistics has already been !-- called once after the current time step IF ( flow_statistics_called ) THEN message_string = 'flow_statistics is called two times within one ' // & 'timestep' CALL message( 'flow_statistics', 'PA0190', 1, 2, 0, 6, 0 ) ENDIF ! !-- Compute statistics for each (sub-)region DO sr = 0, statistic_regions ! !-- Initialize (local) summation array sums_l = 0.0 ! !-- Store sums that have been computed in other subroutines in summation !-- array sums_l(:,11,:) = sums_l_l(:,sr,:) ! mixing length from diffusivities !-- WARNING: next line still has to be adjusted for OpenMP sums_l(:,21,0) = sums_wsts_bc_l(:,sr) ! heat flux from advec_s_bc sums_l(nzb+9,pr_palm,0) = sums_divold_l(sr) ! old divergence from pres sums_l(nzb+10,pr_palm,0) = sums_divnew_l(sr) ! new divergence from pres ! !-- Copy the turbulent quantities, evaluated in the advection routines to !-- the local array sums_l() for further computations IF ( ws_scheme_mom .AND. sr == 0 ) THEN ! !-- According to the Neumann bc for the horizontal velocity components, !-- the corresponding fluxes has to satisfiy the same bc. IF ( ocean ) THEN sums_us2_ws_l(nzt+1,:) = sums_us2_ws_l(nzt,:) sums_vs2_ws_l(nzt+1,:) = sums_vs2_ws_l(nzt,:) ENDIF DO i = 0, threads_per_task-1 ! !-- Swap the turbulent quantities evaluated in advec_ws. sums_l(:,13,i) = sums_wsus_ws_l(:,i) ! w*u* sums_l(:,15,i) = sums_wsvs_ws_l(:,i) ! w*v* sums_l(:,30,i) = sums_us2_ws_l(:,i) ! u*2 sums_l(:,31,i) = sums_vs2_ws_l(:,i) ! v*2 sums_l(:,32,i) = sums_ws2_ws_l(:,i) ! w*2 sums_l(:,34,i) = sums_l(:,34,i) + 0.5 * & ( sums_us2_ws_l(:,i) + sums_vs2_ws_l(:,i) + & sums_ws2_ws_l(:,i) ) ! e* DO k = nzb, nzt sums_l(nzb+5,pr_palm,i) = sums_l(nzb+5,pr_palm,i) + 0.5 * ( & sums_us2_ws_l(k,i) + & sums_vs2_ws_l(k,i) + & sums_ws2_ws_l(k,i) ) ENDDO ENDDO ENDIF IF ( ws_scheme_sca .AND. sr == 0 ) THEN DO i = 0, threads_per_task-1 sums_l(:,17,i) = sums_wspts_ws_l(:,i) ! w*pt* from advec_s_ws IF ( ocean ) sums_l(:,66,i) = sums_wssas_ws_l(:,i) ! w*sa* IF ( humidity .OR. passive_scalar ) sums_l(:,49,i) = & sums_wsqs_ws_l(:,i) !w*q* ENDDO ENDIF ! !-- Horizontally averaged profiles of horizontal velocities and temperature. !-- They must have been computed before, because they are already required !-- for other horizontal averages. tn = 0 !$OMP PARALLEL PRIVATE( i, j, k, tn ) #if defined( __intel_openmp_bug ) tn = omp_get_thread_num() #else !$ tn = omp_get_thread_num() #endif !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i), nzt+1 sums_l(k,1,tn) = sums_l(k,1,tn) + u(k,j,i) * rmask(j,i,sr) sums_l(k,2,tn) = sums_l(k,2,tn) + v(k,j,i) * rmask(j,i,sr) sums_l(k,4,tn) = sums_l(k,4,tn) + pt(k,j,i) * rmask(j,i,sr) ENDDO ENDDO ENDDO ! !-- Horizontally averaged profile of salinity IF ( ocean ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i), nzt+1 sums_l(k,23,tn) = sums_l(k,23,tn) + & sa(k,j,i) * rmask(j,i,sr) ENDDO ENDDO ENDDO ENDIF ! !-- Horizontally averaged profiles of virtual potential temperature, !-- total water content, specific humidity and liquid water potential !-- temperature IF ( humidity ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i), nzt+1 sums_l(k,44,tn) = sums_l(k,44,tn) + & vpt(k,j,i) * rmask(j,i,sr) sums_l(k,41,tn) = sums_l(k,41,tn) + & q(k,j,i) * rmask(j,i,sr) ENDDO ENDDO ENDDO IF ( cloud_physics ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i), nzt+1 sums_l(k,42,tn) = sums_l(k,42,tn) + & ( q(k,j,i) - ql(k,j,i) ) * rmask(j,i,sr) sums_l(k,43,tn) = sums_l(k,43,tn) + ( & pt(k,j,i) + l_d_cp*pt_d_t(k) * ql(k,j,i) & ) * rmask(j,i,sr) ENDDO ENDDO ENDDO ENDIF ENDIF ! !-- Horizontally averaged profiles of passive scalar IF ( passive_scalar ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i), nzt+1 sums_l(k,41,tn) = sums_l(k,41,tn) + q(k,j,i) * rmask(j,i,sr) ENDDO ENDDO ENDDO ENDIF !$OMP END PARALLEL ! !-- Summation of thread sums IF ( threads_per_task > 1 ) THEN DO i = 1, threads_per_task-1 sums_l(:,1,0) = sums_l(:,1,0) + sums_l(:,1,i) sums_l(:,2,0) = sums_l(:,2,0) + sums_l(:,2,i) sums_l(:,4,0) = sums_l(:,4,0) + sums_l(:,4,i) IF ( ocean ) THEN sums_l(:,23,0) = sums_l(:,23,0) + sums_l(:,23,i) ENDIF IF ( humidity ) THEN sums_l(:,41,0) = sums_l(:,41,0) + sums_l(:,41,i) sums_l(:,44,0) = sums_l(:,44,0) + sums_l(:,44,i) IF ( cloud_physics ) THEN sums_l(:,42,0) = sums_l(:,42,0) + sums_l(:,42,i) sums_l(:,43,0) = sums_l(:,43,0) + sums_l(:,43,i) ENDIF ENDIF IF ( passive_scalar ) THEN sums_l(:,41,0) = sums_l(:,41,0) + sums_l(:,41,i) ENDIF ENDDO ENDIF #if defined( __parallel ) ! !-- Compute total sum from local sums IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), nzt+2-nzb, MPI_REAL, & MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,2,0), sums(nzb,2), nzt+2-nzb, MPI_REAL, & MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,4,0), sums(nzb,4), nzt+2-nzb, MPI_REAL, & MPI_SUM, comm2d, ierr ) IF ( ocean ) THEN IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,23,0), sums(nzb,23), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) ENDIF IF ( humidity ) THEN IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,44,0), sums(nzb,44), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,41,0), sums(nzb,41), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) IF ( cloud_physics ) THEN IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,42,0), sums(nzb,42), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,43,0), sums(nzb,43), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) ENDIF ENDIF IF ( passive_scalar ) THEN IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,41,0), sums(nzb,41), nzt+2-nzb, & MPI_REAL, MPI_SUM, comm2d, ierr ) ENDIF #else sums(:,1) = sums_l(:,1,0) sums(:,2) = sums_l(:,2,0) sums(:,4) = sums_l(:,4,0) IF ( ocean ) sums(:,23) = sums_l(:,23,0) IF ( humidity ) THEN sums(:,44) = sums_l(:,44,0) sums(:,41) = sums_l(:,41,0) IF ( cloud_physics ) THEN sums(:,42) = sums_l(:,42,0) sums(:,43) = sums_l(:,43,0) ENDIF ENDIF IF ( passive_scalar ) sums(:,41) = sums_l(:,41,0) #endif ! !-- Final values are obtained by division by the total number of grid points !-- used for summation. After that store profiles. sums(:,1) = sums(:,1) / ngp_2dh(sr) sums(:,2) = sums(:,2) / ngp_2dh(sr) sums(:,4) = sums(:,4) / ngp_2dh_s_inner(:,sr) hom(:,1,1,sr) = sums(:,1) ! u hom(:,1,2,sr) = sums(:,2) ! v hom(:,1,4,sr) = sums(:,4) ! pt ! !-- Salinity IF ( ocean ) THEN sums(:,23) = sums(:,23) / ngp_2dh_s_inner(:,sr) hom(:,1,23,sr) = sums(:,23) ! sa ENDIF ! !-- Humidity and cloud parameters IF ( humidity ) THEN sums(:,44) = sums(:,44) / ngp_2dh_s_inner(:,sr) sums(:,41) = sums(:,41) / ngp_2dh_s_inner(:,sr) hom(:,1,44,sr) = sums(:,44) ! vpt hom(:,1,41,sr) = sums(:,41) ! qv (q) IF ( cloud_physics ) THEN sums(:,42) = sums(:,42) / ngp_2dh_s_inner(:,sr) sums(:,43) = sums(:,43) / ngp_2dh_s_inner(:,sr) hom(:,1,42,sr) = sums(:,42) ! qv hom(:,1,43,sr) = sums(:,43) ! pt ENDIF ENDIF ! !-- Passive scalar IF ( passive_scalar ) hom(:,1,41,sr) = sums(:,41) / & ngp_2dh_s_inner(:,sr) ! s (q) ! !-- Horizontally averaged profiles of the remaining prognostic variables, !-- variances, the total and the perturbation energy (single values in last !-- column of sums_l) and some diagnostic quantities. !-- NOTE: for simplicity, nzb_s_inner is used below, although strictly !-- ---- speaking the following k-loop would have to be split up and !-- rearranged according to the staggered grid. !-- However, this implies no error since staggered velocity components !-- are zero at the walls and inside buildings. tn = 0 #if defined( __intel_openmp_bug ) !$OMP PARALLEL PRIVATE( i, j, k, pts, sums_ll, sums_l_eper, sums_l_etot, & !$OMP tn, ust, ust2, u2, vst, vst2, v2, w2 ) tn = omp_get_thread_num() #else !$OMP PARALLEL PRIVATE( i, j, k, pts, sums_ll, sums_l_eper, sums_l_etot, tn, ust, ust2, u2, vst, vst2, v2, w2 ) !$ tn = omp_get_thread_num() #endif !$OMP DO DO i = nxl, nxr DO j = nys, nyn sums_l_etot = 0.0 DO k = nzb_s_inner(j,i), nzt+1 ! !-- Prognostic and diagnostic variables sums_l(k,3,tn) = sums_l(k,3,tn) + w(k,j,i) * rmask(j,i,sr) sums_l(k,8,tn) = sums_l(k,8,tn) + e(k,j,i) * rmask(j,i,sr) sums_l(k,9,tn) = sums_l(k,9,tn) + km(k,j,i) * rmask(j,i,sr) sums_l(k,10,tn) = sums_l(k,10,tn) + kh(k,j,i) * rmask(j,i,sr) sums_l(k,40,tn) = sums_l(k,40,tn) + p(k,j,i) sums_l(k,33,tn) = sums_l(k,33,tn) + & ( pt(k,j,i)-hom(k,1,4,sr) )**2 * rmask(j,i,sr) IF ( humidity ) THEN sums_l(k,70,tn) = sums_l(k,70,tn) + & ( q(k,j,i)-hom(k,1,41,sr) )**2 * rmask(j,i,sr) ENDIF ! !-- Higher moments !-- (Computation of the skewness of w further below) sums_l(k,38,tn) = sums_l(k,38,tn) + w(k,j,i)**3 * rmask(j,i,sr) sums_l_etot = sums_l_etot + & 0.5 * ( u(k,j,i)**2 + v(k,j,i)**2 + & w(k,j,i)**2 ) * rmask(j,i,sr) ENDDO ! !-- Total and perturbation energy for the total domain (being !-- collected in the last column of sums_l). Summation of these !-- quantities is seperated from the previous loop in order to !-- allow vectorization of that loop. sums_l(nzb+4,pr_palm,tn) = sums_l(nzb+4,pr_palm,tn) + sums_l_etot ! !-- 2D-arrays (being collected in the last column of sums_l) sums_l(nzb,pr_palm,tn) = sums_l(nzb,pr_palm,tn) + & us(j,i) * rmask(j,i,sr) sums_l(nzb+1,pr_palm,tn) = sums_l(nzb+1,pr_palm,tn) + & usws(j,i) * rmask(j,i,sr) sums_l(nzb+2,pr_palm,tn) = sums_l(nzb+2,pr_palm,tn) + & vsws(j,i) * rmask(j,i,sr) sums_l(nzb+3,pr_palm,tn) = sums_l(nzb+3,pr_palm,tn) + & ts(j,i) * rmask(j,i,sr) IF ( humidity ) THEN sums_l(nzb+12,pr_palm,tn) = sums_l(nzb+12,pr_palm,tn) + & qs(j,i) * rmask(j,i,sr) ENDIF ENDDO ENDDO ! !-- Computation of statistics when ws-scheme is not used. Else these !-- quantities are evaluated in the advection routines. IF ( .NOT. ws_scheme_mom .OR. sr /= 0 ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn sums_l_eper = 0.0 DO k = nzb_s_inner(j,i), nzt+1 u2 = u(k,j,i)**2 v2 = v(k,j,i)**2 w2 = w(k,j,i)**2 ust2 = ( u(k,j,i) - hom(k,1,1,sr) )**2 vst2 = ( v(k,j,i) - hom(k,1,2,sr) )**2 sums_l(k,30,tn) = sums_l(k,30,tn) + ust2 * rmask(j,i,sr) sums_l(k,31,tn) = sums_l(k,31,tn) + vst2 * rmask(j,i,sr) sums_l(k,32,tn) = sums_l(k,32,tn) + w2 * rmask(j,i,sr) ! !-- Perturbation energy sums_l(k,34,tn) = sums_l(k,34,tn) + 0.5 * & ( ust2 + vst2 + w2 ) * rmask(j,i,sr) sums_l_eper = sums_l_eper + & 0.5 * ( ust2+vst2+w2 ) * rmask(j,i,sr) ENDDO sums_l(nzb+5,pr_palm,tn) = sums_l(nzb+5,pr_palm,tn) & + sums_l_eper ENDDO ENDDO ENDIF ! !-- Horizontally averaged profiles of the vertical fluxes !$OMP DO DO i = nxl, nxr DO j = nys, nyn ! !-- Subgridscale fluxes (without Prandtl layer from k=nzb, !-- oterwise from k=nzb+1) !-- NOTE: for simplicity, nzb_diff_s_inner is used below, although !-- ---- strictly speaking the following k-loop would have to be !-- split up according to the staggered grid. !-- However, this implies no error since staggered velocity !-- components are zero at the walls and inside buildings. DO k = nzb_diff_s_inner(j,i)-1, nzt_diff ! !-- Momentum flux w"u" sums_l(k,12,tn) = sums_l(k,12,tn) - 0.25 * ( & km(k,j,i)+km(k+1,j,i)+km(k,j,i-1)+km(k+1,j,i-1) & ) * ( & ( u(k+1,j,i) - u(k,j,i) ) * ddzu(k+1) & + ( w(k,j,i) - w(k,j,i-1) ) * ddx & ) * rmask(j,i,sr) ! !-- Momentum flux w"v" sums_l(k,14,tn) = sums_l(k,14,tn) - 0.25 * ( & km(k,j,i)+km(k+1,j,i)+km(k,j-1,i)+km(k+1,j-1,i) & ) * ( & ( v(k+1,j,i) - v(k,j,i) ) * ddzu(k+1) & + ( w(k,j,i) - w(k,j-1,i) ) * ddy & ) * rmask(j,i,sr) ! !-- Heat flux w"pt" sums_l(k,16,tn) = sums_l(k,16,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( pt(k+1,j,i) - pt(k,j,i) ) & * ddzu(k+1) * rmask(j,i,sr) ! !-- Salinity flux w"sa" IF ( ocean ) THEN sums_l(k,65,tn) = sums_l(k,65,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( sa(k+1,j,i) - sa(k,j,i) ) & * ddzu(k+1) * rmask(j,i,sr) ENDIF ! !-- Buoyancy flux, water flux (humidity flux) w"q" IF ( humidity ) THEN sums_l(k,45,tn) = sums_l(k,45,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( vpt(k+1,j,i) - vpt(k,j,i) ) & * ddzu(k+1) * rmask(j,i,sr) sums_l(k,48,tn) = sums_l(k,48,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( q(k+1,j,i) - q(k,j,i) ) & * ddzu(k+1) * rmask(j,i,sr) IF ( cloud_physics ) THEN sums_l(k,51,tn) = sums_l(k,51,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( ( q(k+1,j,i) - ql(k+1,j,i) )& - ( q(k,j,i) - ql(k,j,i) ) ) & * ddzu(k+1) * rmask(j,i,sr) ENDIF ENDIF ! !-- Passive scalar flux IF ( passive_scalar ) THEN sums_l(k,48,tn) = sums_l(k,48,tn) & - 0.5 * ( kh(k,j,i) + kh(k+1,j,i) ) & * ( q(k+1,j,i) - q(k,j,i) ) & * ddzu(k+1) * rmask(j,i,sr) ENDIF ENDDO ! !-- Subgridscale fluxes in the Prandtl layer IF ( use_surface_fluxes ) THEN sums_l(nzb,12,tn) = sums_l(nzb,12,tn) + & usws(j,i) * rmask(j,i,sr) ! w"u" sums_l(nzb,14,tn) = sums_l(nzb,14,tn) + & vsws(j,i) * rmask(j,i,sr) ! w"v" sums_l(nzb,16,tn) = sums_l(nzb,16,tn) + & shf(j,i) * rmask(j,i,sr) ! w"pt" sums_l(nzb,58,tn) = sums_l(nzb,58,tn) + & 0.0 * rmask(j,i,sr) ! u"pt" sums_l(nzb,61,tn) = sums_l(nzb,61,tn) + & 0.0 * rmask(j,i,sr) ! v"pt" IF ( ocean ) THEN sums_l(nzb,65,tn) = sums_l(nzb,65,tn) + & saswsb(j,i) * rmask(j,i,sr) ! w"sa" ENDIF IF ( humidity ) THEN sums_l(nzb,48,tn) = sums_l(nzb,48,tn) + & qsws(j,i) * rmask(j,i,sr) ! w"q" (w"qv") sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & ( 1.0 + 0.61 * q(nzb,j,i) ) * & shf(j,i) + 0.61 * pt(nzb,j,i) * & qsws(j,i) ) IF ( cloud_droplets ) THEN sums_l(nzb,45,tn) = sums_l(nzb,45,tn) + ( & ( 1.0 + 0.61 * q(nzb,j,i) - & ql(nzb,j,i) ) * shf(j,i) + & 0.61 * pt(nzb,j,i) * qsws(j,i) ) ENDIF IF ( cloud_physics ) THEN ! !-- Formula does not work if ql(nzb) /= 0.0 sums_l(nzb,51,tn) = sums_l(nzb,51,tn) + & ! w"q" (w"qv") qsws(j,i) * rmask(j,i,sr) ENDIF ENDIF IF ( passive_scalar ) THEN sums_l(nzb,48,tn) = sums_l(nzb,48,tn) + & qsws(j,i) * rmask(j,i,sr) ! w"q" (w"qv") ENDIF ENDIF ! !-- Subgridscale fluxes at the top surface IF ( use_top_fluxes ) THEN sums_l(nzt:nzt+1,12,tn) = sums_l(nzt:nzt+1,12,tn) + & uswst(j,i) * rmask(j,i,sr) ! w"u" sums_l(nzt:nzt+1,14,tn) = sums_l(nzt:nzt+1,14,tn) + & vswst(j,i) * rmask(j,i,sr) ! w"v" sums_l(nzt:nzt+1,16,tn) = sums_l(nzt:nzt+1,16,tn) + & tswst(j,i) * rmask(j,i,sr) ! w"pt" sums_l(nzt:nzt+1,58,tn) = sums_l(nzt:nzt+1,58,tn) + & 0.0 * rmask(j,i,sr) ! u"pt" sums_l(nzt:nzt+1,61,tn) = sums_l(nzt:nzt+1,61,tn) + & 0.0 * rmask(j,i,sr) ! v"pt" IF ( ocean ) THEN sums_l(nzt,65,tn) = sums_l(nzt,65,tn) + & saswst(j,i) * rmask(j,i,sr) ! w"sa" ENDIF IF ( humidity ) THEN sums_l(nzt,48,tn) = sums_l(nzt,48,tn) + & qswst(j,i) * rmask(j,i,sr) ! w"q" (w"qv") sums_l(nzt,45,tn) = sums_l(nzt,45,tn) + ( & ( 1.0 + 0.61 * q(nzt,j,i) ) * & tswst(j,i) + 0.61 * pt(nzt,j,i) * & qswst(j,i) ) IF ( cloud_droplets ) THEN sums_l(nzt,45,tn) = sums_l(nzt,45,tn) + ( & ( 1.0 + 0.61 * q(nzt,j,i) - & ql(nzt,j,i) ) * tswst(j,i) + & 0.61 * pt(nzt,j,i) * qswst(j,i) ) ENDIF IF ( cloud_physics ) THEN ! !-- Formula does not work if ql(nzb) /= 0.0 sums_l(nzt,51,tn) = sums_l(nzt,51,tn) + & ! w"q" (w"qv") qswst(j,i) * rmask(j,i,sr) ENDIF ENDIF IF ( passive_scalar ) THEN sums_l(nzt,48,tn) = sums_l(nzt,48,tn) + & qswst(j,i) * rmask(j,i,sr) ! w"q" (w"qv") ENDIF ENDIF ! !-- Resolved fluxes (can be computed for all horizontal points) !-- NOTE: for simplicity, nzb_s_inner is used below, although strictly !-- ---- speaking the following k-loop would have to be split up and !-- rearranged according to the staggered grid. DO k = nzb_s_inner(j,i), nzt ust = 0.5 * ( u(k,j,i) - hom(k,1,1,sr) + & u(k+1,j,i) - hom(k+1,1,1,sr) ) vst = 0.5 * ( v(k,j,i) - hom(k,1,2,sr) + & v(k+1,j,i) - hom(k+1,1,2,sr) ) pts = 0.5 * ( pt(k,j,i) - hom(k,1,4,sr) + & pt(k+1,j,i) - hom(k+1,1,4,sr) ) !-- Higher moments sums_l(k,35,tn) = sums_l(k,35,tn) + pts * w(k,j,i)**2 * & rmask(j,i,sr) sums_l(k,36,tn) = sums_l(k,36,tn) + pts**2 * w(k,j,i) * & rmask(j,i,sr) ! !-- Salinity flux and density (density does not belong to here, !-- but so far there is no other suitable place to calculate) IF ( ocean ) THEN IF( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN pts = 0.5 * ( sa(k,j,i) - hom(k,1,23,sr) + & sa(k+1,j,i) - hom(k+1,1,23,sr) ) sums_l(k,66,tn) = sums_l(k,66,tn) + pts * w(k,j,i) * & rmask(j,i,sr) ENDIF sums_l(k,64,tn) = sums_l(k,64,tn) + rho(k,j,i) * & rmask(j,i,sr) sums_l(k,71,tn) = sums_l(k,71,tn) + prho(k,j,i) * & rmask(j,i,sr) ENDIF ! !-- Buoyancy flux, water flux, humidity flux and liquid water !-- content IF ( humidity ) THEN IF ( cloud_physics .OR. cloud_droplets ) THEN pts = 0.5 * ( vpt(k,j,i) - hom(k,1,44,sr) + & vpt(k+1,j,i) - hom(k+1,1,44,sr) ) sums_l(k,46,tn) = sums_l(k,46,tn) + pts * w(k,j,i) * & rmask(j,i,sr) sums_l(k,54,tn) = sums_l(k,54,tn) + ql(k,j,i) * & rmask(j,i,sr) ELSE IF( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN pts = 0.5 * ( vpt(k,j,i) - hom(k,1,44,sr) + & vpt(k+1,j,i) - hom(k+1,1,44,sr) ) sums_l(k,46,tn) = sums_l(k,46,tn) + pts * w(k,j,i) * & rmask(j,i,sr) ELSE IF ( ws_scheme_sca .AND. sr == 0 ) THEN sums_l(k,46,tn) = ( 1.0 + 0.61 * hom(k,1,41,sr) ) * & sums_l(k,17,tn) + & 0.61 * hom(k,1,4,sr) * sums_l(k,49,tn) END IF END IF ENDIF ! !-- Passive scalar flux IF ( passive_scalar .AND. ( .NOT. ws_scheme_sca & .OR. sr /= 0 ) ) THEN pts = 0.5 * ( q(k,j,i) - hom(k,1,41,sr) + & q(k+1,j,i) - hom(k+1,1,41,sr) ) sums_l(k,49,tn) = sums_l(k,49,tn) + pts * w(k,j,i) * & rmask(j,i,sr) ENDIF ! !-- Energy flux w*e* !-- has to be adjusted sums_l(k,37,tn) = sums_l(k,37,tn) + w(k,j,i) * 0.5 * & ( ust**2 + vst**2 + w(k,j,i)**2 )& * rmask(j,i,sr) ENDDO ENDDO ENDDO ! !-- For speed optimization fluxes which have been computed in part directly !-- inside the WS advection routines are treated seperatly !-- Momentum fluxes first: IF ( .NOT. ws_scheme_mom .OR. sr /= 0 ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_diff_s_inner(j,i)-1, nzt_diff ust = 0.5 * ( u(k,j,i) - hom(k,1,1,sr) + & u(k+1,j,i) - hom(k+1,1,1,sr) ) vst = 0.5 * ( v(k,j,i) - hom(k,1,2,sr) + & v(k+1,j,i) - hom(k+1,1,2,sr) ) ! !-- Momentum flux w*u* sums_l(k,13,tn) = sums_l(k,13,tn) + 0.5 * & ( w(k,j,i-1) + w(k,j,i) ) & * ust * rmask(j,i,sr) ! !-- Momentum flux w*v* sums_l(k,15,tn) = sums_l(k,15,tn) + 0.5 * & ( w(k,j-1,i) + w(k,j,i) ) & * vst * rmask(j,i,sr) ENDDO ENDDO ENDDO ENDIF IF ( .NOT. ws_scheme_sca .OR. sr /= 0 ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_diff_s_inner(j,i)-1, nzt_diff ! !-- Vertical heat flux sums_l(k,17,tn) = sums_l(k,17,tn) + 0.5 * & ( pt(k,j,i) - hom(k,1,4,sr) + & pt(k+1,j,i) - hom(k+1,1,4,sr) ) & * w(k,j,i) * rmask(j,i,sr) IF ( humidity ) THEN pts = 0.5 * ( q(k,j,i) - hom(k,1,41,sr) + & q(k+1,j,i) - hom(k+1,1,41,sr) ) sums_l(k,49,tn) = sums_l(k,49,tn) + pts * w(k,j,i) * & rmask(j,i,sr) ENDIF ENDDO ENDDO ENDDO ENDIF ! !-- Density at top follows Neumann condition IF ( ocean ) THEN sums_l(nzt+1,64,tn) = sums_l(nzt,64,tn) sums_l(nzt+1,71,tn) = sums_l(nzt,71,tn) ENDIF ! !-- Divergence of vertical flux of resolved scale energy and pressure !-- fluctuations as well as flux of pressure fluctuation itself (68). !-- First calculate the products, then the divergence. !-- Calculation is time consuming. Do it only, if profiles shall be plotted. IF ( hom(nzb+1,2,55,0) /= 0.0 .OR. hom(nzb+1,2,68,0) /= 0.0 ) THEN sums_ll = 0.0 ! local array !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt sums_ll(k,1) = sums_ll(k,1) + 0.5 * w(k,j,i) * ( & ( 0.25 * ( u(k,j,i)+u(k+1,j,i)+u(k,j,i+1)+u(k+1,j,i+1) & - 0.5 * ( hom(k,1,1,sr) + hom(k+1,1,1,sr) ) & ) )**2 & + ( 0.25 * ( v(k,j,i)+v(k+1,j,i)+v(k,j+1,i)+v(k+1,j+1,i) & - 0.5 * ( hom(k,1,2,sr) + hom(k+1,1,2,sr) ) & ) )**2 & + w(k,j,i)**2 ) sums_ll(k,2) = sums_ll(k,2) + 0.5 * w(k,j,i) & * ( p(k,j,i) + p(k+1,j,i) ) ENDDO ENDDO ENDDO sums_ll(0,1) = 0.0 ! because w is zero at the bottom sums_ll(nzt+1,1) = 0.0 sums_ll(0,2) = 0.0 sums_ll(nzt+1,2) = 0.0 DO k = nzb+1, nzt sums_l(k,55,tn) = ( sums_ll(k,1) - sums_ll(k-1,1) ) * ddzw(k) sums_l(k,56,tn) = ( sums_ll(k,2) - sums_ll(k-1,2) ) * ddzw(k) sums_l(k,68,tn) = sums_ll(k,2) ENDDO sums_l(nzb,55,tn) = sums_l(nzb+1,55,tn) sums_l(nzb,56,tn) = sums_l(nzb+1,56,tn) sums_l(nzb,68,tn) = 0.0 ! because w* = 0 at nzb ENDIF ! !-- Divergence of vertical flux of SGS TKE and the flux itself (69) IF ( hom(nzb+1,2,57,0) /= 0.0 .OR. hom(nzb+1,2,69,0) /= 0.0 ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt sums_l(k,57,tn) = sums_l(k,57,tn) - 0.5 * ( & (km(k,j,i)+km(k+1,j,i)) * (e(k+1,j,i)-e(k,j,i)) * ddzu(k+1) & - (km(k-1,j,i)+km(k,j,i)) * (e(k,j,i)-e(k-1,j,i)) * ddzu(k) & ) * ddzw(k) sums_l(k,69,tn) = sums_l(k,69,tn) - 0.5 * ( & (km(k,j,i)+km(k+1,j,i)) * (e(k+1,j,i)-e(k,j,i)) * ddzu(k+1) & ) ENDDO ENDDO ENDDO sums_l(nzb,57,tn) = sums_l(nzb+1,57,tn) sums_l(nzb,69,tn) = sums_l(nzb+1,69,tn) ENDIF ! !-- Horizontal heat fluxes (subgrid, resolved, total). !-- Do it only, if profiles shall be plotted. IF ( hom(nzb+1,2,58,0) /= 0.0 ) THEN !$OMP DO DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt ! !-- Subgrid horizontal heat fluxes u"pt", v"pt" sums_l(k,58,tn) = sums_l(k,58,tn) - 0.5 * & ( kh(k,j,i) + kh(k,j,i-1) ) & * ( pt(k,j,i-1) - pt(k,j,i) ) & * ddx * rmask(j,i,sr) sums_l(k,61,tn) = sums_l(k,61,tn) - 0.5 * & ( kh(k,j,i) + kh(k,j-1,i) ) & * ( pt(k,j-1,i) - pt(k,j,i) ) & * ddy * rmask(j,i,sr) ! !-- Resolved horizontal heat fluxes u*pt*, v*pt* sums_l(k,59,tn) = sums_l(k,59,tn) + & ( u(k,j,i) - hom(k,1,1,sr) ) & * 0.5 * ( pt(k,j,i-1) - hom(k,1,4,sr) + & pt(k,j,i) - hom(k,1,4,sr) ) pts = 0.5 * ( pt(k,j-1,i) - hom(k,1,4,sr) + & pt(k,j,i) - hom(k,1,4,sr) ) sums_l(k,62,tn) = sums_l(k,62,tn) + & ( v(k,j,i) - hom(k,1,2,sr) ) & * 0.5 * ( pt(k,j-1,i) - hom(k,1,4,sr) + & pt(k,j,i) - hom(k,1,4,sr) ) ENDDO ENDDO ENDDO ! !-- Fluxes at the surface must be zero (e.g. due to the Prandtl-layer) sums_l(nzb,58,tn) = 0.0 sums_l(nzb,59,tn) = 0.0 sums_l(nzb,60,tn) = 0.0 sums_l(nzb,61,tn) = 0.0 sums_l(nzb,62,tn) = 0.0 sums_l(nzb,63,tn) = 0.0 ENDIF ! !-- Calculate the user-defined profiles CALL user_statistics( 'profiles', sr, tn ) !$OMP END PARALLEL ! !-- Summation of thread sums IF ( threads_per_task > 1 ) THEN DO i = 1, threads_per_task-1 sums_l(:,3,0) = sums_l(:,3,0) + sums_l(:,3,i) sums_l(:,4:40,0) = sums_l(:,4:40,0) + sums_l(:,4:40,i) sums_l(:,45:pr_palm,0) = sums_l(:,45:pr_palm,0) + & sums_l(:,45:pr_palm,i) IF ( max_pr_user > 0 ) THEN sums_l(:,pr_palm+1:pr_palm+max_pr_user,0) = & sums_l(:,pr_palm+1:pr_palm+max_pr_user,0) + & sums_l(:,pr_palm+1:pr_palm+max_pr_user,i) ENDIF ENDDO ENDIF #if defined( __parallel ) ! !-- Compute total sum from local sums IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), ngp_sums, MPI_REAL, & MPI_SUM, comm2d, ierr ) #else sums = sums_l(:,:,0) #endif ! !-- Final values are obtained by division by the total number of grid points !-- used for summation. After that store profiles. !-- Profiles: DO k = nzb, nzt+1 sums(k,3) = sums(k,3) / ngp_2dh(sr) sums(k,8:11) = sums(k,8:11) / ngp_2dh_s_inner(k,sr) sums(k,12:22) = sums(k,12:22) / ngp_2dh(sr) sums(k,23:29) = sums(k,23:29) / ngp_2dh_s_inner(k,sr) sums(k,30:32) = sums(k,30:32) / ngp_2dh(sr) sums(k,33:34) = sums(k,33:34) / ngp_2dh_s_inner(k,sr) sums(k,35:39) = sums(k,35:39) / ngp_2dh(sr) sums(k,40) = sums(k,40) / ngp_2dh_s_inner(k,sr) sums(k,45:53) = sums(k,45:53) / ngp_2dh(sr) sums(k,54) = sums(k,54) / ngp_2dh_s_inner(k,sr) sums(k,55:63) = sums(k,55:63) / ngp_2dh(sr) sums(k,64) = sums(k,64) / ngp_2dh_s_inner(k,sr) sums(k,65:69) = sums(k,65:69) / ngp_2dh(sr) sums(k,70:pr_palm-2) = sums(k,70:pr_palm-2)/ ngp_2dh_s_inner(k,sr) ENDDO !-- Upstream-parts sums(nzb:nzb+11,pr_palm-1) = sums(nzb:nzb+11,pr_palm-1) / ngp_3d(sr) !-- u* and so on !-- As sums(nzb:nzb+3,pr_palm) are full 2D arrays (us, usws, vsws, ts) whose !-- size is always ( nx + 1 ) * ( ny + 1 ), defined at the first grid layer !-- above the topography, they are being divided by ngp_2dh(sr) sums(nzb:nzb+3,pr_palm) = sums(nzb:nzb+3,pr_palm) / & ngp_2dh(sr) sums(nzb+12,pr_palm) = sums(nzb+12,pr_palm) / & ! qs ngp_2dh(sr) !-- eges, e* sums(nzb+4:nzb+5,pr_palm) = sums(nzb+4:nzb+5,pr_palm) / & ngp_3d(sr) !-- Old and new divergence sums(nzb+9:nzb+10,pr_palm) = sums(nzb+9:nzb+10,pr_palm) / & ngp_3d_inner(sr) !-- User-defined profiles IF ( max_pr_user > 0 ) THEN DO k = nzb, nzt+1 sums(k,pr_palm+1:pr_palm+max_pr_user) = & sums(k,pr_palm+1:pr_palm+max_pr_user) / & ngp_2dh_s_inner(k,sr) ENDDO ENDIF ! !-- Collect horizontal average in hom. !-- Compute deduced averages (e.g. total heat flux) hom(:,1,3,sr) = sums(:,3) ! w hom(:,1,8,sr) = sums(:,8) ! e profiles 5-7 are initial profiles hom(:,1,9,sr) = sums(:,9) ! km hom(:,1,10,sr) = sums(:,10) ! kh hom(:,1,11,sr) = sums(:,11) ! l hom(:,1,12,sr) = sums(:,12) ! w"u" hom(:,1,13,sr) = sums(:,13) ! w*u* hom(:,1,14,sr) = sums(:,14) ! w"v" hom(:,1,15,sr) = sums(:,15) ! w*v* hom(:,1,16,sr) = sums(:,16) ! w"pt" hom(:,1,17,sr) = sums(:,17) ! w*pt* hom(:,1,18,sr) = sums(:,16) + sums(:,17) ! wpt hom(:,1,19,sr) = sums(:,12) + sums(:,13) ! wu hom(:,1,20,sr) = sums(:,14) + sums(:,15) ! wv hom(:,1,21,sr) = sums(:,21) ! w*pt*BC hom(:,1,22,sr) = sums(:,16) + sums(:,21) ! wptBC ! profile 24 is initial profile (sa) ! profiles 25-29 left empty for initial ! profiles hom(:,1,30,sr) = sums(:,30) ! u*2 hom(:,1,31,sr) = sums(:,31) ! v*2 hom(:,1,32,sr) = sums(:,32) ! w*2 hom(:,1,33,sr) = sums(:,33) ! pt*2 hom(:,1,34,sr) = sums(:,34) ! e* hom(:,1,35,sr) = sums(:,35) ! w*2pt* hom(:,1,36,sr) = sums(:,36) ! w*pt*2 hom(:,1,37,sr) = sums(:,37) ! w*e* hom(:,1,38,sr) = sums(:,38) ! w*3 hom(:,1,39,sr) = sums(:,38) / ( abs( sums(:,32) ) + 1E-20 )**1.5 ! Sw hom(:,1,40,sr) = sums(:,40) ! p hom(:,1,45,sr) = sums(:,45) ! w"vpt" hom(:,1,46,sr) = sums(:,46) ! w*vpt* hom(:,1,47,sr) = sums(:,45) + sums(:,46) ! wvpt hom(:,1,48,sr) = sums(:,48) ! w"q" (w"qv") hom(:,1,49,sr) = sums(:,49) ! w*q* (w*qv*) hom(:,1,50,sr) = sums(:,48) + sums(:,49) ! wq (wqv) hom(:,1,51,sr) = sums(:,51) ! w"qv" hom(:,1,52,sr) = sums(:,52) ! w*qv* hom(:,1,53,sr) = sums(:,52) + sums(:,51) ! wq (wqv) hom(:,1,54,sr) = sums(:,54) ! ql hom(:,1,55,sr) = sums(:,55) ! w*u*u*/dz hom(:,1,56,sr) = sums(:,56) ! w*p*/dz hom(:,1,57,sr) = sums(:,57) ! ( w"e + w"p"/rho )/dz hom(:,1,58,sr) = sums(:,58) ! u"pt" hom(:,1,59,sr) = sums(:,59) ! u*pt* hom(:,1,60,sr) = sums(:,58) + sums(:,59) ! upt_t hom(:,1,61,sr) = sums(:,61) ! v"pt" hom(:,1,62,sr) = sums(:,62) ! v*pt* hom(:,1,63,sr) = sums(:,61) + sums(:,62) ! vpt_t hom(:,1,64,sr) = sums(:,64) ! rho hom(:,1,65,sr) = sums(:,65) ! w"sa" hom(:,1,66,sr) = sums(:,66) ! w*sa* hom(:,1,67,sr) = sums(:,65) + sums(:,66) ! wsa hom(:,1,68,sr) = sums(:,68) ! w*p* hom(:,1,69,sr) = sums(:,69) ! w"e + w"p"/rho hom(:,1,70,sr) = sums(:,70) ! q*2 hom(:,1,71,sr) = sums(:,71) ! prho hom(:,1,72,sr) = hyp * 1E-4 ! hyp in dbar hom(:,1,pr_palm-1,sr) = sums(:,pr_palm-1) ! upstream-parts u_x, u_y, u_z, v_x, ! v_y, usw. (in last but one profile) hom(:,1,pr_palm,sr) = sums(:,pr_palm) ! u*, w'u', w'v', t* (in last profile) IF ( max_pr_user > 0 ) THEN ! user-defined profiles hom(:,1,pr_palm+1:pr_palm+max_pr_user,sr) = & sums(:,pr_palm+1:pr_palm+max_pr_user) ENDIF ! !-- Determine the boundary layer height using two different schemes. !-- First scheme: Starting from the Earth's (Ocean's) surface, look for the !-- first relative minimum (maximum) of the total heat flux. !-- The corresponding height is assumed as the boundary layer height, if it !-- is less than 1.5 times the height where the heat flux becomes negative !-- (positive) for the first time. z_i(1) = 0.0 first = .TRUE. IF ( ocean ) THEN DO k = nzt, nzb+1, -1 IF ( first .AND. hom(k,1,18,sr) < 0.0 & .AND. abs(hom(k,1,18,sr)) > 1.0E-8) THEN first = .FALSE. height = zw(k) ENDIF IF ( hom(k,1,18,sr) < 0.0 .AND. & abs(hom(k,1,18,sr)) > 1.0E-8 .AND. & hom(k-1,1,18,sr) > hom(k,1,18,sr) ) THEN IF ( zw(k) < 1.5 * height ) THEN z_i(1) = zw(k) ELSE z_i(1) = height ENDIF EXIT ENDIF ENDDO ELSE DO k = nzb, nzt-1 IF ( first .AND. hom(k,1,18,sr) < 0.0 & .AND. abs(hom(k,1,18,sr)) > 1.0E-8 ) THEN first = .FALSE. height = zw(k) ENDIF IF ( hom(k,1,18,sr) < 0.0 .AND. & abs(hom(k,1,18,sr)) > 1.0E-8 .AND. & hom(k+1,1,18,sr) > hom(k,1,18,sr) ) THEN IF ( zw(k) < 1.5 * height ) THEN z_i(1) = zw(k) ELSE z_i(1) = height ENDIF EXIT ENDIF ENDDO ENDIF ! !-- Second scheme: Gradient scheme from Sullivan et al. (1998), modified !-- by Uhlenbrock(2006). The boundary layer height is the height with the !-- maximal local temperature gradient: starting from the second (the last !-- but one) vertical gridpoint, the local gradient must be at least !-- 0.2K/100m and greater than the next four gradients. !-- WARNING: The threshold value of 0.2K/100m must be adjusted for the !-- ocean case! z_i(2) = 0.0 DO k = nzb+1, nzt+1 dptdz(k) = ( hom(k,1,4,sr) - hom(k-1,1,4,sr) ) * ddzu(k) ENDDO dptdz_threshold = 0.2 / 100.0 IF ( ocean ) THEN DO k = nzt+1, nzb+5, -1 IF ( dptdz(k) > dptdz_threshold .AND. & dptdz(k) > dptdz(k-1) .AND. dptdz(k) > dptdz(k-2) .AND. & dptdz(k) > dptdz(k-3) .AND. dptdz(k) > dptdz(k-4) ) THEN z_i(2) = zw(k-1) EXIT ENDIF ENDDO ELSE DO k = nzb+1, nzt-3 IF ( dptdz(k) > dptdz_threshold .AND. & dptdz(k) > dptdz(k+1) .AND. dptdz(k) > dptdz(k+2) .AND. & dptdz(k) > dptdz(k+3) .AND. dptdz(k) > dptdz(k+4) ) THEN z_i(2) = zw(k-1) EXIT ENDIF ENDDO ENDIF hom(nzb+6,1,pr_palm,sr) = z_i(1) hom(nzb+7,1,pr_palm,sr) = z_i(2) ! !-- Computation of both the characteristic vertical velocity and !-- the characteristic convective boundary layer temperature. !-- The horizontal average at nzb+1 is input for the average temperature. IF ( hom(nzb,1,18,sr) > 0.0 .AND. abs(hom(nzb,1,18,sr)) > 1.0E-8 & .AND. z_i(1) /= 0.0 ) THEN hom(nzb+8,1,pr_palm,sr) = ( g / hom(nzb+1,1,4,sr) * & hom(nzb,1,18,sr) * & ABS( z_i(1) ) )**0.333333333 !-- so far this only works if Prandtl layer is used hom(nzb+11,1,pr_palm,sr) = hom(nzb,1,16,sr) / hom(nzb+8,1,pr_palm,sr) ELSE hom(nzb+8,1,pr_palm,sr) = 0.0 hom(nzb+11,1,pr_palm,sr) = 0.0 ENDIF ! !-- Collect the time series quantities ts_value(1,sr) = hom(nzb+4,1,pr_palm,sr) ! E ts_value(2,sr) = hom(nzb+5,1,pr_palm,sr) ! E* ts_value(3,sr) = dt_3d ts_value(4,sr) = hom(nzb,1,pr_palm,sr) ! u* ts_value(5,sr) = hom(nzb+3,1,pr_palm,sr) ! th* ts_value(6,sr) = u_max ts_value(7,sr) = v_max ts_value(8,sr) = w_max ts_value(9,sr) = hom(nzb+10,1,pr_palm,sr) ! new divergence ts_value(10,sr) = hom(nzb+9,1,pr_palm,sr) ! old Divergence ts_value(11,sr) = hom(nzb+6,1,pr_palm,sr) ! z_i(1) ts_value(12,sr) = hom(nzb+7,1,pr_palm,sr) ! z_i(2) ts_value(13,sr) = hom(nzb+8,1,pr_palm,sr) ! w* ts_value(14,sr) = hom(nzb,1,16,sr) ! w'pt' at k=0 ts_value(15,sr) = hom(nzb+1,1,16,sr) ! w'pt' at k=1 ts_value(16,sr) = hom(nzb+1,1,18,sr) ! wpt at k=1 ts_value(17,sr) = hom(nzb,1,4,sr) ! pt(0) ts_value(18,sr) = hom(nzb+1,1,4,sr) ! pt(zp) ts_value(19,sr) = hom(nzb+1,1,pr_palm,sr) ! u'w' at k=0 ts_value(20,sr) = hom(nzb+2,1,pr_palm,sr) ! v'w' at k=0 ts_value(21,sr) = hom(nzb,1,48,sr) ! w"q" at k=0 IF ( ts_value(5,sr) /= 0.0 ) THEN ts_value(22,sr) = ts_value(4,sr)**2 / & ( kappa * g * ts_value(5,sr) / ts_value(18,sr) ) ! L ELSE ts_value(22,sr) = 10000.0 ENDIF ts_value(23,sr) = hom(nzb+12,1,pr_palm,sr) ! q* ! !-- Calculate additional statistics provided by the user interface CALL user_statistics( 'time_series', sr, 0 ) ENDDO ! loop of the subregions ! !-- If required, sum up horizontal averages for subsequent time averaging IF ( do_sum ) THEN IF ( average_count_pr == 0 ) hom_sum = 0.0 hom_sum = hom_sum + hom(:,1,:,:) average_count_pr = average_count_pr + 1 do_sum = .FALSE. ENDIF ! !-- Set flag for other UPs (e.g. output routines, but also buoyancy). !-- This flag is reset after each time step in time_integration. flow_statistics_called = .TRUE. CALL cpu_log( log_point(10), 'flow_statistics', 'stop' ) END SUBROUTINE flow_statistics