source: palm/trunk/SOURCE/pres.f90 @ 1683

Last change on this file since 1683 was 1683, checked in by knoop, 9 years ago

last commit documented

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