[1873] | 1 | !> @file coriolis.f90 |
---|
[1036] | 2 | !--------------------------------------------------------------------------------! |
---|
| 3 | ! This file is part of PALM. |
---|
| 4 | ! |
---|
| 5 | ! PALM is free software: you can redistribute it and/or modify it under the terms |
---|
| 6 | ! of the GNU General Public License as published by the Free Software Foundation, |
---|
| 7 | ! either version 3 of the License, or (at your option) any later version. |
---|
| 8 | ! |
---|
| 9 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 10 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 11 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 12 | ! |
---|
| 13 | ! You should have received a copy of the GNU General Public License along with |
---|
| 14 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 15 | ! |
---|
[1818] | 16 | ! Copyright 1997-2016 Leibniz Universitaet Hannover |
---|
[1036] | 17 | !--------------------------------------------------------------------------------! |
---|
| 18 | ! |
---|
[254] | 19 | ! Current revisions: |
---|
[1] | 20 | ! ----------------- |
---|
[1873] | 21 | ! Module renamed (removed _mod) |
---|
[1354] | 22 | ! |
---|
[1683] | 23 | ! |
---|
[1321] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: coriolis.f90 1873 2016-04-18 14:50:06Z maronga $ |
---|
| 27 | ! |
---|
[1851] | 28 | ! 1850 2016-04-08 13:29:27Z maronga |
---|
| 29 | ! Module renamed |
---|
| 30 | ! |
---|
| 31 | ! |
---|
[1683] | 32 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 33 | ! Code annotations made doxygen readable |
---|
| 34 | ! |
---|
[1354] | 35 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
| 36 | ! REAL constants provided with KIND-attribute |
---|
| 37 | ! |
---|
[1321] | 38 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 39 | ! ONLY-attribute added to USE-statements, |
---|
| 40 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 41 | ! kinds are defined in new module kinds, |
---|
| 42 | ! revision history before 2012 removed, |
---|
| 43 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 44 | ! all variable declaration statements |
---|
[1321] | 45 | ! |
---|
[1258] | 46 | ! 1257 2013-11-08 15:18:40Z raasch |
---|
| 47 | ! openacc loop and loop vector clauses removed |
---|
| 48 | ! |
---|
[1132] | 49 | ! 1128 2013-04-12 06:19:32Z raasch |
---|
| 50 | ! loop index bounds in accelerator version replaced by i_left, i_right, j_south, |
---|
| 51 | ! j_north |
---|
| 52 | ! |
---|
[1037] | 53 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 54 | ! code put under GPL (PALM 3.9) |
---|
| 55 | ! |
---|
[1017] | 56 | ! 1015 2012-09-27 09:23:24Z raasch |
---|
| 57 | ! accelerator version (*_acc) added |
---|
| 58 | ! |
---|
[1] | 59 | ! Revision 1.1 1997/08/29 08:57:38 raasch |
---|
| 60 | ! Initial revision |
---|
| 61 | ! |
---|
| 62 | ! |
---|
| 63 | ! Description: |
---|
| 64 | ! ------------ |
---|
[1682] | 65 | !> Computation of all Coriolis terms in the equations of motion. |
---|
[1] | 66 | !------------------------------------------------------------------------------! |
---|
[1682] | 67 | MODULE coriolis_mod |
---|
| 68 | |
---|
[1] | 69 | |
---|
| 70 | PRIVATE |
---|
[1015] | 71 | PUBLIC coriolis, coriolis_acc |
---|
[1] | 72 | |
---|
| 73 | INTERFACE coriolis |
---|
| 74 | MODULE PROCEDURE coriolis |
---|
| 75 | MODULE PROCEDURE coriolis_ij |
---|
| 76 | END INTERFACE coriolis |
---|
| 77 | |
---|
[1015] | 78 | INTERFACE coriolis_acc |
---|
| 79 | MODULE PROCEDURE coriolis_acc |
---|
| 80 | END INTERFACE coriolis_acc |
---|
| 81 | |
---|
[1] | 82 | CONTAINS |
---|
| 83 | |
---|
| 84 | |
---|
| 85 | !------------------------------------------------------------------------------! |
---|
[1682] | 86 | ! Description: |
---|
| 87 | ! ------------ |
---|
| 88 | !> Call for all grid points |
---|
[1] | 89 | !------------------------------------------------------------------------------! |
---|
| 90 | SUBROUTINE coriolis( component ) |
---|
| 91 | |
---|
[1320] | 92 | USE arrays_3d, & |
---|
| 93 | ONLY: tend, u, ug, v, vg, w |
---|
| 94 | |
---|
| 95 | USE control_parameters, & |
---|
| 96 | ONLY: f, fs, message_string |
---|
| 97 | |
---|
| 98 | USE indices, & |
---|
| 99 | ONLY: nxl, nxlu, nxr, nyn, nys, nysv, nzb_u_inner, nzb_v_inner, & |
---|
| 100 | nzb_w_inner, nzt |
---|
| 101 | |
---|
| 102 | USE kinds |
---|
[1] | 103 | |
---|
| 104 | IMPLICIT NONE |
---|
| 105 | |
---|
[1682] | 106 | INTEGER(iwp) :: component !< |
---|
| 107 | INTEGER(iwp) :: i !< |
---|
| 108 | INTEGER(iwp) :: j !< |
---|
| 109 | INTEGER(iwp) :: k !< |
---|
[1] | 110 | |
---|
| 111 | |
---|
| 112 | ! |
---|
| 113 | !-- Compute Coriolis terms for the three velocity components |
---|
| 114 | SELECT CASE ( component ) |
---|
| 115 | |
---|
| 116 | ! |
---|
| 117 | !-- u-component |
---|
| 118 | CASE ( 1 ) |
---|
[106] | 119 | DO i = nxlu, nxr |
---|
[1] | 120 | DO j = nys, nyn |
---|
| 121 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
[1353] | 122 | tend(k,j,i) = tend(k,j,i) + f * ( 0.25_wp * & |
---|
[1320] | 123 | ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & |
---|
| 124 | v(k,j+1,i) ) - vg(k) ) & |
---|
[1353] | 125 | - fs * ( 0.25_wp * & |
---|
[1320] | 126 | ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + & |
---|
[1] | 127 | w(k,j,i) ) & |
---|
| 128 | ) |
---|
| 129 | ENDDO |
---|
| 130 | ENDDO |
---|
| 131 | ENDDO |
---|
| 132 | |
---|
| 133 | ! |
---|
| 134 | !-- v-component |
---|
| 135 | CASE ( 2 ) |
---|
| 136 | DO i = nxl, nxr |
---|
[106] | 137 | DO j = nysv, nyn |
---|
[1] | 138 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
[1353] | 139 | tend(k,j,i) = tend(k,j,i) - f * ( 0.25_wp * & |
---|
[1320] | 140 | ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & |
---|
[1] | 141 | u(k,j,i+1) ) - ug(k) ) |
---|
| 142 | ENDDO |
---|
| 143 | ENDDO |
---|
| 144 | ENDDO |
---|
| 145 | |
---|
| 146 | ! |
---|
| 147 | !-- w-component |
---|
| 148 | CASE ( 3 ) |
---|
| 149 | DO i = nxl, nxr |
---|
| 150 | DO j = nys, nyn |
---|
| 151 | DO k = nzb_w_inner(j,i)+1, nzt |
---|
[1353] | 152 | tend(k,j,i) = tend(k,j,i) + fs * 0.25_wp * & |
---|
[1320] | 153 | ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & |
---|
[1] | 154 | u(k+1,j,i+1) ) |
---|
| 155 | ENDDO |
---|
| 156 | ENDDO |
---|
| 157 | ENDDO |
---|
| 158 | |
---|
| 159 | CASE DEFAULT |
---|
| 160 | |
---|
[254] | 161 | WRITE( message_string, * ) ' wrong component: ', component |
---|
| 162 | CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) |
---|
[1] | 163 | |
---|
| 164 | END SELECT |
---|
| 165 | |
---|
| 166 | END SUBROUTINE coriolis |
---|
| 167 | |
---|
| 168 | |
---|
| 169 | !------------------------------------------------------------------------------! |
---|
[1682] | 170 | ! Description: |
---|
| 171 | ! ------------ |
---|
| 172 | !> Call for all grid points - accelerator version |
---|
[1015] | 173 | !------------------------------------------------------------------------------! |
---|
| 174 | SUBROUTINE coriolis_acc( component ) |
---|
| 175 | |
---|
[1320] | 176 | USE arrays_3d, & |
---|
| 177 | ONLY: tend, u, ug, v, vg, w |
---|
| 178 | |
---|
| 179 | USE control_parameters, & |
---|
| 180 | ONLY: f, fs, message_string |
---|
| 181 | |
---|
| 182 | USE indices, & |
---|
| 183 | ONLY: i_left, i_right, j_north, j_south, nzb_u_inner, & |
---|
| 184 | nzb_v_inner, nzb_w_inner, nzt |
---|
| 185 | |
---|
| 186 | USE kinds |
---|
[1015] | 187 | |
---|
| 188 | IMPLICIT NONE |
---|
| 189 | |
---|
[1682] | 190 | INTEGER(iwp) :: component !< |
---|
| 191 | INTEGER(iwp) :: i !< |
---|
| 192 | INTEGER(iwp) :: j !< |
---|
| 193 | INTEGER(iwp) :: k !< |
---|
[1015] | 194 | |
---|
| 195 | |
---|
| 196 | ! |
---|
| 197 | !-- Compute Coriolis terms for the three velocity components |
---|
| 198 | SELECT CASE ( component ) |
---|
| 199 | |
---|
| 200 | ! |
---|
| 201 | !-- u-component |
---|
| 202 | CASE ( 1 ) |
---|
| 203 | !$acc kernels present( nzb_u_inner, tend, v, vg, w ) |
---|
[1128] | 204 | DO i = i_left, i_right |
---|
| 205 | DO j = j_south, j_north |
---|
[1015] | 206 | DO k = 1, nzt |
---|
| 207 | IF ( k > nzb_u_inner(j,i) ) THEN |
---|
[1353] | 208 | tend(k,j,i) = tend(k,j,i) + f * ( 0.25_wp * & |
---|
[1015] | 209 | ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & |
---|
| 210 | v(k,j+1,i) ) - vg(k) ) & |
---|
[1353] | 211 | - fs * ( 0.25_wp * & |
---|
[1015] | 212 | ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) & |
---|
| 213 | + w(k,j,i) ) & |
---|
| 214 | ) |
---|
| 215 | ENDIF |
---|
| 216 | ENDDO |
---|
| 217 | ENDDO |
---|
| 218 | ENDDO |
---|
| 219 | !$acc end kernels |
---|
| 220 | |
---|
| 221 | ! |
---|
| 222 | !-- v-component |
---|
| 223 | CASE ( 2 ) |
---|
| 224 | !$acc kernels present( nzb_v_inner, tend, u, ug ) |
---|
[1128] | 225 | DO i = i_left, i_right |
---|
| 226 | DO j = j_south, j_north |
---|
[1015] | 227 | DO k = 1, nzt |
---|
| 228 | IF ( k > nzb_v_inner(j,i) ) THEN |
---|
[1353] | 229 | tend(k,j,i) = tend(k,j,i) - f * ( 0.25_wp * & |
---|
[1015] | 230 | ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & |
---|
| 231 | u(k,j,i+1) ) - ug(k) ) |
---|
| 232 | ENDIF |
---|
| 233 | ENDDO |
---|
| 234 | ENDDO |
---|
| 235 | ENDDO |
---|
| 236 | !$acc end kernels |
---|
| 237 | |
---|
| 238 | ! |
---|
| 239 | !-- w-component |
---|
| 240 | CASE ( 3 ) |
---|
| 241 | !$acc kernels present( nzb_w_inner, tend, u ) |
---|
[1128] | 242 | DO i = i_left, i_right |
---|
| 243 | DO j = j_south, j_north |
---|
[1015] | 244 | DO k = 1, nzt |
---|
| 245 | IF ( k > nzb_w_inner(j,i) ) THEN |
---|
[1353] | 246 | tend(k,j,i) = tend(k,j,i) + fs * 0.25_wp * & |
---|
[1320] | 247 | ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & |
---|
[1015] | 248 | u(k+1,j,i+1) ) |
---|
| 249 | ENDIF |
---|
| 250 | ENDDO |
---|
| 251 | ENDDO |
---|
| 252 | ENDDO |
---|
| 253 | !$acc end kernels |
---|
| 254 | |
---|
| 255 | CASE DEFAULT |
---|
| 256 | |
---|
| 257 | WRITE( message_string, * ) ' wrong component: ', component |
---|
| 258 | CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) |
---|
| 259 | |
---|
| 260 | END SELECT |
---|
| 261 | |
---|
| 262 | END SUBROUTINE coriolis_acc |
---|
| 263 | |
---|
| 264 | |
---|
| 265 | !------------------------------------------------------------------------------! |
---|
[1682] | 266 | ! Description: |
---|
| 267 | ! ------------ |
---|
| 268 | !> Call for grid point i,j |
---|
[1] | 269 | !------------------------------------------------------------------------------! |
---|
| 270 | SUBROUTINE coriolis_ij( i, j, component ) |
---|
| 271 | |
---|
[1320] | 272 | USE arrays_3d, & |
---|
| 273 | ONLY: tend, u, ug, v, vg, w |
---|
| 274 | |
---|
| 275 | USE control_parameters, & |
---|
| 276 | ONLY: f, fs, message_string |
---|
| 277 | |
---|
| 278 | USE indices, & |
---|
| 279 | ONLY: nzb_u_inner, nzb_v_inner, nzb_w_inner, nzt |
---|
| 280 | |
---|
| 281 | USE kinds |
---|
| 282 | |
---|
[1] | 283 | IMPLICIT NONE |
---|
| 284 | |
---|
[1682] | 285 | INTEGER(iwp) :: component !< |
---|
| 286 | INTEGER(iwp) :: i !< |
---|
| 287 | INTEGER(iwp) :: j !< |
---|
| 288 | INTEGER(iwp) :: k !< |
---|
[1] | 289 | |
---|
| 290 | ! |
---|
| 291 | !-- Compute Coriolis terms for the three velocity components |
---|
| 292 | SELECT CASE ( component ) |
---|
| 293 | |
---|
| 294 | ! |
---|
| 295 | !-- u-component |
---|
| 296 | CASE ( 1 ) |
---|
| 297 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
[1353] | 298 | tend(k,j,i) = tend(k,j,i) + f * ( 0.25_wp * & |
---|
[1320] | 299 | ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) + & |
---|
| 300 | v(k,j+1,i) ) - vg(k) ) & |
---|
[1353] | 301 | - fs * ( 0.25_wp * & |
---|
[1320] | 302 | ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + & |
---|
| 303 | w(k,j,i) ) ) |
---|
[1] | 304 | ENDDO |
---|
| 305 | |
---|
| 306 | ! |
---|
| 307 | !-- v-component |
---|
| 308 | CASE ( 2 ) |
---|
| 309 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
[1353] | 310 | tend(k,j,i) = tend(k,j,i) - f * ( 0.25_wp * & |
---|
[1320] | 311 | ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + & |
---|
[1] | 312 | u(k,j,i+1) ) - ug(k) ) |
---|
| 313 | ENDDO |
---|
| 314 | |
---|
| 315 | ! |
---|
| 316 | !-- w-component |
---|
| 317 | CASE ( 3 ) |
---|
| 318 | DO k = nzb_w_inner(j,i)+1, nzt |
---|
[1353] | 319 | tend(k,j,i) = tend(k,j,i) + fs * 0.25_wp * & |
---|
[1320] | 320 | ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + & |
---|
[1] | 321 | u(k+1,j,i+1) ) |
---|
| 322 | ENDDO |
---|
| 323 | |
---|
| 324 | CASE DEFAULT |
---|
| 325 | |
---|
[254] | 326 | WRITE( message_string, * ) ' wrong component: ', component |
---|
| 327 | CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 ) |
---|
[1] | 328 | |
---|
| 329 | END SELECT |
---|
| 330 | |
---|
| 331 | END SUBROUTINE coriolis_ij |
---|
| 332 | |
---|
| 333 | END MODULE coriolis_mod |
---|