MODULE coriolis_mod !------------------------------------------------------------------------------! ! Actual revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Log: coriolis.f90,v $ ! Revision 1.12 2006/02/23 10:08:57 raasch ! nzb_2d replaced by nzb_u/v/w_inner ! ! Revision 1.11 2005/06/29 10:30:17 steinfeld ! A potential dependency of the geostrophic wind components ug and vg on height ! is considered ! ! Revision 1.10 2005/03/26 20:02:04 raasch ! Extension of horizontal loop upper bounds for non-cyclic boundary conditions ! ! Revision 1.9 2004/01/30 10:17:34 raasch ! Scalar lower k index nzb replaced by 2d-array nzb_2d ! ! Revision 1.8 2003/03/12 16:23:47 raasch ! Full code replaced in the call for all gridpoints instead of calling the ! _ij version (required by NEC, because otherwise no vectorization) ! ! Revision 1.7 2002/12/19 14:08:14 raasch ! STOP statement replaced by call of subroutine local_stop ! ! Revision 1.6 2002/06/11 12:51:09 raasch ! Former subroutine changed to a module which allows to be called for all grid ! points of a single vertical column with index i,j or for all grid points by ! using function overloading. ! ! Revision 1.5 2001/03/30 06:59:09 raasch ! Translation of remaining German identifiers (variables, subroutines, etc.) ! ! Revision 1.4 2001/01/22 05:44:16 raasch ! Module test_variables removed ! ! Revision 1.3 2000/12/20 09:50:27 letzel ! All comments translated into English. ! ! Revision 1.2 1998/07/06 12:09:54 raasch ! + USE test_variables ! ! 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 INTERFACE coriolis MODULE PROCEDURE coriolis MODULE PROCEDURE coriolis_ij END INTERFACE coriolis 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 = nxl, nxr+uxrp 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 = nys, nyn+vynp 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 IF ( myid == 0 ) PRINT*,'+++ coriolis: wrong component: ', & component CALL local_stop END SELECT END SUBROUTINE coriolis !------------------------------------------------------------------------------! ! 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 IF ( myid == 0 ) PRINT*,'+++ coriolis: wrong component: ', & component CALL local_stop END SELECT END SUBROUTINE coriolis_ij END MODULE coriolis_mod