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

Last change on this file since 4346 was 4346, checked in by motisi, 4 years ago

Introduction of wall_flags_total_0, which currently sets bits based on static topography information used in wall_flags_static_0

  • Property svn:keywords set to Id
File size: 30.1 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 4346 2019-12-18 11:55:56Z motisi $
[4346]27! Introduction of wall_flags_total_0, which currently sets bits based on static
28! topography information used in wall_flags_static_0
29!
30! 4329 2019-12-10 15:46:36Z motisi
[4329]31! Renamed wall_flags_0 to wall_flags_static_0
32!
33! 4182 2019-08-22 15:20:23Z scharf
[4182]34! Corrected "Former revisions" section
35!
36! 4015 2019-06-05 13:25:35Z raasch
[4015]37! variable child_domain_nvn eliminated
38!
39! 3849 2019-04-01 16:35:16Z knoop
[3634]40! OpenACC port for SPEC
[1919]41!
[4182]42! Revision 1.1  1997/07/24 11:24:44  raasch
43! Initial revision
44!
45!
[1]46! Description:
47! ------------
[1682]48!> Compute the divergence of the provisional velocity field. Solve the Poisson
49!> equation for the perturbation pressure. Compute the final velocities using
50!> this perturbation pressure. Compute the remaining divergence.
[1]51!------------------------------------------------------------------------------!
[1682]52 SUBROUTINE pres
53 
[1]54
[1320]55    USE arrays_3d,                                                             &
[2037]56        ONLY:  d, ddzu, ddzu_pres, ddzw, dzw, p, p_loc, rho_air, rho_air_zw,   &
57               tend, u, v, w
[1320]58
59    USE control_parameters,                                                    &
[3182]60        ONLY:  bc_lr_cyc, bc_ns_cyc, bc_radiation_l, bc_radiation_n,           &
61               bc_radiation_r, bc_radiation_s, child_domain,                   &
62               conserve_volume_flow, coupling_mode,                            &
[3012]63               dt_3d, gathered_size, ibc_p_b, ibc_p_t,                         &
64               intermediate_timestep_count, intermediate_timestep_count_max,   &
[3347]65               mg_switch_to_pe0_level, nesting_offline,                        &
66               psolver, subdomain_size,                                        &
[3182]67               topography, volume_flow, volume_flow_area, volume_flow_initial
[1320]68
69    USE cpulog,                                                                &
70        ONLY:  cpu_log, log_point, log_point_s
71
72    USE grid_variables,                                                        &
73        ONLY:  ddx, ddy
74
75    USE indices,                                                               &
76        ONLY:  nbgp, ngp_2dh_outer, nx, nxl, nxlg, nxl_mg, nxr, nxrg, nxr_mg,  &
[2232]77               ny, nys, nysg, nys_mg, nyn, nyng, nyn_mg, nzb, nzt, nzt_mg,     &
[4346]78               wall_flags_total_0
[1320]79
80    USE kinds
81
[1]82    USE pegrid
[1933]83   
84    USE pmc_interface,                                                         &
85        ONLY:  nesting_mode 
[1]86
[1320]87    USE poisfft_mod,                                                           &
88        ONLY:  poisfft
89
[1575]90    USE poismg_mod
91
[2696]92    USE poismg_noopt_mod
93
[1320]94    USE statistics,                                                            &
95        ONLY:  statistic_regions, sums_divnew_l, sums_divold_l, weight_pres,   &
96               weight_substep
97
[2232]98    USE surface_mod,                                                           &
99        ONLY :  bc_h
100
[1]101    IMPLICIT NONE
102
[1682]103    INTEGER(iwp) ::  i              !<
104    INTEGER(iwp) ::  j              !<
105    INTEGER(iwp) ::  k              !<
[2232]106    INTEGER(iwp) ::  m              !<
[1]107
[1682]108    REAL(wp)     ::  ddt_3d         !<
[1918]109    REAL(wp)     ::  d_weight_pres  !<
[1682]110    REAL(wp)     ::  localsum       !<
111    REAL(wp)     ::  threadsum      !<
[1918]112    REAL(wp)     ::  weight_pres_l  !<
[1929]113    REAL(wp)     ::  weight_substep_l !<
[1]114
[1762]115    REAL(wp), DIMENSION(1:3)   ::  volume_flow_l       !<
116    REAL(wp), DIMENSION(1:3)   ::  volume_flow_offset  !<
[1682]117    REAL(wp), DIMENSION(1:nzt) ::  w_l                 !<
118    REAL(wp), DIMENSION(1:nzt) ::  w_l_l               !<
[1]119
120
121    CALL cpu_log( log_point(8), 'pres', 'start' )
122
[1918]123!
124!-- Calculate quantities to be used locally
[1342]125    ddt_3d = 1.0_wp / dt_3d
[1918]126    IF ( intermediate_timestep_count == 0 )  THEN
127!
128!--    If pres is called before initial time step
[1929]129       weight_pres_l    = 1.0_wp
130       d_weight_pres    = 1.0_wp
131       weight_substep_l = 1.0_wp
[1918]132    ELSE
[1929]133       weight_pres_l    = weight_pres(intermediate_timestep_count)
134       d_weight_pres    = 1.0_wp / weight_pres(intermediate_timestep_count)
135       weight_substep_l = weight_substep(intermediate_timestep_count)
[1918]136    ENDIF
[85]137
[1]138!
[707]139!-- Multigrid method expects array d to have one ghost layer.
140!--
[1575]141    IF ( psolver(1:9) == 'multigrid' )  THEN
[667]142     
[1]143       DEALLOCATE( d )
[667]144       ALLOCATE( d(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ) 
[707]145
146!
147!--    Since p is later used to hold the weighted average of the substeps, it
148!--    cannot be used in the iterative solver. Therefore, its initial value is
149!--    stored on p_loc, which is then iteratively advanced in every substep.
[1762]150       IF ( intermediate_timestep_count <= 1 )  THEN
[707]151          DO  i = nxl-1, nxr+1
152             DO  j = nys-1, nyn+1
153                DO  k = nzb, nzt+1
154                   p_loc(k,j,i) = p(k,j,i)
155                ENDDO
156             ENDDO
157          ENDDO
[667]158       ENDIF
159       
[1918]160    ELSEIF ( psolver == 'sor'  .AND.  intermediate_timestep_count <= 1 )  THEN
[707]161
162!
163!--    Since p is later used to hold the weighted average of the substeps, it
164!--    cannot be used in the iterative solver. Therefore, its initial value is
165!--    stored on p_loc, which is then iteratively advanced in every substep.
166       p_loc = p
167
[1]168    ENDIF
169
170!
[75]171!-- Conserve the volume flow at the outflow in case of non-cyclic lateral
172!-- boundary conditions
[106]173!-- WARNING: so far, this conservation does not work at the left/south
174!--          boundary if the topography at the inflow differs from that at the
175!--          outflow! For this case, volume_flow_area needs adjustment!
176!
177!-- Left/right
[3182]178    IF ( conserve_volume_flow  .AND.  ( bc_radiation_l .OR.                    &
179                                        bc_radiation_r ) )  THEN
[680]180
[1342]181       volume_flow(1)   = 0.0_wp
182       volume_flow_l(1) = 0.0_wp
[106]183
[3182]184       IF ( bc_radiation_l )  THEN
[106]185          i = 0
[3182]186       ELSEIF ( bc_radiation_r )  THEN
[106]187          i = nx+1
188       ENDIF
189
190       DO  j = nys, nyn
191!
192!--       Sum up the volume flow through the south/north boundary
[2232]193          DO  k = nzb+1, nzt
194             volume_flow_l(1) = volume_flow_l(1) + u(k,j,i) * dzw(k)           &
195                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]196                                        BTEST( wall_flags_total_0(k,j,i), 1 )  &
[2232]197                                            )
[106]198          ENDDO
199       ENDDO
200
201#if defined( __parallel )   
[680]202       IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
[106]203       CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 1, MPI_REAL, &
204                           MPI_SUM, comm1dy, ierr )   
205#else
206       volume_flow = volume_flow_l 
207#endif
[709]208       volume_flow_offset(1) = ( volume_flow_initial(1) - volume_flow(1) ) &
[106]209                               / volume_flow_area(1)
210
[667]211       DO  j = nysg, nyng
[2232]212          DO  k = nzb+1, nzt
213             u(k,j,i) = u(k,j,i) + volume_flow_offset(1)                       &
214                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]215                                        BTEST( wall_flags_total_0(k,j,i), 1 )  &
[2232]216                                            )
[106]217          ENDDO
218       ENDDO
219
220    ENDIF
221
222!
223!-- South/north
[3182]224    IF ( conserve_volume_flow  .AND.  ( bc_radiation_n .OR. bc_radiation_s ) )  THEN
[106]225
[1342]226       volume_flow(2)   = 0.0_wp
227       volume_flow_l(2) = 0.0_wp
[75]228
[3182]229       IF ( bc_radiation_s )  THEN
[106]230          j = 0
[3182]231       ELSEIF ( bc_radiation_n )  THEN
[75]232          j = ny+1
[106]233       ENDIF
234
235       DO  i = nxl, nxr
[75]236!
[106]237!--       Sum up the volume flow through the south/north boundary
[2232]238          DO  k = nzb+1, nzt
239             volume_flow_l(2) = volume_flow_l(2) + v(k,j,i) * dzw(k)           &
240                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]241                                        BTEST( wall_flags_total_0(k,j,i), 2 )  &
[2232]242                                            )
[75]243          ENDDO
[106]244       ENDDO
245
[75]246#if defined( __parallel )   
[680]247       IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
[75]248       CALL MPI_ALLREDUCE( volume_flow_l(2), volume_flow(2), 1, MPI_REAL, &
249                           MPI_SUM, comm1dx, ierr )   
250#else
251       volume_flow = volume_flow_l 
252#endif
253       volume_flow_offset(2) = ( volume_flow_initial(2) - volume_flow(2) )    &
[106]254                               / volume_flow_area(2)
[75]255
[667]256       DO  i = nxlg, nxrg
[2232]257          DO  k = nzb+1, nzt
258             v(k,j,i) = v(k,j,i) + volume_flow_offset(2)                       &
259                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]260                                        BTEST( wall_flags_total_0(k,j,i), 2 )  &
[2232]261                                            )
[75]262          ENDDO
[106]263       ENDDO
[75]264
265    ENDIF
266
[76]267!
[4015]268!-- Remove mean vertical velocity in case that Neumann conditions are used both at bottom and top
269!-- boundary. With Neumann conditions at both vertical boundaries, the solver cannot remove
270!-- mean vertical velocities. They should be removed, because incompressibility requires that
271!-- the vertical gradient of vertical velocity is zero. Since w=0 at the solid surface, it must be
272!-- zero everywhere.
273!-- This must not be done in case of a 3d-nesting child domain, because a mean vertical velocity
274!-- can physically exist in such a domain.
275!-- Also in case of offline nesting, mean vertical velocities may exist (and must not be removed),
276!-- caused by horizontal divergence/convergence of the large scale flow that is prescribed at the
277!-- side boundaries.
278!-- The removal cannot be done before the first initial time step because ngp_2dh_outer
[1918]279!-- is not yet known then.
[4015]280    IF ( ibc_p_b == 1  .AND.  ibc_p_t == 1  .AND.  .NOT. nesting_offline                           &
281         .AND. .NOT. ( child_domain .AND. nesting_mode /= 'vertical' )                             &
282         .AND. intermediate_timestep_count /= 0 )                                                  &
[1918]283    THEN
284       w_l = 0.0_wp;  w_l_l = 0.0_wp
285       DO  i = nxl, nxr
286          DO  j = nys, nyn
[2232]287             DO  k = nzb+1, nzt
[4015]288                w_l_l(k) = w_l_l(k) + w(k,j,i)                                                     &
[4346]289                                     * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 3 ) )
[76]290             ENDDO
291          ENDDO
[1918]292       ENDDO
[76]293#if defined( __parallel )   
[1918]294       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[4015]295       CALL MPI_ALLREDUCE( w_l_l(1), w_l(1), nzt, MPI_REAL, MPI_SUM, comm2d, ierr )
[76]296#else
[1918]297       w_l = w_l_l
[76]298#endif
[1918]299       DO  k = 1, nzt
300          w_l(k) = w_l(k) / ngp_2dh_outer(k,0)
301       ENDDO
302       DO  i = nxlg, nxrg
303          DO  j = nysg, nyng
[2232]304             DO  k = nzb+1, nzt
[4015]305                w(k,j,i) = w(k,j,i) - w_l(k)                                                       &
[4346]306                                     * MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_total_0(k,j,i), 3 ) )
[76]307             ENDDO
308          ENDDO
[1918]309       ENDDO
[76]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
[2232]317       !$OMP PARALLEL DO SCHEDULE( STATIC ) PRIVATE (i,j,k)
[1]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
[2232]326       !$OMP PARALLEL DO SCHEDULE( STATIC ) PRIVATE (i,j,k)
[3634]327       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
328       !$ACC PRESENT(d)
[1003]329       DO  i = nxl, nxr
330          DO  j = nys, nyn
331             DO  k = nzb+1, nzt
[1342]332                d(k,j,i) = 0.0_wp
[1]333             ENDDO
334          ENDDO
335       ENDDO
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
[2232]346          DO  k = nzb+1, nzt
[2037]347             d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * rho_air(k) * ddx +       &
348                          ( v(k,j+1,i) - v(k,j,i) ) * rho_air(k) * ddy +       &
349                          ( w(k,j,i)   * rho_air_zw(k) -                       &
350                            w(k-1,j,i) * rho_air_zw(k-1) ) * ddzw(k)           &
[2232]351                        ) * ddt_3d * d_weight_pres                             &
352                                   * MERGE( 1.0_wp, 0.0_wp,                    &
[4346]353                                      BTEST( wall_flags_total_0(k,j,i), 0 )    &
[2232]354                                          )
[1]355          ENDDO
356!
357!--       Compute possible PE-sum of divergences for flow_statistics
[2232]358          DO  k = nzb+1, nzt
359             threadsum = threadsum + ABS( d(k,j,i) )                           &
360                                   * MERGE( 1.0_wp, 0.0_wp,                    &
[4346]361                                      BTEST( wall_flags_total_0(k,j,i), 0 )    &
[2232]362                                          )
[1]363          ENDDO
364
365       ENDDO
366    ENDDO
367
[1918]368    IF ( intermediate_timestep_count == intermediate_timestep_count_max  .OR.  &
369         intermediate_timestep_count == 0 )  THEN
370       localsum = localsum + threadsum * dt_3d * weight_pres_l
[1908]371    ENDIF
[1]372    !$OMP END PARALLEL
373#else
374
[1111]375    !$OMP PARALLEL PRIVATE (i,j,k)
376    !$OMP DO SCHEDULE( STATIC )
[3634]377    !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
[4346]378    !$ACC PRESENT(u, v, w, rho_air, rho_air_zw, ddzw, wall_flags_total_0) &
[3634]379    !$ACC PRESENT(d)
[1111]380    DO  i = nxl, nxr
381       DO  j = nys, nyn
382          DO  k = 1, nzt
[2037]383             d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * rho_air(k) * ddx +       &
384                          ( v(k,j+1,i) - v(k,j,i) ) * rho_air(k) * ddy +       &
385                          ( w(k,j,i)   * rho_air_zw(k) -                       &
386                            w(k-1,j,i) * rho_air_zw(k-1) ) * ddzw(k)           &
[2232]387                        ) * ddt_3d * d_weight_pres                             &
388                                   * MERGE( 1.0_wp, 0.0_wp,                    &
[4346]389                                      BTEST( wall_flags_total_0(k,j,i), 0 )    &
[2232]390                                          )     
[1]391          ENDDO
392       ENDDO
[1111]393    ENDDO
394    !$OMP END PARALLEL
[1]395
396!
[1908]397!-- Compute possible PE-sum of divergences for flow_statistics. Carry out
398!-- computation only at last Runge-Kutta substep.
[1918]399    IF ( intermediate_timestep_count == intermediate_timestep_count_max  .OR.  &
400         intermediate_timestep_count == 0 )  THEN
[1908]401       !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum)
402       !$OMP DO SCHEDULE( STATIC )
[3634]403       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i,j,k) &
404       !$ACC REDUCTION(+:threadsum) COPY(threadsum) &
405       !$ACC PRESENT(d)
[1908]406       DO  i = nxl, nxr
407          DO  j = nys, nyn
408             DO  k = nzb+1, nzt
409                threadsum = threadsum + ABS( d(k,j,i) )
410             ENDDO
[1]411          ENDDO
412       ENDDO
[1918]413       localsum = localsum + threadsum * dt_3d * weight_pres_l
[1908]414       !$OMP END PARALLEL
415    ENDIF
[1]416#endif
417
418!
419!-- For completeness, set the divergence sum of all statistic regions to those
420!-- of the total domain
[1918]421    IF ( intermediate_timestep_count == intermediate_timestep_count_max  .OR.  &
422         intermediate_timestep_count == 0 )  THEN
[1908]423       sums_divold_l(0:statistic_regions) = localsum
[1918]424    ENDIF
[1]425
426    CALL cpu_log( log_point_s(1), 'divergence', 'stop' )
427
428!
429!-- Compute the pressure perturbation solving the Poisson equation
[1575]430    IF ( psolver == 'poisfft' )  THEN
[1]431
432!
433!--    Solve Poisson equation via FFT and solution of tridiagonal matrices
[1575]434       CALL poisfft( d )
[1212]435
[1]436!
437!--    Store computed perturbation pressure and set boundary condition in
438!--    z-direction
[2232]439       !$OMP PARALLEL DO PRIVATE (i,j,k)
[3634]440       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
441       !$ACC PRESENT(d, tend)
[1]442       DO  i = nxl, nxr
443          DO  j = nys, nyn
444             DO  k = nzb+1, nzt
445                tend(k,j,i) = d(k,j,i)
446             ENDDO
447          ENDDO
448       ENDDO
449
450!
451!--    Bottom boundary:
452!--    This condition is only required for internal output. The pressure
453!--    gradient (dp(nzb+1)-dp(nzb))/dz is not used anywhere else.
454       IF ( ibc_p_b == 1 )  THEN
455!
[2232]456!--       Neumann (dp/dz = 0). Using surfae data type, first for non-natural
457!--       surfaces, then for natural and urban surfaces
458!--       Upward facing
459          !$OMP PARALLEL DO PRIVATE( i, j, k )
[3634]460          !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
461          !$ACC PRESENT(bc_h, tend)
[2232]462          DO  m = 1, bc_h(0)%ns
[2298]463             i = bc_h(0)%i(m)
464             j = bc_h(0)%j(m)
465             k = bc_h(0)%k(m)
[2232]466             tend(k-1,j,i) = tend(k,j,i)
[1]467          ENDDO
[2232]468!
469!--       Downward facing
470          !$OMP PARALLEL DO PRIVATE( i, j, k )
[3634]471          !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
472          !$ACC PRESENT(bc_h, tend)
[2232]473          DO  m = 1, bc_h(1)%ns
[2298]474             i = bc_h(1)%i(m)
475             j = bc_h(1)%j(m)
476             k = bc_h(1)%k(m)
[2232]477             tend(k+1,j,i) = tend(k,j,i)
478          ENDDO
[1]479
480       ELSE
481!
[2232]482!--       Dirichlet. Using surface data type, first for non-natural
483!--       surfaces, then for natural and urban surfaces
484!--       Upward facing
485          !$OMP PARALLEL DO PRIVATE( i, j, k )
486          DO  m = 1, bc_h(0)%ns
[2298]487             i = bc_h(0)%i(m)
488             j = bc_h(0)%j(m)
489             k = bc_h(0)%k(m)
[2232]490             tend(k-1,j,i) = 0.0_wp
[1]491          ENDDO
[2232]492!
493!--       Downward facing
494          !$OMP PARALLEL DO PRIVATE( i, j, k )
495          DO  m = 1, bc_h(1)%ns
[2298]496             i = bc_h(1)%i(m)
497             j = bc_h(1)%j(m)
498             k = bc_h(1)%k(m)
[2232]499             tend(k+1,j,i) = 0.0_wp
500          ENDDO
[1]501
502       ENDIF
503
504!
505!--    Top boundary
506       IF ( ibc_p_t == 1 )  THEN
507!
508!--       Neumann
[2232]509          !$OMP PARALLEL DO PRIVATE (i,j,k)
[667]510          DO  i = nxlg, nxrg
511             DO  j = nysg, nyng
[1]512                tend(nzt+1,j,i) = tend(nzt,j,i)
513             ENDDO
514          ENDDO
515
516       ELSE
517!
518!--       Dirichlet
[2232]519          !$OMP PARALLEL DO PRIVATE (i,j,k)
[3634]520          !$ACC PARALLEL LOOP COLLAPSE(2) PRIVATE(i, j) &
521          !$ACC PRESENT(tend)
[667]522          DO  i = nxlg, nxrg
523             DO  j = nysg, nyng
[1342]524                tend(nzt+1,j,i) = 0.0_wp
[1]525             ENDDO
526          ENDDO
527
528       ENDIF
529
530!
531!--    Exchange boundaries for p
[667]532       CALL exchange_horiz( tend, nbgp )
[1]533     
534    ELSEIF ( psolver == 'sor' )  THEN
535
536!
537!--    Solve Poisson equation for perturbation pressure using SOR-Red/Black
538!--    scheme
[707]539       CALL sor( d, ddzu_pres, ddzw, p_loc )
540       tend = p_loc
[1]541
[1575]542    ELSEIF ( psolver(1:9) == 'multigrid' )  THEN
[1]543
544!
545!--    Solve Poisson equation for perturbation pressure using Multigrid scheme,
[2298]546!--    array tend is used to store the residuals.
[778]547
548!--    If the number of grid points of the gathered grid, which is collected
549!--    on PE0, is larger than the number of grid points of an PE, than array
550!--    tend will be enlarged.
551       IF ( gathered_size > subdomain_size )  THEN
552          DEALLOCATE( tend )
553          ALLOCATE( tend(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg(          &
554                    mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,&
555                    nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg(                   &
556                    mg_switch_to_pe0_level)+1) )
557       ENDIF
558
[1575]559       IF ( psolver == 'multigrid' )  THEN
560          CALL poismg( tend )
561       ELSE
[1931]562          CALL poismg_noopt( tend )
[1575]563       ENDIF
[707]564
[778]565       IF ( gathered_size > subdomain_size )  THEN
566          DEALLOCATE( tend )
567          ALLOCATE( tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
568       ENDIF
569
[1]570!
571!--    Restore perturbation pressure on tend because this array is used
572!--    further below to correct the velocity fields
[707]573       DO  i = nxl-1, nxr+1
574          DO  j = nys-1, nyn+1
575             DO  k = nzb, nzt+1
576                tend(k,j,i) = p_loc(k,j,i)
577             ENDDO
578          ENDDO
579       ENDDO
[667]580
[1]581    ENDIF
582
583!
[707]584!-- Store perturbation pressure on array p, used for pressure data output.
585!-- Ghost layers are added in the output routines (except sor-method: see below)
[1918]586    IF ( intermediate_timestep_count <= 1 )  THEN
[707]587       !$OMP PARALLEL PRIVATE (i,j,k)
588       !$OMP DO
[3634]589       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
590       !$ACC PRESENT(p, tend)
[707]591       DO  i = nxl-1, nxr+1
592          DO  j = nys-1, nyn+1
593             DO  k = nzb, nzt+1
594                p(k,j,i) = tend(k,j,i) * &
[1929]595                           weight_substep_l
[673]596             ENDDO
[1]597          ENDDO
[707]598       ENDDO
599       !$OMP END PARALLEL
600
[1762]601    ELSEIF ( intermediate_timestep_count > 1 )  THEN
[707]602       !$OMP PARALLEL PRIVATE (i,j,k)
603       !$OMP DO
[3634]604       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
605       !$ACC PRESENT(p, tend)
[707]606       DO  i = nxl-1, nxr+1
607          DO  j = nys-1, nyn+1
608             DO  k = nzb, nzt+1
609                p(k,j,i) = p(k,j,i) + tend(k,j,i) * &
[1929]610                           weight_substep_l
[673]611             ENDDO
612          ENDDO
[707]613       ENDDO
614       !$OMP END PARALLEL
615
616    ENDIF
[673]617       
[707]618!
619!-- SOR-method needs ghost layers for the next timestep
620    IF ( psolver == 'sor' )  CALL exchange_horiz( p, nbgp )
[682]621
[1]622!
623!-- Correction of the provisional velocities with the current perturbation
624!-- pressure just computed
[709]625    IF ( conserve_volume_flow  .AND.  ( bc_lr_cyc .OR. bc_ns_cyc ) )  THEN
[1342]626       volume_flow_l(1) = 0.0_wp
627       volume_flow_l(2) = 0.0_wp
[1]628    ENDIF
[3347]629!
630!-- Add pressure gradients to the velocity components. Note, the loops are
631!-- running over the entire model domain, even though, in case of non-cyclic
632!-- boundaries u- and v-component are not prognostic at i=0 and j=0,
633!-- respectiveley. However, in case of Dirichlet boundary conditions for the
634!-- velocities, zero-gradient conditions for the pressure are set, so that
635!-- no modification is imposed at the boundaries.
[1]636    !$OMP PARALLEL PRIVATE (i,j,k)
637    !$OMP DO
[3634]638    !$ACC PARALLEL LOOP COLLAPSE(2) PRIVATE(i, j, k) &
[4346]639    !$ACC PRESENT(u, v, w, tend, ddzu, wall_flags_total_0)
[673]640    DO  i = nxl, nxr   
[1]641       DO  j = nys, nyn
[2118]642
[2232]643          DO  k = nzb+1, nzt
644             w(k,j,i) = w(k,j,i) - dt_3d *                                     &
645                           ( tend(k+1,j,i) - tend(k,j,i) ) * ddzu(k+1)         &
646                                     * weight_pres_l                           &
647                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]648                                        BTEST( wall_flags_total_0(k,j,i), 3 )  &
[2232]649                                            )
[1]650          ENDDO
[2118]651
[2232]652          DO  k = nzb+1, nzt
653             u(k,j,i) = u(k,j,i) - dt_3d *                                     &
654                           ( tend(k,j,i) - tend(k,j,i-1) ) * ddx               &
655                                     * weight_pres_l                           &
656                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]657                                        BTEST( wall_flags_total_0(k,j,i), 1 )  &
[2232]658                                            )
[1]659          ENDDO
[2118]660
[2232]661          DO  k = nzb+1, nzt
662             v(k,j,i) = v(k,j,i) - dt_3d *                                     &
663                           ( tend(k,j,i) - tend(k,j-1,i) ) * ddy               &
664                                     * weight_pres_l                           &
665                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]666                                        BTEST( wall_flags_total_0(k,j,i), 2 )  &
[2232]667                                            )
[673]668          ENDDO                                                         
[1]669
670       ENDDO
671    ENDDO
672    !$OMP END PARALLEL
[1113]673
674!
[3012]675!-- The vertical velocity is not set to zero at nzt + 1 for nested domains
676!-- Instead it is set to the values of nzt (see routine vnest_boundary_conds
677!-- or pmci_interp_tril_t) BEFORE calling the pressure solver. To avoid jumps
678!-- while plotting profiles w at the top has to be set to the values in the
679!-- height nzt after above modifications. Hint: w level nzt+1 does not impact
680!-- results.
[3182]681    IF ( child_domain  .OR.  coupling_mode == 'vnested_fine' ) THEN
[3012]682       w(nzt+1,:,:) = w(nzt,:,:)
683    ENDIF
684
685!
[1113]686!-- Sum up the volume flow through the right and north boundary
[2232]687    IF ( conserve_volume_flow  .AND.  bc_lr_cyc  .AND.  bc_ns_cyc  .AND.       &
[1113]688         nxr == nx )  THEN
689
690       !$OMP PARALLEL PRIVATE (j,k)
691       !$OMP DO
692       DO  j = nys, nyn
693          !$OMP CRITICAL
[2232]694          DO  k = nzb+1, nzt
695             volume_flow_l(1) = volume_flow_l(1) + u(k,j,nxr) * dzw(k)         &
696                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]697                                        BTEST( wall_flags_total_0(k,j,nxr), 1 )&
[2232]698                                            )
[1113]699          ENDDO
700          !$OMP END CRITICAL
701       ENDDO
702       !$OMP END PARALLEL
703
704    ENDIF
705
[2232]706    IF ( conserve_volume_flow  .AND.  bc_ns_cyc  .AND.  bc_lr_cyc  .AND.       &
[1113]707         nyn == ny )  THEN
708
709       !$OMP PARALLEL PRIVATE (i,k)
710       !$OMP DO
711       DO  i = nxl, nxr
712          !$OMP CRITICAL
[2232]713          DO  k = nzb+1, nzt
714             volume_flow_l(2) = volume_flow_l(2) + v(k,nyn,i) * dzw(k)         &
715                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]716                                        BTEST( wall_flags_total_0(k,nyn,i), 2 )&
[2232]717                                            )
[1113]718           ENDDO
719          !$OMP END CRITICAL
720       ENDDO
721       !$OMP END PARALLEL
722
723    ENDIF
[673]724   
[1]725!
726!-- Conserve the volume flow
[707]727    IF ( conserve_volume_flow  .AND.  ( bc_lr_cyc  .AND.  bc_ns_cyc ) )  THEN
[1]728
729#if defined( __parallel )   
[622]730       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[1]731       CALL MPI_ALLREDUCE( volume_flow_l(1), volume_flow(1), 2, MPI_REAL, &
732                           MPI_SUM, comm2d, ierr ) 
733#else
734       volume_flow = volume_flow_l 
735#endif   
736
[1799]737       volume_flow_offset(1:2) = ( volume_flow_initial(1:2) - volume_flow(1:2) ) / &
738                            volume_flow_area(1:2)
[1]739
740       !$OMP PARALLEL PRIVATE (i,j,k)
741       !$OMP DO
742       DO  i = nxl, nxr
743          DO  j = nys, nyn
[2232]744             DO  k = nzb+1, nzt
745                u(k,j,i) = u(k,j,i) + volume_flow_offset(1)                    &
746                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]747                                        BTEST( wall_flags_total_0(k,j,i), 1 )  &
[2232]748                                            )
[719]749             ENDDO
[2232]750             DO  k = nzb+1, nzt
751                v(k,j,i) = v(k,j,i) + volume_flow_offset(2)                    &
752                                     * MERGE( 1.0_wp, 0.0_wp,                  &
[4346]753                                        BTEST( wall_flags_total_0(k,j,i), 2 )  &
[2232]754                                            )
[667]755             ENDDO
[1]756          ENDDO
757       ENDDO
[667]758
[1]759       !$OMP END PARALLEL
760
761    ENDIF
762
763!
764!-- Exchange of boundaries for the velocities
[667]765    CALL exchange_horiz( u, nbgp )
766    CALL exchange_horiz( v, nbgp )
767    CALL exchange_horiz( w, nbgp )
[1]768
769!
770!-- Compute the divergence of the corrected velocity field,
[1908]771!-- A possible PE-sum is computed in flow_statistics. Carry out computation
772!-- only at last Runge-Kutta step.
[1918]773    IF ( intermediate_timestep_count == intermediate_timestep_count_max  .OR.  &
774         intermediate_timestep_count == 0 )  THEN
[1908]775       CALL cpu_log( log_point_s(1), 'divergence', 'start' )
776       sums_divnew_l = 0.0_wp
[1]777
778!
[1908]779!--    d must be reset to zero because it can contain nonzero values below the
780!--    topography
781       IF ( topography /= 'flat' )  d = 0.0_wp
[1]782
[1908]783       localsum  = 0.0_wp
784       threadsum = 0.0_wp
[1]785
[1908]786       !$OMP PARALLEL PRIVATE (i,j,k) FIRSTPRIVATE(threadsum) REDUCTION(+:localsum)
[2073]787#if defined( __ibm )
[1908]788       !$OMP DO SCHEDULE( STATIC )
789       DO  i = nxl, nxr
790          DO  j = nys, nyn
[2232]791             DO  k = nzb+1, nzt
792             d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * rho_air(k) * ddx +       &
793                          ( v(k,j+1,i) - v(k,j,i) ) * rho_air(k) * ddy +       &
794                          ( w(k,j,i)   * rho_air_zw(k) -                       &
795                            w(k-1,j,i) * rho_air_zw(k-1) ) * ddzw(k)           &
796                        ) * MERGE( 1.0_wp, 0.0_wp,                             &
[4346]797                             BTEST( wall_flags_total_0(k,j,i), 0 )             &
[2232]798                                 )
[1908]799             ENDDO
800             DO  k = nzb+1, nzt
801                threadsum = threadsum + ABS( d(k,j,i) )
802             ENDDO
[1]803          ENDDO
804       ENDDO
805#else
[2073]806       !$OMP DO SCHEDULE( STATIC )
[3634]807       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
[4346]808       !$ACC PRESENT(u, v, w, rho_air, rho_air_zw, ddzw, wall_flags_total_0) &
[3634]809       !$ACC PRESENT(d)
[1908]810       DO  i = nxl, nxr
811          DO  j = nys, nyn
[2232]812             DO  k = nzb+1, nzt
[2037]813                d(k,j,i) = ( ( u(k,j,i+1) - u(k,j,i) ) * rho_air(k) * ddx +    &
814                             ( v(k,j+1,i) - v(k,j,i) ) * rho_air(k) * ddy +    &
815                             ( w(k,j,i)   * rho_air_zw(k) -                    &
816                               w(k-1,j,i) * rho_air_zw(k-1) ) * ddzw(k)        &
[2232]817                           ) * MERGE( 1.0_wp, 0.0_wp,                          &
[4346]818                             BTEST( wall_flags_total_0(k,j,i), 0 )             &
[2232]819                                    )
[1908]820             ENDDO
[1113]821          ENDDO
822       ENDDO
823!
[1908]824!--    Compute possible PE-sum of divergences for flow_statistics
[2073]825       !$OMP DO SCHEDULE( STATIC )
[3634]826       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i, j, k) &
827       !$ACC REDUCTION(+:threadsum) COPY(threadsum) &
828       !$ACC PRESENT(d)
[1908]829       DO  i = nxl, nxr
830          DO  j = nys, nyn
831             DO  k = nzb+1, nzt
832                threadsum = threadsum + ABS( d(k,j,i) )
833             ENDDO
[1]834          ENDDO
835       ENDDO
836#endif
[667]837
[1908]838       localsum = localsum + threadsum
839       !$OMP END PARALLEL
[1]840
841!
[1908]842!--    For completeness, set the divergence sum of all statistic regions to those
843!--    of the total domain
844       sums_divnew_l(0:statistic_regions) = localsum
[1]845
[1908]846       CALL cpu_log( log_point_s(1), 'divergence', 'stop' )
[1]847
[1908]848    ENDIF
849
[1]850    CALL cpu_log( log_point(8), 'pres', 'stop' )
851
852 END SUBROUTINE pres
Note: See TracBrowser for help on using the repository browser.