MODULE coriolis_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: coriolis.f90 1132 2013-04-12 14:35:30Z witha $ ! ! 1128 2013-04-12 06:19:32Z raasch ! loop index bounds in accelerator version replaced by i_left, i_right, j_south, ! j_north ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1015 2012-09-27 09:23:24Z raasch ! accelerator version (*_acc) added ! ! 254 2009-03-05 15:33:42Z heinze ! Output of messages replaced by message handling routine. ! ! 106 2007-08-16 14:30:26Z raasch ! loops for u and v are starting from index nxlu, nysv, respectively (needed ! for non-cyclic boundary conditions) ! ! 75 2007-03-22 09:54:05Z raasch ! uxrp, vynp eliminated ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.12 2006/02/23 10:08:57 raasch ! nzb_2d replaced by nzb_u/v/w_inner ! ! Revision 1.1 1997/08/29 08:57:38 raasch ! Initial revision ! ! ! Description: ! ------------ ! Computation of all Coriolis terms in the equations of motion. !------------------------------------------------------------------------------! PRIVATE PUBLIC coriolis, coriolis_acc INTERFACE coriolis MODULE PROCEDURE coriolis MODULE PROCEDURE coriolis_ij END INTERFACE coriolis INTERFACE coriolis_acc MODULE PROCEDURE coriolis_acc END INTERFACE coriolis_acc CONTAINS !------------------------------------------------------------------------------! ! Call for all grid points !------------------------------------------------------------------------------! SUBROUTINE coriolis( component ) USE arrays_3d USE control_parameters USE indices USE pegrid IMPLICIT NONE INTEGER :: component, i, j, k ! !-- Compute Coriolis terms for the three velocity components SELECT CASE ( component ) ! !-- u-component CASE ( 1 ) DO i = nxlu, nxr DO j = nys, nyn DO k = nzb_u_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) + f * ( 0.25 * & ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & v(k,j+1,i) ) - vg(k) ) & - fs * ( 0.25 * & ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + & w(k,j,i) ) & ) ENDDO ENDDO ENDDO ! !-- v-component CASE ( 2 ) DO i = nxl, nxr DO j = nysv, nyn DO k = nzb_v_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) - f * ( 0.25 * & ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & u(k,j,i+1) ) - ug(k) ) ENDDO ENDDO ENDDO ! !-- w-component CASE ( 3 ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb_w_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) + fs * 0.25 * & ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & u(k+1,j,i+1) ) ENDDO ENDDO ENDDO CASE DEFAULT WRITE( message_string, * ) ' wrong component: ', component CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) END SELECT END SUBROUTINE coriolis !------------------------------------------------------------------------------! ! Call for all grid points - accelerator version !------------------------------------------------------------------------------! SUBROUTINE coriolis_acc( component ) USE arrays_3d USE control_parameters USE indices USE pegrid IMPLICIT NONE INTEGER :: component, i, j, k ! !-- Compute Coriolis terms for the three velocity components SELECT CASE ( component ) ! !-- u-component CASE ( 1 ) !$acc kernels present( nzb_u_inner, tend, v, vg, w ) !$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 tend(k,j,i) = tend(k,j,i) + f * ( 0.25 * & ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & v(k,j+1,i) ) - vg(k) ) & - fs * ( 0.25 * & ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) & + w(k,j,i) ) & ) ENDIF ENDDO ENDDO ENDDO !$acc end kernels ! !-- v-component CASE ( 2 ) !$acc kernels present( nzb_v_inner, tend, u, ug ) !$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 tend(k,j,i) = tend(k,j,i) - f * ( 0.25 * & ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & u(k,j,i+1) ) - ug(k) ) ENDIF ENDDO ENDDO ENDDO !$acc end kernels ! !-- w-component CASE ( 3 ) !$acc kernels present( nzb_w_inner, tend, u ) !$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_w_inner(j,i) ) THEN tend(k,j,i) = tend(k,j,i) + fs * 0.25 * & ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & u(k+1,j,i+1) ) ENDIF ENDDO ENDDO ENDDO !$acc end kernels CASE DEFAULT WRITE( message_string, * ) ' wrong component: ', component CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) END SELECT END SUBROUTINE coriolis_acc !------------------------------------------------------------------------------! ! Call for grid point i,j !------------------------------------------------------------------------------! SUBROUTINE coriolis_ij( i, j, component ) USE arrays_3d USE control_parameters USE indices USE pegrid IMPLICIT NONE INTEGER :: component, i, j, k ! !-- Compute Coriolis terms for the three velocity components SELECT CASE ( component ) ! !-- u-component CASE ( 1 ) DO k = nzb_u_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) + f * ( 0.25 * & ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & v(k,j+1,i) ) - vg(k) ) & - fs * ( 0.25 * & ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + & w(k,j,i) ) & ) ENDDO ! !-- v-component CASE ( 2 ) DO k = nzb_v_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) - f * ( 0.25 * & ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & u(k,j,i+1) ) - ug(k) ) ENDDO ! !-- w-component CASE ( 3 ) DO k = nzb_w_inner(j,i)+1, nzt tend(k,j,i) = tend(k,j,i) + fs * 0.25 * & ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & u(k+1,j,i+1) ) ENDDO CASE DEFAULT WRITE( message_string, * ) ' wrong component: ', component CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) END SELECT END SUBROUTINE coriolis_ij END MODULE coriolis_mod