[1682] | 1 | !> @file pres.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[1036] | 3 | ! This file is part of PALM. |
---|
| 4 | ! |
---|
[2000] | 5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
| 6 | ! terms of the GNU General Public License as published by the Free Software |
---|
| 7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
| 8 | ! version. |
---|
[1036] | 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
[1818] | 17 | ! Copyright 1997-2016 Leibniz Universitaet Hannover |
---|
[2000] | 18 | !------------------------------------------------------------------------------! |
---|
[1036] | 19 | ! |
---|
[484] | 20 | ! Current revisions: |
---|
[1212] | 21 | ! ------------------ |
---|
[2000] | 22 | ! Forced header and separation lines into 80 columns |
---|
[1932] | 23 | ! |
---|
[1919] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: pres.f90 2000 2016-08-20 18:09:15Z knoop $ |
---|
| 27 | ! |
---|
[1933] | 28 | ! 1932 2016-06-10 12:09:21Z suehring |
---|
| 29 | ! Initial version of purely vertical nesting introduced. |
---|
| 30 | ! |
---|
[1932] | 31 | ! 1931 2016-06-10 12:06:59Z suehring |
---|
| 32 | ! Rename multigrid into multigrid_noopt and multigrid_fast into multigrid |
---|
| 33 | ! |
---|
[1930] | 34 | ! 1929 2016-06-09 16:25:25Z suehring |
---|
| 35 | ! Bugfix: weight_substep for initial call, replace by local variable |
---|
| 36 | ! |
---|
[1919] | 37 | ! 1918 2016-05-27 14:35:57Z raasch |
---|
[1918] | 38 | ! sum of divergence is also calculated when pres is called before the initial |
---|
| 39 | ! first time step, |
---|
| 40 | ! stearing is modified by using intermediate_timestep_count = 0 in order to |
---|
| 41 | ! determine if pres is called before the first initial timestep, |
---|
| 42 | ! bugfix: in case of Neumann conditions for pressure both at bottom and top, |
---|
| 43 | ! mean vertical velocity is also removed for the first time step |
---|
| 44 | ! bugfix for calculating divergences |
---|
[1321] | 45 | ! |
---|
[1909] | 46 | ! 1908 2016-05-25 17:22:32Z suehring |
---|
| 47 | ! New divergence for output into RUN_CONTROL file is calculated only at last |
---|
| 48 | ! Runge-Kutta step |
---|
| 49 | ! |
---|
[1846] | 50 | ! 1845 2016-04-08 08:29:13Z raasch |
---|
| 51 | ! nzb_2d replace by nzb_u|v_inner |
---|
| 52 | ! |
---|
[1800] | 53 | ! 1799 2016-04-05 08:35:55Z gronemeier |
---|
| 54 | ! Bugfix: excluded third dimension from horizontal volume flow calculation |
---|
| 55 | ! |
---|
[1763] | 56 | ! 1762 2016-02-25 12:31:13Z hellstea |
---|
| 57 | ! Introduction of nested domain feature |
---|
| 58 | ! |
---|
[1683] | 59 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 60 | ! Code annotations made doxygen readable |
---|
| 61 | ! |
---|
[1576] | 62 | ! 1575 2015-03-27 09:56:27Z raasch |
---|
| 63 | ! poismg_fast + respective module added, adjustments for psolver-queries |
---|
| 64 | ! |
---|
[1343] | 65 | ! 1342 2014-03-26 17:04:47Z kanani |
---|
| 66 | ! REAL constants defined as wp-kind |
---|
| 67 | ! |
---|
[1321] | 68 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 69 | ! ONLY-attribute added to USE-statements, |
---|
| 70 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 71 | ! kinds are defined in new module kinds, |
---|
| 72 | ! old module precision_kind is removed, |
---|
| 73 | ! revision history before 2012 removed, |
---|
| 74 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 75 | ! all variable declaration statements |
---|
[708] | 76 | ! |
---|
[1319] | 77 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
| 78 | ! module interfaces removed |
---|
| 79 | ! |
---|
[1307] | 80 | ! 1306 2014-03-13 14:30:59Z raasch |
---|
| 81 | ! second argument removed from routine poisfft |
---|
| 82 | ! |
---|
[1258] | 83 | ! 1257 2013-11-08 15:18:40Z raasch |
---|
| 84 | ! openacc loop and loop vector clauses removed, independent clauses added, |
---|
| 85 | ! end parallel replaced by end parallel loop |
---|
| 86 | ! |
---|
[1222] | 87 | ! 1221 2013-09-10 08:59:13Z raasch |
---|
| 88 | ! openACC porting of reduction operations, loops for calculating d are |
---|
| 89 | ! using the rflags_s_inner multiply flag instead of the nzb_s_inner loop index |
---|
| 90 | ! |
---|
[1213] | 91 | ! 1212 2013-08-15 08:46:27Z raasch |
---|
| 92 | ! call of poisfft_hybrid removed |
---|
| 93 | ! |
---|
[1118] | 94 | ! 1117 2013-03-27 11:15:36Z suehring |
---|
| 95 | ! Bugfix in OpenMP parallelization. |
---|
| 96 | ! |
---|
[1114] | 97 | ! 1113 2013-03-10 02:48:14Z raasch |
---|
| 98 | ! GPU-porting of several loops, some loops rearranged |
---|
| 99 | ! |
---|
[1112] | 100 | ! 1111 2013-03-08 23:54:10Z |
---|
| 101 | ! openACC statements added, |
---|
| 102 | ! ibc_p_b = 2 removed |
---|
| 103 | ! |
---|
[1093] | 104 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 105 | ! unused variables removed |
---|
| 106 | ! |
---|
[1037] | 107 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 108 | ! code put under GPL (PALM 3.9) |
---|
| 109 | ! |
---|
[1004] | 110 | ! 1003 2012-09-14 14:35:53Z raasch |
---|
| 111 | ! adjustment of array tend for cases with unequal subdomain sizes removed |
---|
| 112 | ! |
---|
[1] | 113 | ! Revision 1.1 1997/07/24 11:24:44 raasch |
---|
| 114 | ! Initial revision |
---|
| 115 | ! |
---|
| 116 | ! |
---|
| 117 | ! Description: |
---|
| 118 | ! ------------ |
---|
[1682] | 119 | !> Compute the divergence of the provisional velocity field. Solve the Poisson |
---|
| 120 | !> equation for the perturbation pressure. Compute the final velocities using |
---|
| 121 | !> this perturbation pressure. Compute the remaining divergence. |
---|
[1] | 122 | !------------------------------------------------------------------------------! |
---|
[1682] | 123 | SUBROUTINE pres |
---|
| 124 | |
---|
[1] | 125 | |
---|
[1320] | 126 | USE arrays_3d, & |
---|
| 127 | ONLY: d, ddzu, ddzu_pres, ddzw, dzw, p, p_loc, tend, u, v, w |
---|
| 128 | |
---|
| 129 | USE control_parameters, & |
---|
| 130 | ONLY: bc_lr_cyc, bc_ns_cyc, conserve_volume_flow, dt_3d, & |
---|
| 131 | gathered_size, ibc_p_b, ibc_p_t, intermediate_timestep_count, & |
---|
[1908] | 132 | intermediate_timestep_count_max, mg_switch_to_pe0_level, & |
---|
[1933] | 133 | nest_domain, on_device, outflow_l, outflow_n, outflow_r, & |
---|
[1918] | 134 | outflow_s, psolver, subdomain_size, topography, volume_flow, & |
---|
| 135 | volume_flow_area, volume_flow_initial |
---|
[1320] | 136 | |
---|
| 137 | USE cpulog, & |
---|
| 138 | ONLY: cpu_log, log_point, log_point_s |
---|
| 139 | |
---|
| 140 | USE grid_variables, & |
---|
| 141 | ONLY: ddx, ddy |
---|
| 142 | |
---|
| 143 | USE indices, & |
---|
| 144 | ONLY: nbgp, ngp_2dh_outer, nx, nxl, nxlg, nxl_mg, nxr, nxrg, nxr_mg, & |
---|
| 145 | ny, nys, nysg, nys_mg, nyn, nyng, nyn_mg, nzb, nzb_s_inner, & |
---|
[1845] | 146 | nzb_u_inner, nzb_v_inner, nzb_w_inner, nzt, nzt_mg, & |
---|
[1320] | 147 | rflags_s_inner |
---|
| 148 | |
---|
| 149 | USE kinds |
---|
| 150 | |
---|
[1] | 151 | USE pegrid |
---|
[1933] | 152 | |
---|
| 153 | USE pmc_interface, & |
---|
| 154 | ONLY: nesting_mode |
---|
[1] | 155 | |
---|
[1320] | 156 | USE poisfft_mod, & |
---|
| 157 | ONLY: poisfft |
---|
| 158 | |
---|
[1575] | 159 | USE poismg_mod |
---|
| 160 | |
---|
[1320] | 161 | USE statistics, & |
---|
| 162 | ONLY: statistic_regions, sums_divnew_l, sums_divold_l, weight_pres, & |
---|
| 163 | weight_substep |
---|
| 164 | |
---|
[1] | 165 | IMPLICIT NONE |
---|
| 166 | |
---|
[1682] | 167 | INTEGER(iwp) :: i !< |
---|
| 168 | INTEGER(iwp) :: j !< |
---|
| 169 | INTEGER(iwp) :: k !< |
---|
[1] | 170 | |
---|
[1682] | 171 | REAL(wp) :: ddt_3d !< |
---|
[1918] | 172 | REAL(wp) :: d_weight_pres !< |
---|
[1682] | 173 | REAL(wp) :: localsum !< |
---|
| 174 | REAL(wp) :: threadsum !< |
---|
[1918] | 175 | REAL(wp) :: weight_pres_l !< |
---|
[1929] | 176 | REAL(wp) :: weight_substep_l !< |
---|
[1] | 177 | |
---|
[1762] | 178 | REAL(wp), DIMENSION(1:3) :: volume_flow_l !< |
---|
| 179 | REAL(wp), DIMENSION(1:3) :: volume_flow_offset !< |
---|
[1682] | 180 | REAL(wp), DIMENSION(1:nzt) :: w_l !< |
---|
| 181 | REAL(wp), DIMENSION(1:nzt) :: w_l_l !< |
---|
[1] | 182 | |
---|
[1933] | 183 | LOGICAL :: nest_domain_nvn !< |
---|
[1] | 184 | |
---|
[1933] | 185 | |
---|
[1] | 186 | CALL cpu_log( log_point(8), 'pres', 'start' ) |
---|
| 187 | |
---|
[85] | 188 | |
---|
[1918] | 189 | ! |
---|
| 190 | !-- Calculate quantities to be used locally |
---|
[1342] | 191 | ddt_3d = 1.0_wp / dt_3d |
---|
[1918] | 192 | IF ( intermediate_timestep_count == 0 ) THEN |
---|
| 193 | ! |
---|
| 194 | !-- If pres is called before initial time step |
---|
[1929] | 195 | weight_pres_l = 1.0_wp |
---|
| 196 | d_weight_pres = 1.0_wp |
---|
| 197 | weight_substep_l = 1.0_wp |
---|
[1918] | 198 | ELSE |
---|
[1929] | 199 | weight_pres_l = weight_pres(intermediate_timestep_count) |
---|
| 200 | d_weight_pres = 1.0_wp / weight_pres(intermediate_timestep_count) |
---|
| 201 | weight_substep_l = weight_substep(intermediate_timestep_count) |
---|
[1918] | 202 | ENDIF |
---|
[85] | 203 | |
---|
[1] | 204 | ! |
---|
[707] | 205 | !-- Multigrid method expects array d to have one ghost layer. |
---|
| 206 | !-- |
---|
[1575] | 207 | IF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[667] | 208 | |
---|
[1] | 209 | DEALLOCATE( d ) |
---|
[667] | 210 | ALLOCATE( d(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) |
---|
[707] | 211 | |
---|
| 212 | ! |
---|
| 213 | !-- Since p is later used to hold the weighted average of the substeps, it |
---|
| 214 | !-- cannot be used in the iterative solver. Therefore, its initial value is |
---|
| 215 | !-- stored on p_loc, which is then iteratively advanced in every substep. |
---|
[1762] | 216 | IF ( intermediate_timestep_count <= 1 ) THEN |
---|
[707] | 217 | DO i = nxl-1, nxr+1 |
---|
| 218 | DO j = nys-1, nyn+1 |
---|
| 219 | DO k = nzb, nzt+1 |
---|
| 220 | p_loc(k,j,i) = p(k,j,i) |
---|
| 221 | ENDDO |
---|
| 222 | ENDDO |
---|
| 223 | ENDDO |
---|
[667] | 224 | ENDIF |
---|
| 225 | |
---|
[1918] | 226 | ELSEIF ( psolver == 'sor' .AND. intermediate_timestep_count <= 1 ) THEN |
---|
[707] | 227 | |
---|
| 228 | ! |
---|
| 229 | !-- Since p is later used to hold the weighted average of the substeps, it |
---|
| 230 | !-- cannot be used in the iterative solver. Therefore, its initial value is |
---|
| 231 | !-- stored on p_loc, which is then iteratively advanced in every substep. |
---|
| 232 | p_loc = p |
---|
| 233 | |
---|
[1] | 234 | ENDIF |
---|
| 235 | |
---|
| 236 | ! |
---|
[75] | 237 | !-- Conserve the volume flow at the outflow in case of non-cyclic lateral |
---|
| 238 | !-- boundary conditions |
---|
[106] | 239 | !-- WARNING: so far, this conservation does not work at the left/south |
---|
| 240 | !-- boundary if the topography at the inflow differs from that at the |
---|
| 241 | !-- outflow! For this case, volume_flow_area needs adjustment! |
---|
| 242 | ! |
---|
| 243 | !-- Left/right |
---|
[709] | 244 | IF ( conserve_volume_flow .AND. ( outflow_l .OR. outflow_r ) ) THEN |
---|
[680] | 245 | |
---|
[1342] | 246 | volume_flow(1) = 0.0_wp |
---|
| 247 | volume_flow_l(1) = 0.0_wp |
---|
[106] | 248 | |
---|
| 249 | IF ( outflow_l ) THEN |
---|
| 250 | i = 0 |
---|
| 251 | ELSEIF ( outflow_r ) THEN |
---|
| 252 | i = nx+1 |
---|
| 253 | ENDIF |
---|
| 254 | |
---|
| 255 | DO j = nys, nyn |
---|
| 256 | ! |
---|
| 257 | !-- Sum up the volume flow through the south/north boundary |
---|
[1845] | 258 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
[667] | 259 | volume_flow_l(1) = volume_flow_l(1) + u(k,j,i) * dzw(k) |
---|
[106] | 260 | ENDDO |
---|
| 261 | ENDDO |
---|
| 262 | |
---|
| 263 | #if defined( __parallel ) |
---|
[680] | 264 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
[106] | 265 | CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 1, MPI_REAL, & |
---|
| 266 | MPI_SUM, comm1dy, ierr ) |
---|
| 267 | #else |
---|
| 268 | volume_flow = volume_flow_l |
---|
| 269 | #endif |
---|
[709] | 270 | volume_flow_offset(1) = ( volume_flow_initial(1) - volume_flow(1) ) & |
---|
[106] | 271 | / volume_flow_area(1) |
---|
| 272 | |
---|
[667] | 273 | DO j = nysg, nyng |
---|
[1845] | 274 | DO k = nzb_u_inner(j,i)+1, nzt |
---|
[106] | 275 | u(k,j,i) = u(k,j,i) + volume_flow_offset(1) |
---|
| 276 | ENDDO |
---|
| 277 | ENDDO |
---|
| 278 | |
---|
| 279 | ENDIF |
---|
| 280 | |
---|
| 281 | ! |
---|
| 282 | !-- South/north |
---|
[709] | 283 | IF ( conserve_volume_flow .AND. ( outflow_n .OR. outflow_s ) ) THEN |
---|
[106] | 284 | |
---|
[1342] | 285 | volume_flow(2) = 0.0_wp |
---|
| 286 | volume_flow_l(2) = 0.0_wp |
---|
[75] | 287 | |
---|
[106] | 288 | IF ( outflow_s ) THEN |
---|
| 289 | j = 0 |
---|
| 290 | ELSEIF ( outflow_n ) THEN |
---|
[75] | 291 | j = ny+1 |
---|
[106] | 292 | ENDIF |
---|
| 293 | |
---|
| 294 | DO i = nxl, nxr |
---|
[75] | 295 | ! |
---|
[106] | 296 | !-- Sum up the volume flow through the south/north boundary |
---|
[1845] | 297 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
[667] | 298 | volume_flow_l(2) = volume_flow_l(2) + v(k,j,i) * dzw(k) |
---|
[75] | 299 | ENDDO |
---|
[106] | 300 | ENDDO |
---|
| 301 | |
---|
[75] | 302 | #if defined( __parallel ) |
---|
[680] | 303 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
[75] | 304 | CALL MPI_ALLREDUCE( volume_flow_l(2), volume_flow(2), 1, MPI_REAL, & |
---|
| 305 | MPI_SUM, comm1dx, ierr ) |
---|
| 306 | #else |
---|
| 307 | volume_flow = volume_flow_l |
---|
| 308 | #endif |
---|
| 309 | volume_flow_offset(2) = ( volume_flow_initial(2) - volume_flow(2) ) & |
---|
[106] | 310 | / volume_flow_area(2) |
---|
[75] | 311 | |
---|
[667] | 312 | DO i = nxlg, nxrg |
---|
[709] | 313 | DO k = nzb_v_inner(j,i)+1, nzt |
---|
[106] | 314 | v(k,j,i) = v(k,j,i) + volume_flow_offset(2) |
---|
[75] | 315 | ENDDO |
---|
[106] | 316 | ENDDO |
---|
[75] | 317 | |
---|
| 318 | ENDIF |
---|
| 319 | |
---|
[76] | 320 | ! |
---|
[1762] | 321 | !-- Remove mean vertical velocity in case that Neumann conditions are |
---|
[1933] | 322 | !-- used both at bottom and top boundary, and if not a nested domain in a |
---|
| 323 | !-- normal nesting run. In case of vertical nesting, this must be done. |
---|
| 324 | !-- Therefore an auxiliary logical variable nest_domain_nvn is used here, and |
---|
| 325 | !-- nvn stands for non-vertical nesting. |
---|
[1918] | 326 | !-- This cannot be done before the first initial time step because ngp_2dh_outer |
---|
| 327 | !-- is not yet known then. |
---|
[1933] | 328 | nest_domain_nvn = nest_domain |
---|
| 329 | IF ( nest_domain .AND. nesting_mode == 'vertical' ) THEN |
---|
| 330 | nest_domain_nvn = .FALSE. |
---|
| 331 | ENDIF |
---|
| 332 | |
---|
| 333 | IF ( ibc_p_b == 1 .AND. ibc_p_t == 1 .AND. & |
---|
| 334 | .NOT. nest_domain_nvn .AND. intermediate_timestep_count /= 0 ) & |
---|
[1918] | 335 | THEN |
---|
| 336 | w_l = 0.0_wp; w_l_l = 0.0_wp |
---|
| 337 | DO i = nxl, nxr |
---|
| 338 | DO j = nys, nyn |
---|
| 339 | DO k = nzb_w_inner(j,i)+1, nzt |
---|
| 340 | w_l_l(k) = w_l_l(k) + w(k,j,i) |
---|
[76] | 341 | ENDDO |
---|
| 342 | ENDDO |
---|
[1918] | 343 | ENDDO |
---|
[76] | 344 | #if defined( __parallel ) |
---|
[1918] | 345 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 346 | CALL MPI_ALLREDUCE( w_l_l(1), w_l(1), nzt, MPI_REAL, MPI_SUM, & |
---|
| 347 | comm2d, ierr ) |
---|
[76] | 348 | #else |
---|
[1918] | 349 | w_l = w_l_l |
---|
[76] | 350 | #endif |
---|
[1918] | 351 | DO k = 1, nzt |
---|
| 352 | w_l(k) = w_l(k) / ngp_2dh_outer(k,0) |
---|
| 353 | ENDDO |
---|
| 354 | DO i = nxlg, nxrg |
---|
| 355 | DO j = nysg, nyng |
---|
| 356 | DO k = nzb_w_inner(j,i)+1, nzt |
---|
| 357 | w(k,j,i) = w(k,j,i) - w_l(k) |
---|
[76] | 358 | ENDDO |
---|
| 359 | ENDDO |
---|
[1918] | 360 | ENDDO |
---|
[76] | 361 | ENDIF |
---|
[75] | 362 | |
---|
| 363 | ! |
---|
[1] | 364 | !-- Compute the divergence of the provisional velocity field. |
---|
| 365 | CALL cpu_log( log_point_s(1), 'divergence', 'start' ) |
---|
| 366 | |
---|
[1575] | 367 | IF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[1] | 368 | !$OMP PARALLEL DO SCHEDULE( STATIC ) |
---|
| 369 | DO i = nxl-1, nxr+1 |
---|
| 370 | DO j = nys-1, nyn+1 |
---|
| 371 | DO k = nzb, nzt+1 |
---|
[1342] | 372 | d(k,j,i) = 0.0_wp |
---|
[1] | 373 | ENDDO |
---|
| 374 | ENDDO |
---|
| 375 | ENDDO |
---|
| 376 | ELSE |
---|
| 377 | !$OMP PARALLEL DO SCHEDULE( STATIC ) |
---|
[1111] | 378 | !$acc kernels present( d ) |
---|
[1003] | 379 | DO i = nxl, nxr |
---|
| 380 | DO j = nys, nyn |
---|
| 381 | DO k = nzb+1, nzt |
---|
[1342] | 382 | d(k,j,i) = 0.0_wp |
---|
[1] | 383 | ENDDO |
---|
| 384 | ENDDO |
---|
| 385 | ENDDO |
---|
[1111] | 386 | !$acc end kernels |
---|
[1] | 387 | ENDIF |
---|
| 388 | |
---|
[1342] | 389 | localsum = 0.0_wp |
---|
| 390 | threadsum = 0.0_wp |
---|
[1] | 391 | |
---|
| 392 | #if defined( __ibm ) |
---|
| 393 | !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) |
---|
| 394 | !$OMP DO SCHEDULE( STATIC ) |
---|
| 395 | DO i = nxl, nxr |
---|
| 396 | DO j = nys, nyn |
---|
| 397 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
[85] | 398 | d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & |
---|
| 399 | ( v(k,j+1,i) - v(k,j,i) ) * ddy + & |
---|
[673] | 400 | ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ) * ddt_3d & |
---|
| 401 | * d_weight_pres |
---|
[1] | 402 | ENDDO |
---|
| 403 | ! |
---|
| 404 | !-- Compute possible PE-sum of divergences for flow_statistics |
---|
| 405 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
| 406 | threadsum = threadsum + ABS( d(k,j,i) ) |
---|
| 407 | ENDDO |
---|
| 408 | |
---|
| 409 | ENDDO |
---|
| 410 | ENDDO |
---|
| 411 | |
---|
[1918] | 412 | IF ( intermediate_timestep_count == intermediate_timestep_count_max .OR. & |
---|
| 413 | intermediate_timestep_count == 0 ) THEN |
---|
| 414 | localsum = localsum + threadsum * dt_3d * weight_pres_l |
---|
[1908] | 415 | ENDIF |
---|
[1] | 416 | !$OMP END PARALLEL |
---|
| 417 | #else |
---|
| 418 | |
---|
[1111] | 419 | !$OMP PARALLEL PRIVATE (i,j,k) |
---|
| 420 | !$OMP DO SCHEDULE( STATIC ) |
---|
[1221] | 421 | !$acc kernels present( d, ddzw, rflags_s_inner, u, v, w ) |
---|
| 422 | !$acc loop collapse( 3 ) |
---|
[1111] | 423 | DO i = nxl, nxr |
---|
| 424 | DO j = nys, nyn |
---|
| 425 | DO k = 1, nzt |
---|
[1221] | 426 | d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & |
---|
| 427 | ( v(k,j+1,i) - v(k,j,i) ) * ddy + & |
---|
| 428 | ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) ) * ddt_3d & |
---|
| 429 | * d_weight_pres * rflags_s_inner(k,j,i) |
---|
[1] | 430 | ENDDO |
---|
| 431 | ENDDO |
---|
[1111] | 432 | ENDDO |
---|
| 433 | !$acc end kernels |
---|
| 434 | !$OMP END PARALLEL |
---|
[1] | 435 | |
---|
| 436 | ! |
---|
[1908] | 437 | !-- Compute possible PE-sum of divergences for flow_statistics. Carry out |
---|
| 438 | !-- computation only at last Runge-Kutta substep. |
---|
[1918] | 439 | IF ( intermediate_timestep_count == intermediate_timestep_count_max .OR. & |
---|
| 440 | intermediate_timestep_count == 0 ) THEN |
---|
[1908] | 441 | !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) |
---|
| 442 | !$OMP DO SCHEDULE( STATIC ) |
---|
| 443 | !$acc parallel loop collapse(3) present( d ) reduction(+:threadsum) |
---|
| 444 | DO i = nxl, nxr |
---|
| 445 | DO j = nys, nyn |
---|
| 446 | DO k = nzb+1, nzt |
---|
| 447 | threadsum = threadsum + ABS( d(k,j,i) ) |
---|
| 448 | ENDDO |
---|
[1] | 449 | ENDDO |
---|
| 450 | ENDDO |
---|
[1908] | 451 | !$acc end parallel loop |
---|
[1918] | 452 | localsum = localsum + threadsum * dt_3d * weight_pres_l |
---|
[1908] | 453 | !$OMP END PARALLEL |
---|
| 454 | ENDIF |
---|
[1] | 455 | #endif |
---|
| 456 | |
---|
| 457 | ! |
---|
| 458 | !-- For completeness, set the divergence sum of all statistic regions to those |
---|
| 459 | !-- of the total domain |
---|
[1918] | 460 | IF ( intermediate_timestep_count == intermediate_timestep_count_max .OR. & |
---|
| 461 | intermediate_timestep_count == 0 ) THEN |
---|
[1908] | 462 | sums_divold_l(0:statistic_regions) = localsum |
---|
[1918] | 463 | ENDIF |
---|
[1] | 464 | |
---|
| 465 | CALL cpu_log( log_point_s(1), 'divergence', 'stop' ) |
---|
| 466 | |
---|
| 467 | ! |
---|
| 468 | !-- Compute the pressure perturbation solving the Poisson equation |
---|
[1575] | 469 | IF ( psolver == 'poisfft' ) THEN |
---|
[1] | 470 | |
---|
| 471 | ! |
---|
| 472 | !-- Solve Poisson equation via FFT and solution of tridiagonal matrices |
---|
[1575] | 473 | CALL poisfft( d ) |
---|
[1212] | 474 | |
---|
[1] | 475 | ! |
---|
| 476 | !-- Store computed perturbation pressure and set boundary condition in |
---|
| 477 | !-- z-direction |
---|
| 478 | !$OMP PARALLEL DO |
---|
[1113] | 479 | !$acc kernels present( d, tend ) |
---|
[1] | 480 | DO i = nxl, nxr |
---|
| 481 | DO j = nys, nyn |
---|
| 482 | DO k = nzb+1, nzt |
---|
| 483 | tend(k,j,i) = d(k,j,i) |
---|
| 484 | ENDDO |
---|
| 485 | ENDDO |
---|
| 486 | ENDDO |
---|
[1113] | 487 | !$acc end kernels |
---|
[1] | 488 | |
---|
| 489 | ! |
---|
| 490 | !-- Bottom boundary: |
---|
| 491 | !-- This condition is only required for internal output. The pressure |
---|
| 492 | !-- gradient (dp(nzb+1)-dp(nzb))/dz is not used anywhere else. |
---|
| 493 | IF ( ibc_p_b == 1 ) THEN |
---|
| 494 | ! |
---|
| 495 | !-- Neumann (dp/dz = 0) |
---|
| 496 | !$OMP PARALLEL DO |
---|
[1113] | 497 | !$acc kernels present( nzb_s_inner, tend ) |
---|
[667] | 498 | DO i = nxlg, nxrg |
---|
| 499 | DO j = nysg, nyng |
---|
[1] | 500 | tend(nzb_s_inner(j,i),j,i) = tend(nzb_s_inner(j,i)+1,j,i) |
---|
| 501 | ENDDO |
---|
| 502 | ENDDO |
---|
[1113] | 503 | !$acc end kernels |
---|
[1] | 504 | |
---|
| 505 | ELSE |
---|
| 506 | ! |
---|
| 507 | !-- Dirichlet |
---|
| 508 | !$OMP PARALLEL DO |
---|
[1113] | 509 | !$acc kernels present( tend ) |
---|
[667] | 510 | DO i = nxlg, nxrg |
---|
| 511 | DO j = nysg, nyng |
---|
[1342] | 512 | tend(nzb_s_inner(j,i),j,i) = 0.0_wp |
---|
[1] | 513 | ENDDO |
---|
| 514 | ENDDO |
---|
[1113] | 515 | !$acc end kernels |
---|
[1] | 516 | |
---|
| 517 | ENDIF |
---|
| 518 | |
---|
| 519 | ! |
---|
| 520 | !-- Top boundary |
---|
| 521 | IF ( ibc_p_t == 1 ) THEN |
---|
| 522 | ! |
---|
| 523 | !-- Neumann |
---|
| 524 | !$OMP PARALLEL DO |
---|
[1113] | 525 | !$acc kernels present( tend ) |
---|
[667] | 526 | DO i = nxlg, nxrg |
---|
| 527 | DO j = nysg, nyng |
---|
[1] | 528 | tend(nzt+1,j,i) = tend(nzt,j,i) |
---|
| 529 | ENDDO |
---|
| 530 | ENDDO |
---|
[1113] | 531 | !$acc end kernels |
---|
[1] | 532 | |
---|
| 533 | ELSE |
---|
| 534 | ! |
---|
| 535 | !-- Dirichlet |
---|
| 536 | !$OMP PARALLEL DO |
---|
[1113] | 537 | !$acc kernels present( tend ) |
---|
[667] | 538 | DO i = nxlg, nxrg |
---|
| 539 | DO j = nysg, nyng |
---|
[1342] | 540 | tend(nzt+1,j,i) = 0.0_wp |
---|
[1] | 541 | ENDDO |
---|
| 542 | ENDDO |
---|
[1113] | 543 | !$acc end kernels |
---|
[1] | 544 | |
---|
| 545 | ENDIF |
---|
| 546 | |
---|
| 547 | ! |
---|
| 548 | !-- Exchange boundaries for p |
---|
[1113] | 549 | IF ( numprocs == 1 ) THEN ! workaround for single-core GPU runs |
---|
| 550 | on_device = .TRUE. ! to be removed after complete porting |
---|
| 551 | ELSE ! of ghost point exchange |
---|
| 552 | !$acc update host( tend ) |
---|
| 553 | ENDIF |
---|
[667] | 554 | CALL exchange_horiz( tend, nbgp ) |
---|
[1113] | 555 | IF ( numprocs == 1 ) THEN ! workaround for single-core GPU runs |
---|
| 556 | on_device = .FALSE. ! to be removed after complete porting |
---|
| 557 | ELSE ! of ghost point exchange |
---|
| 558 | !$acc update device( tend ) |
---|
| 559 | ENDIF |
---|
[1] | 560 | |
---|
| 561 | ELSEIF ( psolver == 'sor' ) THEN |
---|
| 562 | |
---|
| 563 | ! |
---|
| 564 | !-- Solve Poisson equation for perturbation pressure using SOR-Red/Black |
---|
| 565 | !-- scheme |
---|
[707] | 566 | CALL sor( d, ddzu_pres, ddzw, p_loc ) |
---|
| 567 | tend = p_loc |
---|
[1] | 568 | |
---|
[1575] | 569 | ELSEIF ( psolver(1:9) == 'multigrid' ) THEN |
---|
[1] | 570 | |
---|
| 571 | ! |
---|
| 572 | !-- Solve Poisson equation for perturbation pressure using Multigrid scheme, |
---|
[667] | 573 | !-- array tend is used to store the residuals, logical exchange_mg is used |
---|
| 574 | !-- to discern data exchange in multigrid ( 1 ghostpoint ) and normal grid |
---|
| 575 | !-- ( nbgp ghost points ). |
---|
[778] | 576 | |
---|
| 577 | !-- If the number of grid points of the gathered grid, which is collected |
---|
| 578 | !-- on PE0, is larger than the number of grid points of an PE, than array |
---|
| 579 | !-- tend will be enlarged. |
---|
| 580 | IF ( gathered_size > subdomain_size ) THEN |
---|
| 581 | DEALLOCATE( tend ) |
---|
| 582 | ALLOCATE( tend(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg( & |
---|
| 583 | mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,& |
---|
| 584 | nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg( & |
---|
| 585 | mg_switch_to_pe0_level)+1) ) |
---|
| 586 | ENDIF |
---|
| 587 | |
---|
[1575] | 588 | IF ( psolver == 'multigrid' ) THEN |
---|
| 589 | CALL poismg( tend ) |
---|
| 590 | ELSE |
---|
[1931] | 591 | CALL poismg_noopt( tend ) |
---|
[1575] | 592 | ENDIF |
---|
[707] | 593 | |
---|
[778] | 594 | IF ( gathered_size > subdomain_size ) THEN |
---|
| 595 | DEALLOCATE( tend ) |
---|
| 596 | ALLOCATE( tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) |
---|
| 597 | ENDIF |
---|
| 598 | |
---|
[1] | 599 | ! |
---|
| 600 | !-- Restore perturbation pressure on tend because this array is used |
---|
| 601 | !-- further below to correct the velocity fields |
---|
[707] | 602 | DO i = nxl-1, nxr+1 |
---|
| 603 | DO j = nys-1, nyn+1 |
---|
| 604 | DO k = nzb, nzt+1 |
---|
| 605 | tend(k,j,i) = p_loc(k,j,i) |
---|
| 606 | ENDDO |
---|
| 607 | ENDDO |
---|
| 608 | ENDDO |
---|
[667] | 609 | |
---|
[1] | 610 | ENDIF |
---|
| 611 | |
---|
| 612 | ! |
---|
[707] | 613 | !-- Store perturbation pressure on array p, used for pressure data output. |
---|
| 614 | !-- Ghost layers are added in the output routines (except sor-method: see below) |
---|
[1918] | 615 | IF ( intermediate_timestep_count <= 1 ) THEN |
---|
[707] | 616 | !$OMP PARALLEL PRIVATE (i,j,k) |
---|
| 617 | !$OMP DO |
---|
[1929] | 618 | !$acc kernels present( p, tend, weight_substep_l ) |
---|
[1257] | 619 | !$acc loop independent |
---|
[707] | 620 | DO i = nxl-1, nxr+1 |
---|
[1257] | 621 | !$acc loop independent |
---|
[707] | 622 | DO j = nys-1, nyn+1 |
---|
[1257] | 623 | !$acc loop independent |
---|
[707] | 624 | DO k = nzb, nzt+1 |
---|
| 625 | p(k,j,i) = tend(k,j,i) * & |
---|
[1929] | 626 | weight_substep_l |
---|
[673] | 627 | ENDDO |
---|
[1] | 628 | ENDDO |
---|
[707] | 629 | ENDDO |
---|
[1113] | 630 | !$acc end kernels |
---|
[707] | 631 | !$OMP END PARALLEL |
---|
| 632 | |
---|
[1762] | 633 | ELSEIF ( intermediate_timestep_count > 1 ) THEN |
---|
[707] | 634 | !$OMP PARALLEL PRIVATE (i,j,k) |
---|
| 635 | !$OMP DO |
---|
[1929] | 636 | !$acc kernels present( p, tend, weight_substep_l ) |
---|
[1257] | 637 | !$acc loop independent |
---|
[707] | 638 | DO i = nxl-1, nxr+1 |
---|
[1257] | 639 | !$acc loop independent |
---|
[707] | 640 | DO j = nys-1, nyn+1 |
---|
[1257] | 641 | !$acc loop independent |
---|
[707] | 642 | DO k = nzb, nzt+1 |
---|
| 643 | p(k,j,i) = p(k,j,i) + tend(k,j,i) * & |
---|
[1929] | 644 | weight_substep_l |
---|
[673] | 645 | ENDDO |
---|
| 646 | ENDDO |
---|
[707] | 647 | ENDDO |
---|
[1113] | 648 | !$acc end kernels |
---|
[707] | 649 | !$OMP END PARALLEL |
---|
| 650 | |
---|
| 651 | ENDIF |
---|
[673] | 652 | |
---|
[707] | 653 | ! |
---|
| 654 | !-- SOR-method needs ghost layers for the next timestep |
---|
| 655 | IF ( psolver == 'sor' ) CALL exchange_horiz( p, nbgp ) |
---|
[682] | 656 | |
---|
[1] | 657 | ! |
---|
| 658 | !-- Correction of the provisional velocities with the current perturbation |
---|
| 659 | !-- pressure just computed |
---|
[709] | 660 | IF ( conserve_volume_flow .AND. ( bc_lr_cyc .OR. bc_ns_cyc ) ) THEN |
---|
[1342] | 661 | volume_flow_l(1) = 0.0_wp |
---|
| 662 | volume_flow_l(2) = 0.0_wp |
---|
[1] | 663 | ENDIF |
---|
[707] | 664 | |
---|
[1] | 665 | !$OMP PARALLEL PRIVATE (i,j,k) |
---|
| 666 | !$OMP DO |
---|
[1918] | 667 | !$acc kernels present( ddzu, nzb_u_inner, nzb_v_inner, nzb_w_inner, tend, u, v, w ) |
---|
[1257] | 668 | !$acc loop independent |
---|
[673] | 669 | DO i = nxl, nxr |
---|
[1257] | 670 | !$acc loop independent |
---|
[1] | 671 | DO j = nys, nyn |
---|
[1257] | 672 | !$acc loop independent |
---|
[1113] | 673 | DO k = 1, nzt |
---|
| 674 | IF ( k > nzb_w_inner(j,i) ) THEN |
---|
| 675 | w(k,j,i) = w(k,j,i) - dt_3d * & |
---|
| 676 | ( tend(k+1,j,i) - tend(k,j,i) ) * ddzu(k+1) * & |
---|
[1918] | 677 | weight_pres_l |
---|
[1113] | 678 | ENDIF |
---|
[1] | 679 | ENDDO |
---|
[1257] | 680 | !$acc loop independent |
---|
[1113] | 681 | DO k = 1, nzt |
---|
| 682 | IF ( k > nzb_u_inner(j,i) ) THEN |
---|
| 683 | u(k,j,i) = u(k,j,i) - dt_3d * & |
---|
| 684 | ( tend(k,j,i) - tend(k,j,i-1) ) * ddx * & |
---|
[1918] | 685 | weight_pres_l |
---|
[1113] | 686 | ENDIF |
---|
[1] | 687 | ENDDO |
---|
[1257] | 688 | !$acc loop independent |
---|
[1113] | 689 | DO k = 1, nzt |
---|
| 690 | IF ( k > nzb_v_inner(j,i) ) THEN |
---|
| 691 | v(k,j,i) = v(k,j,i) - dt_3d * & |
---|
| 692 | ( tend(k,j,i) - tend(k,j-1,i) ) * ddy * & |
---|
[1918] | 693 | weight_pres_l |
---|
[1113] | 694 | ENDIF |
---|
[673] | 695 | ENDDO |
---|
[1] | 696 | |
---|
| 697 | ENDDO |
---|
| 698 | ENDDO |
---|
[1113] | 699 | !$acc end kernels |
---|
[1] | 700 | !$OMP END PARALLEL |
---|
[1113] | 701 | |
---|
| 702 | ! |
---|
| 703 | !-- Sum up the volume flow through the right and north boundary |
---|
| 704 | IF ( conserve_volume_flow .AND. bc_lr_cyc .AND. bc_ns_cyc .AND. & |
---|
| 705 | nxr == nx ) THEN |
---|
| 706 | |
---|
| 707 | !$OMP PARALLEL PRIVATE (j,k) |
---|
| 708 | !$OMP DO |
---|
| 709 | DO j = nys, nyn |
---|
| 710 | !$OMP CRITICAL |
---|
[1845] | 711 | DO k = nzb_u_inner(j,nx) + 1, nzt |
---|
[1113] | 712 | volume_flow_l(1) = volume_flow_l(1) + u(k,j,nx) * dzw(k) |
---|
| 713 | ENDDO |
---|
| 714 | !$OMP END CRITICAL |
---|
| 715 | ENDDO |
---|
| 716 | !$OMP END PARALLEL |
---|
| 717 | |
---|
| 718 | ENDIF |
---|
| 719 | |
---|
| 720 | IF ( conserve_volume_flow .AND. bc_ns_cyc .AND. bc_lr_cyc .AND. & |
---|
| 721 | nyn == ny ) THEN |
---|
| 722 | |
---|
| 723 | !$OMP PARALLEL PRIVATE (i,k) |
---|
| 724 | !$OMP DO |
---|
| 725 | DO i = nxl, nxr |
---|
| 726 | !$OMP CRITICAL |
---|
[1845] | 727 | DO k = nzb_v_inner(ny,i) + 1, nzt |
---|
[1113] | 728 | volume_flow_l(2) = volume_flow_l(2) + v(k,ny,i) * dzw(k) |
---|
| 729 | ENDDO |
---|
| 730 | !$OMP END CRITICAL |
---|
| 731 | ENDDO |
---|
| 732 | !$OMP END PARALLEL |
---|
| 733 | |
---|
| 734 | ENDIF |
---|
[673] | 735 | |
---|
[1] | 736 | ! |
---|
| 737 | !-- Conserve the volume flow |
---|
[707] | 738 | IF ( conserve_volume_flow .AND. ( bc_lr_cyc .AND. bc_ns_cyc ) ) THEN |
---|
[1] | 739 | |
---|
| 740 | #if defined( __parallel ) |
---|
[622] | 741 | IF ( collective_wait ) CALL MPI_BARRIER( comm2d, ierr ) |
---|
[1] | 742 | CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 2, MPI_REAL, & |
---|
| 743 | MPI_SUM, comm2d, ierr ) |
---|
| 744 | #else |
---|
| 745 | volume_flow = volume_flow_l |
---|
| 746 | #endif |
---|
| 747 | |
---|
[1799] | 748 | volume_flow_offset(1:2) = ( volume_flow_initial(1:2) - volume_flow(1:2) ) / & |
---|
| 749 | volume_flow_area(1:2) |
---|
[1] | 750 | |
---|
| 751 | !$OMP PARALLEL PRIVATE (i,j,k) |
---|
| 752 | !$OMP DO |
---|
| 753 | DO i = nxl, nxr |
---|
| 754 | DO j = nys, nyn |
---|
[667] | 755 | DO k = nzb_u_inner(j,i) + 1, nzt |
---|
| 756 | u(k,j,i) = u(k,j,i) + volume_flow_offset(1) |
---|
[719] | 757 | ENDDO |
---|
[1845] | 758 | DO k = nzb_v_inner(j,i) + 1, nzt |
---|
[667] | 759 | v(k,j,i) = v(k,j,i) + volume_flow_offset(2) |
---|
| 760 | ENDDO |
---|
[1] | 761 | ENDDO |
---|
| 762 | ENDDO |
---|
[667] | 763 | |
---|
[1] | 764 | !$OMP END PARALLEL |
---|
| 765 | |
---|
| 766 | ENDIF |
---|
| 767 | |
---|
| 768 | ! |
---|
| 769 | !-- Exchange of boundaries for the velocities |
---|
[1113] | 770 | IF ( numprocs == 1 ) THEN ! workaround for single-core GPU runs |
---|
| 771 | on_device = .TRUE. ! to be removed after complete porting |
---|
| 772 | ELSE ! of ghost point exchange |
---|
| 773 | !$acc update host( u, v, w ) |
---|
| 774 | ENDIF |
---|
[667] | 775 | CALL exchange_horiz( u, nbgp ) |
---|
| 776 | CALL exchange_horiz( v, nbgp ) |
---|
| 777 | CALL exchange_horiz( w, nbgp ) |
---|
[1113] | 778 | IF ( numprocs == 1 ) THEN ! workaround for single-core GPU runs |
---|
| 779 | on_device = .FALSE. ! to be removed after complete porting |
---|
| 780 | ELSE ! of ghost point exchange |
---|
| 781 | !$acc update device( u, v, w ) |
---|
| 782 | ENDIF |
---|
[1] | 783 | |
---|
| 784 | ! |
---|
| 785 | !-- Compute the divergence of the corrected velocity field, |
---|
[1908] | 786 | !-- A possible PE-sum is computed in flow_statistics. Carry out computation |
---|
| 787 | !-- only at last Runge-Kutta step. |
---|
[1918] | 788 | IF ( intermediate_timestep_count == intermediate_timestep_count_max .OR. & |
---|
| 789 | intermediate_timestep_count == 0 ) THEN |
---|
[1908] | 790 | CALL cpu_log( log_point_s(1), 'divergence', 'start' ) |
---|
| 791 | sums_divnew_l = 0.0_wp |
---|
[1] | 792 | |
---|
| 793 | ! |
---|
[1908] | 794 | !-- d must be reset to zero because it can contain nonzero values below the |
---|
| 795 | !-- topography |
---|
| 796 | IF ( topography /= 'flat' ) d = 0.0_wp |
---|
[1] | 797 | |
---|
[1908] | 798 | localsum = 0.0_wp |
---|
| 799 | threadsum = 0.0_wp |
---|
[1] | 800 | |
---|
[1908] | 801 | !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum) |
---|
| 802 | !$OMP DO SCHEDULE( STATIC ) |
---|
[1] | 803 | #if defined( __ibm ) |
---|
[1908] | 804 | DO i = nxl, nxr |
---|
| 805 | DO j = nys, nyn |
---|
| 806 | DO k = nzb_s_inner(j,i)+1, nzt |
---|
| 807 | d(k,j,i) = ( u(k,j,i+1) - u(k,j,i) ) * ddx + & |
---|
| 808 | ( v(k,j+1,i) - v(k,j,i) ) * ddy + & |
---|
| 809 | ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) |
---|
| 810 | ENDDO |
---|
| 811 | DO k = nzb+1, nzt |
---|
| 812 | threadsum = threadsum + ABS( d(k,j,i) ) |
---|
| 813 | ENDDO |
---|
[1] | 814 | ENDDO |
---|
| 815 | ENDDO |
---|
| 816 | #else |
---|
[1908] | 817 | !$acc kernels present( d, ddzw, rflags_s_inner, u, v, w ) |
---|
| 818 | !$acc loop collapse( 3 ) |
---|
| 819 | DO i = nxl, nxr |
---|
| 820 | DO j = nys, nyn |
---|
| 821 | DO k = 1, nzt |
---|
| 822 | d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * ddx + & |
---|
| 823 | ( v(k,j+1,i) - v(k,j,i) ) * ddy + & |
---|
| 824 | ( w(k,j,i) - w(k-1,j,i) ) * ddzw(k) & |
---|
| 825 | ) * rflags_s_inner(k,j,i) |
---|
| 826 | ENDDO |
---|
[1113] | 827 | ENDDO |
---|
| 828 | ENDDO |
---|
[1908] | 829 | !$acc end kernels |
---|
[1113] | 830 | ! |
---|
[1908] | 831 | !-- Compute possible PE-sum of divergences for flow_statistics |
---|
| 832 | !$acc parallel loop collapse(3) present( d ) reduction(+:threadsum) |
---|
| 833 | DO i = nxl, nxr |
---|
| 834 | DO j = nys, nyn |
---|
| 835 | DO k = nzb+1, nzt |
---|
| 836 | threadsum = threadsum + ABS( d(k,j,i) ) |
---|
| 837 | ENDDO |
---|
[1] | 838 | ENDDO |
---|
| 839 | ENDDO |
---|
[1908] | 840 | !$acc end parallel loop |
---|
[1] | 841 | #endif |
---|
[667] | 842 | |
---|
[1908] | 843 | localsum = localsum + threadsum |
---|
| 844 | !$OMP END PARALLEL |
---|
[1] | 845 | |
---|
| 846 | ! |
---|
[1908] | 847 | !-- For completeness, set the divergence sum of all statistic regions to those |
---|
| 848 | !-- of the total domain |
---|
| 849 | sums_divnew_l(0:statistic_regions) = localsum |
---|
[1] | 850 | |
---|
[1908] | 851 | CALL cpu_log( log_point_s(1), 'divergence', 'stop' ) |
---|
[1] | 852 | |
---|
[1908] | 853 | ENDIF |
---|
| 854 | |
---|
[1] | 855 | CALL cpu_log( log_point(8), 'pres', 'stop' ) |
---|
| 856 | |
---|
| 857 | |
---|
| 858 | END SUBROUTINE pres |
---|