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

Last change on this file since 1863 was 1846, checked in by raasch, 8 years ago

last commit documented

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