MODULE prognostic_equations_mod !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! PALM is free software: you can redistribute it and/or modify it under the terms ! of the GNU General Public License as published by the Free Software Foundation, ! either version 3 of the License, or (at your option) any later version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ------------------ ! ! ! Former revisions: ! ----------------- ! $Id: prognostic_equations.f90 1182 2013-06-14 09:07:24Z raasch $ ! ! 1179 2013-06-14 05:57:58Z raasch ! two arguments removed from routine buoyancy, ref_state updated on device ! ! 1128 2013-04-12 06:19:32Z raasch ! those parts requiring global communication moved to time_integration, ! loop index bounds in accelerator version replaced by i_left, i_right, j_south, ! j_north ! ! 1115 2013-03-26 18:16:16Z hoffmann ! optimized cloud physics: calculation of microphysical tendencies transfered ! to microphysics.f90; qr and nr are only calculated if precipitation is required ! ! 1111 2013-03-08 23:54:10Z raasch ! update directives for prognostic quantities removed ! ! 1106 2013-03-04 05:31:38Z raasch ! small changes in code formatting ! ! 1092 2013-02-02 11:24:22Z raasch ! unused variables removed ! ! 1053 2012-11-13 17:11:03Z hoffmann ! implementation of two new prognostic equations for rain drop concentration (nr) ! and rain water content (qr) ! ! currently, only available for cache loop optimization ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1019 2012-09-28 06:46:45Z raasch ! non-optimized version of prognostic_equations removed ! ! 1015 2012-09-27 09:23:24Z raasch ! new branch prognostic_equations_acc ! OpenACC statements added + code changes required for GPU optimization ! ! 1001 2012-09-13 14:08:46Z raasch ! all actions concerning leapfrog- and upstream-spline-scheme removed ! ! 978 2012-08-09 08:28:32Z fricke ! km_damp_x and km_damp_y removed in calls of diffusion_u and diffusion_v ! add ptdf_x, ptdf_y for damping the potential temperature at the inflow ! boundary in case of non-cyclic lateral boundaries ! Bugfix: first thread index changes for WS-scheme at the inflow ! ! 940 2012-07-09 14:31:00Z raasch ! temperature equation can be switched off ! ! 785 2011-11-28 09:47:19Z raasch ! new factor rdf_sc allows separate Rayleigh damping of scalars ! ! 736 2011-08-17 14:13:26Z suehring ! Bugfix: determination of first thread index i for WS-scheme ! ! 709 2011-03-30 09:31:40Z raasch ! formatting adjustments ! ! 673 2011-01-18 16:19:48Z suehring ! Consideration of the pressure gradient (steered by tsc(4)) during the time ! integration removed. ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! Calls of the advection routines with WS5 added. ! Calls of ws_statistics added to set the statistical arrays to zero after each ! time step. ! ! 531 2010-04-21 06:47:21Z heinze ! add call of subsidence in the equation for humidity / passive scalar ! ! 411 2009-12-11 14:15:58Z heinze ! add call of subsidence in the equation for potential temperature ! ! 388 2009-09-23 09:40:33Z raasch ! prho is used instead of rho in diffusion_e, ! external pressure gradient ! ! 153 2008-03-19 09:41:30Z steinfeld ! add call of plant_canopy_model in the prognostic equation for ! the potential temperature and for the passive scalar ! ! 138 2007-11-28 10:03:58Z letzel ! add call of subroutines that evaluate the canopy drag terms, ! add wall_*flux to parameter list of calls of diffusion_s ! ! 106 2007-08-16 14:30:26Z raasch ! +uswst, vswst as arguments in calls of diffusion_u|v, ! loops for u and v are starting from index nxlu, nysv, respectively (needed ! for non-cyclic boundary conditions) ! ! 97 2007-06-21 08:23:15Z raasch ! prognostic equation for salinity, density is calculated from equation of ! state for seawater and is used for calculation of buoyancy, ! +eqn_state_seawater_mod ! diffusion_e is called with argument rho in case of ocean runs, ! new argument zw in calls of diffusion_e, new argument pt_/prho_reference ! in calls of buoyancy and diffusion_e, calc_mean_pt_profile renamed ! calc_mean_profile ! ! 75 2007-03-22 09:54:05Z raasch ! checking for negative q and limiting for positive values, ! z0 removed from arguments in calls of diffusion_u/v/w, uxrp, vynp eliminated, ! subroutine names changed to .._noopt, .._cache, and .._vector, ! moisture renamed humidity, Bott-Chlond-scheme can be used in the ! _vector-version ! ! 19 2007-02-23 04:53:48Z raasch ! Calculation of e, q, and pt extended for gridpoint nzt, ! handling of given temperature/humidity/scalar fluxes at top surface ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.21 2006/08/04 15:01:07 raasch ! upstream scheme can be forced to be used for tke (use_upstream_for_tke) ! regardless of the timestep scheme used for the other quantities, ! new argument diss in call of diffusion_e ! ! Revision 1.1 2000/04/13 14:56:27 schroeter ! Initial revision ! ! ! Description: ! ------------ ! Solving the prognostic equations. !------------------------------------------------------------------------------! USE arrays_3d USE control_parameters USE cpulog USE eqn_state_seawater_mod USE grid_variables USE indices USE interfaces USE pegrid USE pointer_interfaces USE statistics USE advec_ws USE advec_s_pw_mod USE advec_s_up_mod USE advec_u_pw_mod USE advec_u_up_mod USE advec_v_pw_mod USE advec_v_up_mod USE advec_w_pw_mod USE advec_w_up_mod USE buoyancy_mod USE calc_precipitation_mod USE calc_radiation_mod USE coriolis_mod USE diffusion_e_mod USE diffusion_s_mod USE diffusion_u_mod USE diffusion_v_mod USE diffusion_w_mod USE impact_of_latent_heat_mod USE microphysics_mod USE plant_canopy_model_mod USE production_e_mod USE subsidence_mod USE user_actions_mod PRIVATE PUBLIC prognostic_equations_cache, prognostic_equations_vector, & prognostic_equations_acc INTERFACE prognostic_equations_cache MODULE PROCEDURE prognostic_equations_cache END INTERFACE prognostic_equations_cache INTERFACE prognostic_equations_vector MODULE PROCEDURE prognostic_equations_vector END INTERFACE prognostic_equations_vector INTERFACE prognostic_equations_acc MODULE PROCEDURE prognostic_equations_acc END INTERFACE prognostic_equations_acc CONTAINS SUBROUTINE prognostic_equations_cache !------------------------------------------------------------------------------! ! Version with one optimized loop over all equations. It is only allowed to ! be called for the Wicker and Skamarock or Piascek-Williams advection scheme. ! ! Here the calls of most subroutines are embedded in two DO loops over i and j, ! so communication between CPUs is not allowed (does not make sense) within ! these loops. ! ! (Optimized to avoid cache missings, i.e. for Power4/5-architectures.) !------------------------------------------------------------------------------! IMPLICIT NONE INTEGER :: i, i_omp_start, j, k, omp_get_thread_num, tn = 0 LOGICAL :: loop_start ! !-- Time measurement can only be performed for the whole set of equations CALL cpu_log( log_point(32), 'all progn.equations', 'start' ) ! !-- Loop over all prognostic equations !$OMP PARALLEL private (i,i_omp_start,j,k,loop_start,tn) !$ tn = omp_get_thread_num() loop_start = .TRUE. !$OMP DO DO i = nxl, nxr ! !-- Store the first loop index. It differs for each thread and is required !-- later in advec_ws IF ( loop_start ) THEN loop_start = .FALSE. i_omp_start = i ENDIF DO j = nys, nyn ! !-- Tendency terms for u-velocity component IF ( .NOT. outflow_l .OR. i > nxl ) THEN tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN IF ( ( inflow_l .OR. outflow_l ) .AND. i_omp_start == nxl ) THEN CALL advec_u_ws( i, j, i_omp_start + 1, tn ) ELSE CALL advec_u_ws( i, j, i_omp_start, tn ) ENDIF ELSE CALL advec_u_pw( i, j ) ENDIF ELSE CALL advec_u_up( i, j ) ENDIF CALL diffusion_u( i, j ) CALL coriolis( i, j, 1 ) IF ( sloping_surface .AND. .NOT. neutral ) THEN CALL buoyancy( i, j, pt, 1 ) ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( i, j, 1 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) ENDDO ENDIF CALL user_actions( i, j, 'u-tendency' ) ! !-- Prognostic equation for u-velocity component DO k = nzb_u_inner(j,i)+1, nzt u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tu_m(k,j,i) ) & - tsc(5) * rdf(k) * ( u(k,j,i) - ug(k) ) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_u_inner(j,i)+1, nzt tu_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_u_inner(j,i)+1, nzt tu_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tu_m(k,j,i) ENDDO ENDIF ENDIF ENDIF ! !-- Tendency terms for v-velocity component IF ( .NOT. outflow_s .OR. j > nys ) THEN tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_v_ws( i, j, i_omp_start, tn ) ELSE CALL advec_v_pw( i, j ) ENDIF ELSE CALL advec_v_up( i, j ) ENDIF CALL diffusion_v( i, j ) CALL coriolis( i, j, 2 ) ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( i, j, 2 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) ENDDO ENDIF CALL user_actions( i, j, 'v-tendency' ) ! !-- Prognostic equation for v-velocity component DO k = nzb_v_inner(j,i)+1, nzt v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tv_m(k,j,i) ) & - tsc(5) * rdf(k) * ( v(k,j,i) - vg(k) ) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_v_inner(j,i)+1, nzt tv_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_v_inner(j,i)+1, nzt tv_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tv_m(k,j,i) ENDDO ENDIF ENDIF ENDIF ! !-- Tendency terms for w-velocity component tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_w_ws( i, j, i_omp_start, tn ) ELSE CALL advec_w_pw( i, j ) END IF ELSE CALL advec_w_up( i, j ) ENDIF CALL diffusion_w( i, j ) CALL coriolis( i, j, 3 ) IF ( .NOT. neutral ) THEN IF ( ocean ) THEN CALL buoyancy( i, j, rho, 3 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy( i, j, pt, 3 ) ELSE CALL buoyancy( i, j, vpt, 3 ) ENDIF ENDIF ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( i, j, 3 ) CALL user_actions( i, j, 'w-tendency' ) ! !-- Prognostic equation for w-velocity component DO k = nzb_w_inner(j,i)+1, nzt-1 w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tw_m(k,j,i) ) & - tsc(5) * rdf(k) * w(k,j,i) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_w_inner(j,i)+1, nzt-1 tw_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_w_inner(j,i)+1, nzt-1 tw_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tw_m(k,j,i) ENDDO ENDIF ENDIF ! !-- If required, calculate tendencies for total water content, liquid water !-- potential temperature, rain water content and rain drop concentration IF ( cloud_physics .AND. icloud_scheme == 0 ) CALL microphysics_control( i, j ) ! !-- If required, compute prognostic equation for potential temperature IF ( .NOT. neutral ) THEN ! !-- Tendency terms for potential temperature tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, pt, 'pt', flux_s_pt, diss_s_pt, & flux_l_pt, diss_l_pt, i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, pt ) ENDIF ELSE CALL advec_s_up( i, j, pt ) ENDIF CALL diffusion_s( i, j, pt, shf, tswst, wall_heatflux ) ! !-- If required compute heating/cooling due to long wave radiation !-- processes IF ( radiation ) THEN CALL calc_radiation( i, j ) ENDIF ! !-- Using microphysical tendencies (latent heat) IF ( cloud_physics ) THEN IF ( icloud_scheme == 0 ) THEN tend(:,j,i) = tend(:,j,i) + tend_pt(:,j,i) ELSEIF ( icloud_scheme == 1 .AND. precipitation ) THEN CALL impact_of_latent_heat( i, j ) ENDIF ENDIF ! !-- Consideration of heat sources within the plant canopy IF ( plant_canopy .AND. cthf /= 0.0 ) THEN CALL plant_canopy_model( i, j, 4 ) ENDIF ! !-- If required, compute effect of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( i, j, tend, pt, pt_init ) ENDIF CALL user_actions( i, j, 'pt-tendency' ) ! !-- Prognostic equation for potential temperature DO k = nzb_s_inner(j,i)+1, nzt pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tpt_m(k,j,i) ) & - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt tpt_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt tpt_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tpt_m(k,j,i) ENDDO ENDIF ENDIF ENDIF ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN ! !-- Tendency-terms for salinity tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) & THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, sa, 'sa', flux_s_sa, & diss_s_sa, flux_l_sa, diss_l_sa, i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, sa ) ENDIF ELSE CALL advec_s_up( i, j, sa ) ENDIF CALL diffusion_s( i, j, sa, saswsb, saswst, wall_salinityflux ) CALL user_actions( i, j, 'sa-tendency' ) ! !-- Prognostic equation for salinity DO k = nzb_s_inner(j,i)+1, nzt sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tsa_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( sa(k,j,i) - sa_init(k) ) IF ( sa_p(k,j,i) < 0.0 ) sa_p(k,j,i) = 0.1 * sa(k,j,i) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt tsa_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt tsa_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tsa_m(k,j,i) ENDDO ENDIF ENDIF ! !-- Calculate density by the equation of state for seawater CALL eqn_state_seawater( i, j ) ENDIF ! !-- If required, compute prognostic equation for total water content / !-- scalar IF ( humidity .OR. passive_scalar ) THEN ! !-- Tendency-terms for total water content / scalar tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) & THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, q, 'q', flux_s_q, & diss_s_q, flux_l_q, diss_l_q, i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, q ) ENDIF ELSE CALL advec_s_up( i, j, q ) ENDIF CALL diffusion_s( i, j, q, qsws, qswst, wall_qflux ) ! !-- Using microphysical tendencies IF ( cloud_physics ) THEN IF ( icloud_scheme == 0 ) THEN tend(:,j,i) = tend(:,j,i) + tend_q(:,j,i) ELSEIF ( icloud_scheme == 1 .AND. precipitation ) THEN CALL calc_precipitation( i, j ) ENDIF ENDIF ! !-- Sink or source of scalar concentration due to canopy elements IF ( plant_canopy ) CALL plant_canopy_model( i, j, 5 ) ! !-- If required compute influence of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( i, j, tend, q, q_init ) ENDIF CALL user_actions( i, j, 'q-tendency' ) ! !-- Prognostic equation for total water content / scalar DO k = nzb_s_inner(j,i)+1, nzt q_p(k,j,i) = q(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tq_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( q(k,j,i) - q_init(k) ) IF ( q_p(k,j,i) < 0.0 ) q_p(k,j,i) = 0.1 * q(k,j,i) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt tq_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt tq_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tq_m(k,j,i) ENDDO ENDIF ENDIF ! !-- If required, calculate prognostic equations for rain water content !-- and rain drop concentration IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & precipitation ) THEN ! !-- Calculate prognostic equation for rain water content tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) & THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, qr, 'qr', flux_s_qr, & diss_s_qr, flux_l_qr, diss_l_qr, & i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, qr ) ENDIF ELSE CALL advec_s_up( i, j, qr ) ENDIF CALL diffusion_s( i, j, qr, qrsws, qrswst, wall_qrflux ) ! !-- Using microphysical tendencies (autoconversion, accretion, !-- evaporation; if required: sedimentation) tend(:,j,i) = tend(:,j,i) + tend_qr(:,j,i) CALL user_actions( i, j, 'qr-tendency' ) ! !-- Prognostic equation for rain water content DO k = nzb_s_inner(j,i)+1, nzt qr_p(k,j,i) = qr(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tqr_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * qr(k,j,i) IF ( qr_p(k,j,i) < 0.0 ) qr_p(k,j,i) = 0.0 ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt tqr_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt tqr_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tqr_m(k,j,i) ENDDO ENDIF ENDIF ! !-- Calculate prognostic equation for rain drop concentration. tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, nr, 'nr', flux_s_nr, & diss_s_nr, flux_l_nr, diss_l_nr, & i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, nr ) ENDIF ELSE CALL advec_s_up( i, j, nr ) ENDIF CALL diffusion_s( i, j, nr, nrsws, nrswst, wall_nrflux ) ! !-- Using microphysical tendencies (autoconversion, accretion, !-- selfcollection, breakup, evaporation; !-- if required: sedimentation) tend(:,j,i) = tend(:,j,i) + tend_nr(:,j,i) CALL user_actions( i, j, 'nr-tendency' ) ! !-- Prognostic equation for rain drop concentration DO k = nzb_s_inner(j,i)+1, nzt nr_p(k,j,i) = nr(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tnr_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * nr(k,j,i) IF ( nr_p(k,j,i) < 0.0 ) nr_p(k,j,i) = 0.0 ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt tnr_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt tnr_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tnr_m(k,j,i) ENDDO ENDIF ENDIF ENDIF ENDIF ! !-- If required, compute prognostic equation for turbulent kinetic !-- energy (TKE) IF ( .NOT. constant_diffusion ) THEN ! !-- Tendency-terms for TKE tend(:,j,i) = 0.0 IF ( timestep_scheme(1:5) == 'runge' & .AND. .NOT. use_upstream_for_tke ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, e, 'e', flux_s_e, diss_s_e, & flux_l_e, diss_l_e , i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, e ) ENDIF ELSE CALL advec_s_up( i, j, e ) ENDIF IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( i, j, prho, prho_reference ) ELSE CALL diffusion_e( i, j, pt, pt_reference ) ENDIF ELSE CALL diffusion_e( i, j, vpt, pt_reference ) ENDIF CALL production_e( i, j ) ! !-- Additional sink term for flows through plant canopies IF ( plant_canopy ) CALL plant_canopy_model( i, j, 6 ) CALL user_actions( i, j, 'e-tendency' ) ! !-- Prognostic equation for TKE. !-- Eliminate negative TKE values, which can occur due to numerical !-- reasons in the course of the integration. In such cases the old !-- TKE value is reduced by 90%. DO k = nzb_s_inner(j,i)+1, nzt e_p(k,j,i) = e(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * te_m(k,j,i) ) IF ( e_p(k,j,i) < 0.0 ) e_p(k,j,i) = 0.1 * e(k,j,i) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb_s_inner(j,i)+1, nzt te_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb_s_inner(j,i)+1, nzt te_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * te_m(k,j,i) ENDDO ENDIF ENDIF ENDIF ! TKE equation ENDDO ENDDO !$OMP END PARALLEL CALL cpu_log( log_point(32), 'all progn.equations', 'stop' ) END SUBROUTINE prognostic_equations_cache SUBROUTINE prognostic_equations_vector !------------------------------------------------------------------------------! ! Version for vector machines !------------------------------------------------------------------------------! IMPLICIT NONE INTEGER :: i, j, k REAL :: sbt ! !-- u-velocity component CALL cpu_log( log_point(5), 'u-equation', 'start' ) tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_u_ws ELSE CALL advec_u_pw ENDIF ELSE CALL advec_u_up ENDIF CALL diffusion_u CALL coriolis( 1 ) IF ( sloping_surface .AND. .NOT. neutral ) THEN CALL buoyancy( pt, 1 ) ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 1 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO i = nxlu, nxr DO j = nys, nyn DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) ENDDO ENDDO ENDDO ENDIF CALL user_actions( 'u-tendency' ) ! !-- Prognostic equation for u-velocity component DO i = nxlu, nxr DO j = nys, nyn DO k = nzb_u_inner(j,i)+1, nzt u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tu_m(k,j,i) ) & - tsc(5) * rdf(k) * ( u(k,j,i) - ug(k) ) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxlu, nxr DO j = nys, nyn DO k = nzb_u_inner(j,i)+1, nzt tu_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxlu, nxr DO j = nys, nyn DO k = nzb_u_inner(j,i)+1, nzt tu_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tu_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(5), 'u-equation', 'stop' ) ! !-- v-velocity component CALL cpu_log( log_point(6), 'v-equation', 'start' ) tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_v_ws ELSE CALL advec_v_pw END IF ELSE CALL advec_v_up ENDIF CALL diffusion_v CALL coriolis( 2 ) ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 2 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO i = nxl, nxr DO j = nysv, nyn DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) ENDDO ENDDO ENDDO ENDIF CALL user_actions( 'v-tendency' ) ! !-- Prognostic equation for v-velocity component DO i = nxl, nxr DO j = nysv, nyn DO k = nzb_v_inner(j,i)+1, nzt v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tv_m(k,j,i) ) & - tsc(5) * rdf(k) * ( v(k,j,i) - vg(k) ) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nysv, nyn DO k = nzb_v_inner(j,i)+1, nzt tv_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nysv, nyn DO k = nzb_v_inner(j,i)+1, nzt tv_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tv_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(6), 'v-equation', 'stop' ) ! !-- w-velocity component CALL cpu_log( log_point(7), 'w-equation', 'start' ) tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_w_ws ELSE CALL advec_w_pw ENDIF ELSE CALL advec_w_up ENDIF CALL diffusion_w CALL coriolis( 3 ) IF ( .NOT. neutral ) THEN IF ( ocean ) THEN CALL buoyancy( rho, 3 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy( pt, 3 ) ELSE CALL buoyancy( vpt, 3 ) ENDIF ENDIF ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 3 ) CALL user_actions( 'w-tendency' ) ! !-- Prognostic equation for w-velocity component DO i = nxl, nxr DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt-1 w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tw_m(k,j,i) ) & - tsc(5) * rdf(k) * w(k,j,i) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt-1 tw_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt-1 tw_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tw_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(7), 'w-equation', 'stop' ) ! !-- If required, compute prognostic equation for potential temperature IF ( .NOT. neutral ) THEN CALL cpu_log( log_point(13), 'pt-equation', 'start' ) ! !-- pt-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( pt, 'pt' ) ENDIF ! !-- pt-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( pt, 'pt' ) ELSE CALL advec_s_pw( pt ) ENDIF ELSE CALL advec_s_up( pt ) ENDIF ENDIF CALL diffusion_s( pt, shf, tswst, wall_heatflux ) ! !-- If required compute heating/cooling due to long wave radiation processes IF ( radiation ) THEN CALL calc_radiation ENDIF ! !-- If required compute impact of latent heat due to precipitation IF ( precipitation ) THEN CALL impact_of_latent_heat ENDIF ! !-- Consideration of heat sources within the plant canopy IF ( plant_canopy .AND. ( cthf /= 0.0 ) ) THEN CALL plant_canopy_model( 4 ) ENDIF ! !-- If required compute influence of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( tend, pt, pt_init ) ENDIF CALL user_actions( 'pt-tendency' ) ! !-- Prognostic equation for potential temperature DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tpt_m(k,j,i) ) & - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tpt_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tpt_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tpt_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(13), 'pt-equation', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN CALL cpu_log( log_point(37), 'sa-equation', 'start' ) ! !-- sa-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( sa, 'sa' ) ENDIF ! !-- sa-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( sa, 'sa' ) ELSE CALL advec_s_pw( sa ) ENDIF ELSE CALL advec_s_up( sa ) ENDIF ENDIF CALL diffusion_s( sa, saswsb, saswst, wall_salinityflux ) CALL user_actions( 'sa-tendency' ) ! !-- Prognostic equation for salinity DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tsa_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( sa(k,j,i) - sa_init(k) ) IF ( sa_p(k,j,i) < 0.0 ) sa_p(k,j,i) = 0.1 * sa(k,j,i) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tsa_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tsa_m(k,j,i) = -9.5625 * tend(k,j,i) + & 5.3125 * tsa_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(37), 'sa-equation', 'stop' ) ! !-- Calculate density by the equation of state for seawater CALL cpu_log( log_point(38), 'eqns-seawater', 'start' ) CALL eqn_state_seawater CALL cpu_log( log_point(38), 'eqns-seawater', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for total water content / scalar IF ( humidity .OR. passive_scalar ) THEN CALL cpu_log( log_point(29), 'q/s-equation', 'start' ) ! !-- Scalar/q-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( q, 'q' ) ENDIF ! !-- Scalar/q-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( q, 'q' ) ELSE CALL advec_s_pw( q ) ENDIF ELSE CALL advec_s_up( q ) ENDIF ENDIF CALL diffusion_s( q, qsws, qswst, wall_qflux ) ! !-- If required compute decrease of total water content due to !-- precipitation IF ( precipitation ) THEN CALL calc_precipitation ENDIF ! !-- Sink or source of scalar concentration due to canopy elements IF ( plant_canopy ) CALL plant_canopy_model( 5 ) ! !-- If required compute influence of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( tend, q, q_init ) ENDIF CALL user_actions( 'q-tendency' ) ! !-- Prognostic equation for total water content / scalar DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt q_p(k,j,i) = q(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tq_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( q(k,j,i) - q_init(k) ) IF ( q_p(k,j,i) < 0.0 ) q_p(k,j,i) = 0.1 * q(k,j,i) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tq_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt tq_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tq_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(29), 'q/s-equation', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for turbulent kinetic !-- energy (TKE) IF ( .NOT. constant_diffusion ) THEN CALL cpu_log( log_point(16), 'tke-equation', 'start' ) sbt = tsc(2) IF ( .NOT. use_upstream_for_tke ) THEN IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( e, 'e' ) ENDIF ENDIF ! !-- TKE-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' .OR. use_upstream_for_tke ) THEN IF ( use_upstream_for_tke ) THEN tend = 0.0 CALL advec_s_up( e ) ELSE tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( e, 'e' ) ELSE CALL advec_s_pw( e ) ENDIF ELSE CALL advec_s_up( e ) ENDIF ENDIF ENDIF IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( prho, prho_reference ) ELSE CALL diffusion_e( pt, pt_reference ) ENDIF ELSE CALL diffusion_e( vpt, pt_reference ) ENDIF CALL production_e ! !-- Additional sink term for flows through plant canopies IF ( plant_canopy ) CALL plant_canopy_model( 6 ) CALL user_actions( 'e-tendency' ) ! !-- Prognostic equation for TKE. !-- Eliminate negative TKE values, which can occur due to numerical !-- reasons in the course of the integration. In such cases the old TKE !-- value is reduced by 90%. DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt e_p(k,j,i) = e(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * te_m(k,j,i) ) IF ( e_p(k,j,i) < 0.0 ) e_p(k,j,i) = 0.1 * e(k,j,i) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt te_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_s_inner(j,i)+1, nzt te_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * te_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF CALL cpu_log( log_point(16), 'tke-equation', 'stop' ) ENDIF END SUBROUTINE prognostic_equations_vector SUBROUTINE prognostic_equations_acc !------------------------------------------------------------------------------! ! Version for accelerator boards !------------------------------------------------------------------------------! IMPLICIT NONE INTEGER :: i, j, k, runge_step REAL :: sbt ! !-- Set switch for intermediate Runge-Kutta step runge_step = 0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN runge_step = 1 ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN runge_step = 2 ENDIF ENDIF ! !-- u-velocity component !++ Statistics still not ported to accelerators !$acc update device( hom, ref_state ) CALL cpu_log( log_point(5), 'u-equation', 'start' ) IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_u_ws_acc ELSE tend = 0.0 ! to be removed later?? CALL advec_u_pw ENDIF ELSE CALL advec_u_up ENDIF CALL diffusion_u_acc CALL coriolis_acc( 1 ) IF ( sloping_surface .AND. .NOT. neutral ) THEN CALL buoyancy( pt, 1 ) ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 1 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO i = i_left, i_right DO j = j_south, j_north DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(1) * dp_smooth_factor(k) ENDDO ENDDO ENDDO ENDIF CALL user_actions( 'u-tendency' ) ! !-- Prognostic equation for u-velocity component !$acc kernels present( nzb_u_inner, rdf, tend, tu_m, u, ug, u_p ) !$acc loop DO i = i_left, i_right DO j = j_south, j_north !$acc loop vector( 32 ) DO k = 1, nzt IF ( k > nzb_u_inner(j,i) ) THEN u_p(k,j,i) = u(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tu_m(k,j,i) ) & - tsc(5) * rdf(k) * ( u(k,j,i) - ug(k) ) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tu_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tu_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tu_m(k,j,i) ENDIF ENDIF ENDDO ENDDO ENDDO !$acc end kernels CALL cpu_log( log_point(5), 'u-equation', 'stop' ) ! !-- v-velocity component CALL cpu_log( log_point(6), 'v-equation', 'start' ) IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_v_ws_acc ELSE tend = 0.0 ! to be removed later?? CALL advec_v_pw END IF ELSE CALL advec_v_up ENDIF CALL diffusion_v_acc CALL coriolis_acc( 2 ) ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 2 ) ! !-- External pressure gradient IF ( dp_external ) THEN DO i = i_left, i_right DO j = j_south, j_north DO k = dp_level_ind_b+1, nzt tend(k,j,i) = tend(k,j,i) - dpdxy(2) * dp_smooth_factor(k) ENDDO ENDDO ENDDO ENDIF CALL user_actions( 'v-tendency' ) ! !-- Prognostic equation for v-velocity component !$acc kernels present( nzb_v_inner, rdf, tend, tv_m, v, vg, v_p ) !$acc loop DO i = i_left, i_right DO j = j_south, j_north !$acc loop vector( 32 ) DO k = 1, nzt IF ( k > nzb_v_inner(j,i) ) THEN v_p(k,j,i) = v(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tv_m(k,j,i) ) & - tsc(5) * rdf(k) * ( v(k,j,i) - vg(k) ) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tv_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tv_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tv_m(k,j,i) ENDIF ENDIF ENDDO ENDDO ENDDO !$acc end kernels CALL cpu_log( log_point(6), 'v-equation', 'stop' ) ! !-- w-velocity component CALL cpu_log( log_point(7), 'w-equation', 'start' ) IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_mom ) THEN CALL advec_w_ws_acc ELSE tend = 0.0 ! to be removed later?? CALL advec_w_pw ENDIF ELSE CALL advec_w_up ENDIF CALL diffusion_w_acc CALL coriolis_acc( 3 ) IF ( .NOT. neutral ) THEN IF ( ocean ) THEN CALL buoyancy( rho, 3 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy_acc( pt, 3 ) ELSE CALL buoyancy( vpt, 3 ) ENDIF ENDIF ENDIF ! !-- Drag by plant canopy IF ( plant_canopy ) CALL plant_canopy_model( 3 ) CALL user_actions( 'w-tendency' ) ! !-- Prognostic equation for w-velocity component !$acc kernels present( nzb_w_inner, rdf, tend, tw_m, w, w_p ) !$acc loop DO i = i_left, i_right DO j = j_south, j_north !$acc loop vector( 32 ) DO k = 1, nzt-1 IF ( k > nzb_w_inner(j,i) ) THEN w_p(k,j,i) = w(k,j,i) + dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * tw_m(k,j,i) ) & - tsc(5) * rdf(k) * w(k,j,i) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tw_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tw_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tw_m(k,j,i) ENDIF ENDIF ENDDO ENDDO ENDDO !$acc end kernels CALL cpu_log( log_point(7), 'w-equation', 'stop' ) ! !-- If required, compute prognostic equation for potential temperature IF ( .NOT. neutral ) THEN CALL cpu_log( log_point(13), 'pt-equation', 'start' ) ! !-- pt-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( pt, 'pt' ) ENDIF ! !-- pt-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws_acc( pt, 'pt' ) ELSE tend = 0.0 ! to be removed later?? CALL advec_s_pw( pt ) ENDIF ELSE CALL advec_s_up( pt ) ENDIF ENDIF CALL diffusion_s_acc( pt, shf, tswst, wall_heatflux ) ! !-- If required compute heating/cooling due to long wave radiation processes IF ( radiation ) THEN CALL calc_radiation ENDIF ! !-- If required compute impact of latent heat due to precipitation IF ( precipitation ) THEN CALL impact_of_latent_heat ENDIF ! !-- Consideration of heat sources within the plant canopy IF ( plant_canopy .AND. ( cthf /= 0.0 ) ) THEN CALL plant_canopy_model( 4 ) ENDIF ! !-- If required compute influence of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( tend, pt, pt_init ) ENDIF CALL user_actions( 'pt-tendency' ) ! !-- Prognostic equation for potential temperature !$acc kernels present( nzb_s_inner, rdf_sc, ptdf_x, ptdf_y, pt_init ) & !$acc present( tend, tpt_m, pt, pt_p ) !$acc loop DO i = i_left, i_right DO j = j_south, j_north !$acc loop vector( 32 ) DO k = 1, nzt IF ( k > nzb_s_inner(j,i) ) THEN pt_p(k,j,i) = pt(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tpt_m(k,j,i) ) & - tsc(5) * ( pt(k,j,i) - pt_init(k) ) *& ( rdf_sc(k) + ptdf_x(i) + ptdf_y(j) ) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tpt_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tpt_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tpt_m(k,j,i) ENDIF ENDIF ENDDO ENDDO ENDDO !$acc end kernels CALL cpu_log( log_point(13), 'pt-equation', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN CALL cpu_log( log_point(37), 'sa-equation', 'start' ) ! !-- sa-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( sa, 'sa' ) ENDIF ! !-- sa-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( sa, 'sa' ) ELSE CALL advec_s_pw( sa ) ENDIF ELSE CALL advec_s_up( sa ) ENDIF ENDIF CALL diffusion_s( sa, saswsb, saswst, wall_salinityflux ) CALL user_actions( 'sa-tendency' ) ! !-- Prognostic equation for salinity DO i = i_left, i_right DO j = j_south, j_north DO k = nzb_s_inner(j,i)+1, nzt sa_p(k,j,i) = sa(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tsa_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( sa(k,j,i) - sa_init(k) ) IF ( sa_p(k,j,i) < 0.0 ) sa_p(k,j,i) = 0.1 * sa(k,j,i) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tsa_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tsa_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tsa_m(k,j,i) ENDIF ENDDO ENDDO ENDDO CALL cpu_log( log_point(37), 'sa-equation', 'stop' ) ! !-- Calculate density by the equation of state for seawater CALL cpu_log( log_point(38), 'eqns-seawater', 'start' ) CALL eqn_state_seawater CALL cpu_log( log_point(38), 'eqns-seawater', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for total water content / scalar IF ( humidity .OR. passive_scalar ) THEN CALL cpu_log( log_point(29), 'q/s-equation', 'start' ) ! !-- Scalar/q-tendency terms with communication sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( q, 'q' ) ENDIF ! !-- Scalar/q-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' ) THEN tend = 0.0 IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( q, 'q' ) ELSE CALL advec_s_pw( q ) ENDIF ELSE CALL advec_s_up( q ) ENDIF ENDIF CALL diffusion_s( q, qsws, qswst, wall_qflux ) ! !-- If required compute decrease of total water content due to !-- precipitation IF ( precipitation ) THEN CALL calc_precipitation ENDIF ! !-- Sink or source of scalar concentration due to canopy elements IF ( plant_canopy ) CALL plant_canopy_model( 5 ) ! !-- If required compute influence of large-scale subsidence/ascent IF ( large_scale_subsidence ) THEN CALL subsidence( tend, q, q_init ) ENDIF CALL user_actions( 'q-tendency' ) ! !-- Prognostic equation for total water content / scalar DO i = i_left, i_right DO j = j_south, j_north DO k = nzb_s_inner(j,i)+1, nzt q_p(k,j,i) = q(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * tq_m(k,j,i) ) & - tsc(5) * rdf_sc(k) * & ( q(k,j,i) - q_init(k) ) IF ( q_p(k,j,i) < 0.0 ) q_p(k,j,i) = 0.1 * q(k,j,i) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN tq_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN tq_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * tq_m(k,j,i) ENDIF ENDDO ENDDO ENDDO CALL cpu_log( log_point(29), 'q/s-equation', 'stop' ) ENDIF ! !-- If required, compute prognostic equation for turbulent kinetic !-- energy (TKE) IF ( .NOT. constant_diffusion ) THEN CALL cpu_log( log_point(16), 'tke-equation', 'start' ) sbt = tsc(2) IF ( .NOT. use_upstream_for_tke ) THEN IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step. Thus: sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( e, 'e' ) ENDIF ENDIF ! !-- TKE-tendency terms with no communication IF ( scalar_advec /= 'bc-scheme' .OR. use_upstream_for_tke ) THEN IF ( use_upstream_for_tke ) THEN tend = 0.0 CALL advec_s_up( e ) ELSE IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws_acc( e, 'e' ) ELSE tend = 0.0 ! to be removed later?? CALL advec_s_pw( e ) ENDIF ELSE tend = 0.0 ! to be removed later?? CALL advec_s_up( e ) ENDIF ENDIF ENDIF IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( prho, prho_reference ) ELSE CALL diffusion_e_acc( pt, pt_reference ) ENDIF ELSE CALL diffusion_e( vpt, pt_reference ) ENDIF CALL production_e_acc ! !-- Additional sink term for flows through plant canopies IF ( plant_canopy ) CALL plant_canopy_model( 6 ) CALL user_actions( 'e-tendency' ) ! !-- Prognostic equation for TKE. !-- Eliminate negative TKE values, which can occur due to numerical !-- reasons in the course of the integration. In such cases the old TKE !-- value is reduced by 90%. !$acc kernels present( e, e_p, nzb_s_inner, tend, te_m ) !$acc loop DO i = i_left, i_right DO j = j_south, j_north !$acc loop vector( 32 ) DO k = 1, nzt IF ( k > nzb_s_inner(j,i) ) THEN e_p(k,j,i) = e(k,j,i) + dt_3d * ( sbt * tend(k,j,i) + & tsc(3) * te_m(k,j,i) ) IF ( e_p(k,j,i) < 0.0 ) e_p(k,j,i) = 0.1 * e(k,j,i) ! !-- Tendencies for the next Runge-Kutta step IF ( runge_step == 1 ) THEN te_m(k,j,i) = tend(k,j,i) ELSEIF ( runge_step == 2 ) THEN te_m(k,j,i) = -9.5625 * tend(k,j,i) + 5.3125 * te_m(k,j,i) ENDIF ENDIF ENDDO ENDDO ENDDO !$acc end kernels CALL cpu_log( log_point(16), 'tke-equation', 'stop' ) ENDIF END SUBROUTINE prognostic_equations_acc END MODULE prognostic_equations_mod