MODULE prognostic_equations_mod !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: prognostic_equations.f90 532 2010-04-21 13:33:38Z maronga $ ! ! 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_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 plant_canopy_model_mod USE production_e_mod USE subsidence_mod USE user_actions_mod PRIVATE PUBLIC prognostic_equations_noopt, prognostic_equations_cache, & prognostic_equations_vector INTERFACE prognostic_equations_noopt MODULE PROCEDURE prognostic_equations_noopt END INTERFACE prognostic_equations_noopt 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 CONTAINS SUBROUTINE prognostic_equations_noopt !------------------------------------------------------------------------------! ! Version with single loop optimization ! ! (Optimized over each single prognostic equation.) !------------------------------------------------------------------------------! IMPLICIT NONE CHARACTER (LEN=9) :: time_to_string INTEGER :: i, j, k REAL :: sat, sbt ! !-- Calculate those variables needed in the tendency terms which need !-- global communication CALL calc_mean_profile( pt, 4 ) IF ( ocean ) CALL calc_mean_profile( rho, 64 ) IF ( humidity ) CALL calc_mean_profile( vpt, 44 ) ! !-- u-velocity component CALL cpu_log( log_point(5), 'u-equation', 'start' ) ! !-- u-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_u_ups ENDIF ! !-- u-tendency terms with no communication DO i = nxlu, nxr DO j = nys, nyn ! !-- Tendency terms IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_u_pw( i, j ) ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_u_up( i, j ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_u( i, j, ddzu, ddzw, km_m, km_damp_y, tend, u_m, & usws_m, uswst_m, v_m, w_m ) ELSE CALL diffusion_u( i, j, ddzu, ddzw, km, km_damp_y, tend, u, usws, & uswst, v, w ) ENDIF CALL coriolis( i, j, 1 ) IF ( sloping_surface ) CALL buoyancy( i, j, pt, pt_reference, 1, 4 ) ! !-- 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) = ( 1.0-tsc(1) ) * u_m(k,j,i) + tsc(1) * u(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tu_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j,i-1) ) * ddx & ) - & 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 ENDDO ENDDO CALL cpu_log( log_point(5), 'u-equation', 'stop' ) ! !-- v-velocity component CALL cpu_log( log_point(6), 'v-equation', 'start' ) ! !-- v-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_v_ups ENDIF ! !-- v-tendency terms with no communication DO i = nxl, nxr DO j = nysv, nyn ! !-- Tendency terms IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_v_pw( i, j ) ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_v_up( i, j ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_v( i, j, ddzu, ddzw, km_m, km_damp_x, tend, u_m, & v_m, vsws_m, vswst_m, w_m ) ELSE CALL diffusion_v( i, j, ddzu, ddzw, km, km_damp_x, tend, u, v, & vsws, vswst, w ) ENDIF 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) = ( 1.0-tsc(1) ) * v_m(k,j,i) + tsc(1) * v(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tv_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j-1,i) ) * ddy & ) - & 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 ENDDO ENDDO CALL cpu_log( log_point(6), 'v-equation', 'stop' ) ! !-- w-velocity component CALL cpu_log( log_point(7), 'w-equation', 'start' ) ! !-- w-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_w_ups ENDIF ! !-- w-tendency terms with no communication DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency terms IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_w_pw( i, j ) ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_w_up( i, j ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_w( i, j, ddzu, ddzw, km_m, km_damp_x, km_damp_y, & tend, u_m, v_m, w_m ) ELSE CALL diffusion_w( i, j, ddzu, ddzw, km, km_damp_x, km_damp_y, & tend, u, v, w ) ENDIF CALL coriolis( i, j, 3 ) IF ( ocean ) THEN CALL buoyancy( i, j, rho, rho_reference, 3, 64 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy( i, j, pt, pt_reference, 3, 4 ) ELSE CALL buoyancy( i, j, vpt, pt_reference, 3, 44 ) 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) = ( 1.0-tsc(1) ) * w_m(k,j,i) + tsc(1) * w(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tw_m(k,j,i) & - tsc(4) * ( p(k+1,j,i) - p(k,j,i) ) * ddzu(k+1) & ) - & 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 ENDDO ENDDO CALL cpu_log( log_point(7), 'w-equation', 'stop' ) ! !-- potential temperature CALL cpu_log( log_point(13), 'pt-equation', 'start' ) ! !-- pt-tendency terms with communication sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( pt, 'pt' ) ELSE IF ( tsc(2) /= 2.0 .AND. scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( pt, 'pt' ) ENDIF ENDIF ! !-- pt-tendency terms with no communication DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency terms IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( i, j, ddzu, ddzw, kh, pt, shf, tswst, & wall_heatflux, tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_s_pw( i, j, pt ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_s_up( i, j, pt ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) & THEN CALL diffusion_s( i, j, ddzu, ddzw, kh_m, pt_m, shf_m, & tswst_m, wall_heatflux, tend ) ELSE CALL diffusion_s( i, j, ddzu, ddzw, kh, pt, shf, tswst, & wall_heatflux, tend ) ENDIF ENDIF ! !-- If required compute heating/cooling due to long wave radiation !-- processes IF ( radiation ) THEN CALL calc_radiation( i, j ) ENDIF ! !-- If required compute impact of latent heat due to precipitation IF ( precipitation ) THEN CALL impact_of_latent_heat( i, j ) 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 influence 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) = ( 1 - sat ) * pt_m(k,j,i) + sat * pt(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tpt_m(k,j,i) & ) - & tsc(5) * rdf(k) * ( pt(k,j,i) - pt_init(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_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 ENDDO ENDDO CALL cpu_log( log_point(13), 'pt-equation', 'stop' ) ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN CALL cpu_log( log_point(37), 'sa-equation', 'start' ) ! !-- sa-tendency terms with communication sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( sa, 'sa' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( sa, 'sa' ) ENDIF ENDIF ENDIF ! !-- sa terms with no communication DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency-terms IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( i, j, ddzu, ddzw, kh, sa, saswsb, saswst, & wall_salinityflux, tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_s_pw( i, j, sa ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_s_up( i, j, sa ) ENDIF ENDIF CALL diffusion_s( i, j, ddzu, ddzw, kh, sa, saswsb, saswst, & wall_salinityflux, tend ) ENDIF 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) = sat * sa(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tsa_m(k,j,i) & ) - & tsc(5) * rdf(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 ) ENDDO ENDDO CALL cpu_log( log_point(37), 'sa-equation', '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 sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( q, 'q' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( q, 'q' ) ENDIF ENDIF ENDIF ! !-- Scalar/q-tendency terms with no communication DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency-terms IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( i, j, ddzu, ddzw, kh, q, qsws, qswst, & wall_qflux, tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend(:,j,i) = 0.0 CALL advec_s_pw( i, j, q ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_s_up( i, j, q ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' )& THEN CALL diffusion_s( i, j, ddzu, ddzw, kh_m, q_m, qsws_m, & qswst_m, wall_qflux, tend ) ELSE CALL diffusion_s( i, j, ddzu, ddzw, kh, q, qsws, qswst, & wall_qflux, tend ) ENDIF ENDIF ! !-- If required compute decrease of total water content due to !-- precipitation IF ( precipitation ) THEN CALL calc_precipitation( i, j ) 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) = ( 1 - sat ) * q_m(k,j,i) + sat * q(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tq_m(k,j,i) & ) - & tsc(5) * rdf(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 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' ) ! !-- TKE-tendency terms with communication CALL production_e_init sat = tsc(1) 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 when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( e, 'e' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( e, 'e' ) ENDIF ENDIF ENDIF ENDIF ! !-- TKE-tendency terms with no communication DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency-terms IF ( scalar_advec == 'bc-scheme' .AND. & .NOT. use_upstream_for_tke ) THEN IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, km, & l_grid, prho, prho_reference, rif, & tend, zu, zw ) ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, km, & l_grid, pt, pt_reference, rif, tend, & zu, zw ) ENDIF ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, km, & l_grid, vpt, pt_reference, rif, tend, zu, & zw ) ENDIF ELSE IF ( use_upstream_for_tke ) THEN tend(:,j,i) = 0.0 CALL advec_s_up( i, j, e ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) & THEN tend(:,j,i) = 0.0 CALL advec_s_pw( i, j, e ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend(:,j,i) = 0.0 CALL advec_s_up( i, j, e ) ENDIF ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' )& THEN IF ( .NOT. humidity ) THEN CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e_m, & km_m, l_grid, pt_m, pt_reference, & rif_m, tend, zu, zw ) ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e_m, & km_m, l_grid, vpt_m, pt_reference, & rif_m, tend, zu, zw ) ENDIF ELSE IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, & km, l_grid, prho, prho_reference, & rif, tend, zu, zw ) ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, & km, l_grid, pt, pt_reference, rif, & tend, zu, zw ) ENDIF ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, km, & l_grid, vpt, pt_reference, rif, tend, & zu, zw ) ENDIF ENDIF 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) = ( 1 - sat ) * e_m(k,j,i) + sat * 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 ! !-- 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 ENDDO ENDDO CALL cpu_log( log_point(16), 'tke-equation', 'stop' ) ENDIF END SUBROUTINE prognostic_equations_noopt SUBROUTINE prognostic_equations_cache !------------------------------------------------------------------------------! ! Version with one optimized loop over all equations. It is only allowed to ! be called for the standard 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 CHARACTER (LEN=9) :: time_to_string INTEGER :: i, j, k ! !-- Time measurement can only be performed for the whole set of equations CALL cpu_log( log_point(32), 'all progn.equations', 'start' ) ! !-- Calculate those variables needed in the tendency terms which need !-- global communication CALL calc_mean_profile( pt, 4 ) IF ( ocean ) CALL calc_mean_profile( rho, 64 ) IF ( humidity ) CALL calc_mean_profile( vpt, 44 ) IF ( .NOT. constant_diffusion ) CALL production_e_init ! !-- Loop over all prognostic equations !$OMP PARALLEL private (i,j,k) !$OMP DO DO i = nxl, nxr DO j = nys, nyn ! !-- Tendency terms for u-velocity component IF ( .NOT. outflow_l .OR. i > nxl ) THEN tend(:,j,i) = 0.0 IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN CALL advec_u_pw( i, j ) ELSE CALL advec_u_up( i, j ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) & THEN CALL diffusion_u( i, j, ddzu, ddzw, km_m, km_damp_y, tend, & u_m, usws_m, uswst_m, v_m, w_m ) ELSE CALL diffusion_u( i, j, ddzu, ddzw, km, km_damp_y, tend, u, & usws, uswst, v, w ) ENDIF CALL coriolis( i, j, 1 ) IF ( sloping_surface ) CALL buoyancy( i, j, pt, pt_reference, 1, & 4 ) ! !-- 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) = ( 1.0-tsc(1) ) * u_m(k,j,i) + tsc(1) * u(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tu_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j,i-1) ) * ddx & ) - & 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 ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN CALL advec_v_pw( i, j ) ELSE CALL advec_v_up( i, j ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) & THEN CALL diffusion_v( i, j, ddzu, ddzw, km_m, km_damp_x, tend, & u_m, v_m, vsws_m, vswst_m, w_m ) ELSE CALL diffusion_v( i, j, ddzu, ddzw, km, km_damp_x, tend, u, v, & vsws, vswst, w ) ENDIF 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) = ( 1.0-tsc(1) ) * v_m(k,j,i) + tsc(1) * v(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tv_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j-1,i) ) * ddy & ) - & 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 ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN CALL advec_w_pw( i, j ) ELSE CALL advec_w_up( i, j ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) & THEN CALL diffusion_w( i, j, ddzu, ddzw, km_m, km_damp_x, & km_damp_y, tend, u_m, v_m, w_m ) ELSE CALL diffusion_w( i, j, ddzu, ddzw, km, km_damp_x, km_damp_y, & tend, u, v, w ) ENDIF CALL coriolis( i, j, 3 ) IF ( ocean ) THEN CALL buoyancy( i, j, rho, rho_reference, 3, 64 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy( i, j, pt, pt_reference, 3, 4 ) ELSE CALL buoyancy( i, j, vpt, pt_reference, 3, 44 ) 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) = ( 1.0-tsc(1) ) * w_m(k,j,i) + tsc(1) * w(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tw_m(k,j,i) & - tsc(4) * ( p(k+1,j,i) - p(k,j,i) ) * ddzu(k+1) & ) - & 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 ! !-- Tendency terms for potential temperature tend(:,j,i) = 0.0 IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN CALL advec_s_pw( i, j, pt ) ELSE CALL advec_s_up( i, j, pt ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) & THEN CALL diffusion_s( i, j, ddzu, ddzw, kh_m, pt_m, shf_m, & tswst_m, wall_heatflux, tend ) ELSE CALL diffusion_s( i, j, ddzu, ddzw, kh, pt, shf, tswst, & wall_heatflux, tend ) ENDIF ! !-- If required compute heating/cooling due to long wave radiation !-- processes IF ( radiation ) THEN CALL calc_radiation( i, j ) ENDIF ! !-- If required compute impact of latent heat due to precipitation IF ( precipitation ) THEN CALL impact_of_latent_heat( i, j ) 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 influence 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) = ( 1.0-tsc(1) ) * pt_m(k,j,i) + tsc(1)*pt(k,j,i) +& dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tpt_m(k,j,i) & ) - & tsc(5) * rdf(k) * ( pt(k,j,i) - pt_init(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_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 ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN ! !-- Tendency-terms for salinity tend(:,j,i) = 0.0 IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) & THEN CALL advec_s_pw( i, j, sa ) ELSE CALL advec_s_up( i, j, sa ) ENDIF CALL diffusion_s( i, j, ddzu, ddzw, kh, sa, saswsb, saswst, & wall_salinityflux, tend ) 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) = tsc(1) * sa(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tsa_m(k,j,i) & ) - & tsc(5) * rdf(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 ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) & THEN CALL advec_s_pw( i, j, q ) ELSE CALL advec_s_up( i, j, q ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' )& THEN CALL diffusion_s( i, j, ddzu, ddzw, kh_m, q_m, qsws_m, & qswst_m, wall_qflux, tend ) ELSE CALL diffusion_s( i, j, ddzu, ddzw, kh, q, qsws, qswst, & wall_qflux, tend ) ENDIF ! !-- If required compute decrease of total water content due to !-- precipitation IF ( precipitation ) THEN CALL calc_precipitation( i, j ) 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) = ( 1.0-tsc(1) ) * q_m(k,j,i) + tsc(1)*q(k,j,i) +& dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tq_m(k,j,i) & ) - & tsc(5) * rdf(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 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 ( ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) & .AND. .NOT. use_upstream_for_tke ) THEN CALL advec_s_pw( i, j, e ) ELSE CALL advec_s_up( i, j, e ) ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' )& THEN IF ( .NOT. humidity ) THEN CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e_m, & km_m, l_grid, pt_m, pt_reference, & rif_m, tend, zu, zw ) ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e_m, & km_m, l_grid, vpt_m, pt_reference, & rif_m, tend, zu, zw ) ENDIF ELSE IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, & km, l_grid, prho, prho_reference, & rif, tend, zu, zw ) ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, & km, l_grid, pt, pt_reference, rif, & tend, zu, zw ) ENDIF ELSE CALL diffusion_e( i, j, ddzu, dd2zu, ddzw, diss, e, km, & l_grid, vpt, pt_reference, rif, tend, & zu, zw ) ENDIF 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) = ( 1.0-tsc(1) ) * e_m(k,j,i) + tsc(1)*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 CHARACTER (LEN=9) :: time_to_string INTEGER :: i, j, k REAL :: sat, sbt ! !-- Calculate those variables needed in the tendency terms which need !-- global communication CALL calc_mean_profile( pt, 4 ) IF ( ocean ) CALL calc_mean_profile( rho, 64 ) IF ( humidity ) CALL calc_mean_profile( vpt, 44 ) ! !-- u-velocity component CALL cpu_log( log_point(5), 'u-equation', 'start' ) ! !-- u-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_u_ups ENDIF ! !-- u-tendency terms with no communication IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_u_pw ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_u_up ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_u( ddzu, ddzw, km_m, km_damp_y, tend, u_m, usws_m, & uswst_m, v_m, w_m ) ELSE CALL diffusion_u( ddzu, ddzw, km, km_damp_y, tend, u, usws, uswst, v, w ) ENDIF CALL coriolis( 1 ) IF ( sloping_surface ) CALL buoyancy( pt, pt_reference, 1, 4 ) ! !-- 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) = ( 1.0-tsc(1) ) * u_m(k,j,i) + tsc(1) * u(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tu_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j,i-1) ) * ddx & ) - & 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' ) ! !-- v-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_v_ups ENDIF ! !-- v-tendency terms with no communication IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_v_pw ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_v_up ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_v( ddzu, ddzw, km_m, km_damp_x, tend, u_m, v_m, vsws_m, & vswst_m, w_m ) ELSE CALL diffusion_v( ddzu, ddzw, km, km_damp_x, tend, u, v, vsws, vswst, w ) ENDIF 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) = ( 1.0-tsc(1) ) * v_m(k,j,i) + tsc(1) * v(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tv_m(k,j,i) & - tsc(4) * ( p(k,j,i) - p(k,j-1,i) ) * ddy & ) - & 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' ) ! !-- w-tendency terms with communication IF ( momentum_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_w_ups ENDIF ! !-- w-tendency terms with no communication IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_w_pw ELSE IF ( momentum_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_w_up ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_w( ddzu, ddzw, km_m, km_damp_x, km_damp_y, tend, u_m, & v_m, w_m ) ELSE CALL diffusion_w( ddzu, ddzw, km, km_damp_x, km_damp_y, tend, u, v, w ) ENDIF CALL coriolis( 3 ) IF ( ocean ) THEN CALL buoyancy( rho, rho_reference, 3, 64 ) ELSE IF ( .NOT. humidity ) THEN CALL buoyancy( pt, pt_reference, 3, 4 ) ELSE CALL buoyancy( vpt, pt_reference, 3, 44 ) 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) = ( 1-tsc(1) ) * w_m(k,j,i) + tsc(1) * w(k,j,i) + & dt_3d * ( & tsc(2) * tend(k,j,i) + tsc(3) * tw_m(k,j,i) & - tsc(4) * ( p(k+1,j,i) - p(k,j,i) ) * ddzu(k+1) & ) - & 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' ) ! !-- potential temperature CALL cpu_log( log_point(13), 'pt-equation', 'start' ) ! !-- pt-tendency terms with communication sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( pt, 'pt' ) ELSE IF ( tsc(2) /= 2.0 .AND. scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( pt, 'pt' ) ENDIF ENDIF ! !-- pt-tendency terms with no communication IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( ddzu, ddzw, kh, pt, shf, tswst, wall_heatflux, & tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_s_pw( pt ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_up( pt ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_s( ddzu, ddzw, kh_m, pt_m, shf_m, tswst_m, & wall_heatflux, tend ) ELSE CALL diffusion_s( ddzu, ddzw, kh, pt, shf, tswst, wall_heatflux, & tend ) ENDIF ENDIF ! !-- 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) = ( 1 - sat ) * pt_m(k,j,i) + sat * pt(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tpt_m(k,j,i) & ) - & tsc(5) * rdf(k) * ( pt(k,j,i) - pt_init(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 = 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' ) ! !-- If required, compute prognostic equation for salinity IF ( ocean ) THEN CALL cpu_log( log_point(37), 'sa-equation', 'start' ) ! !-- sa-tendency terms with communication sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( sa, 'sa' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( sa, 'sa' ) ENDIF ENDIF ENDIF ! !-- sa-tendency terms with no communication IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( ddzu, ddzw, kh, sa, saswsb, saswst, & wall_salinityflux, tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_s_pw( sa ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_up( sa ) ENDIF ENDIF CALL diffusion_s( ddzu, ddzw, kh, sa, saswsb, saswst, & wall_salinityflux, tend ) ENDIF 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) = sat * sa(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tsa_m(k,j,i) & ) - & tsc(5) * rdf(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 sat = tsc(1) sbt = tsc(2) IF ( scalar_advec == 'bc-scheme' ) THEN IF ( timestep_scheme(1:5) /= 'runge' ) THEN ! !-- Bott-Chlond scheme always uses Euler time step when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( q, 'q' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( q, 'q' ) ENDIF ENDIF ENDIF ! !-- Scalar/q-tendency terms with no communication IF ( scalar_advec == 'bc-scheme' ) THEN CALL diffusion_s( ddzu, ddzw, kh, q, qsws, qswst, wall_qflux, tend ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_s_pw( q ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_up( q ) ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN CALL diffusion_s( ddzu, ddzw, kh_m, q_m, qsws_m, qswst_m, & wall_qflux, tend ) ELSE CALL diffusion_s( ddzu, ddzw, kh, q, qsws, qswst, & wall_qflux, tend ) ENDIF ENDIF ! !-- 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) = ( 1 - sat ) * q_m(k,j,i) + sat * q(k,j,i) + & dt_3d * ( & sbt * tend(k,j,i) + tsc(3) * tq_m(k,j,i) & ) - & tsc(5) * rdf(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' ) ! !-- TKE-tendency terms with communication CALL production_e_init sat = tsc(1) 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 when leapfrog is !-- switched on. Thus: sat = 1.0 sbt = 1.0 ENDIF tend = 0.0 CALL advec_s_bc( e, 'e' ) ELSE IF ( tsc(2) /= 2.0 ) THEN IF ( scalar_advec == 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_ups( e, 'e' ) ENDIF ENDIF ENDIF ENDIF ! !-- TKE-tendency terms with no communication IF ( scalar_advec == 'bc-scheme' .AND. .NOT. use_upstream_for_tke ) & THEN IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, & prho, prho_reference, rif, tend, zu, zw ) ELSE CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, pt, & pt_reference, rif, tend, zu, zw ) ENDIF ELSE CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, vpt, & pt_reference, rif, tend, zu, zw ) ENDIF ELSE IF ( use_upstream_for_tke ) THEN tend = 0.0 CALL advec_s_up( e ) ELSE IF ( tsc(2) == 2.0 .OR. timestep_scheme(1:5) == 'runge' ) THEN tend = 0.0 CALL advec_s_pw( e ) ELSE IF ( scalar_advec /= 'ups-scheme' ) THEN tend = 0.0 CALL advec_s_up( e ) ENDIF ENDIF ENDIF IF ( tsc(2) == 2.0 .AND. timestep_scheme(1:8) == 'leapfrog' ) THEN IF ( .NOT. humidity ) THEN CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e_m, km_m, l_grid, & pt_m, pt_reference, rif_m, tend, zu, zw ) ELSE CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e_m, km_m, l_grid, & vpt_m, pt_reference, rif_m, tend, zu, zw ) ENDIF ELSE IF ( .NOT. humidity ) THEN IF ( ocean ) THEN CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, & prho, prho_reference, rif, tend, zu, zw ) ELSE CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, & pt, pt_reference, rif, tend, zu, zw ) ENDIF ELSE CALL diffusion_e( ddzu, dd2zu, ddzw, diss, e, km, l_grid, vpt, & pt_reference, rif, tend, zu, zw ) ENDIF ENDIF 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) = ( 1 - sat ) * e_m(k,j,i) + sat * 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 END MODULE prognostic_equations_mod