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

Last change on this file since 3849 was 3849, checked in by knoop, 5 years ago

Bugfix: added proper OpenACC support to disturb_field

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