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

Last change on this file since 1927 was 1919, checked in by raasch, 8 years ago

last commit documented

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