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

Last change on this file since 4104 was 4015, checked in by raasch, 5 years ago

all reals changed to double precision in order to work with 32-bit working precision, otherwise calculated time intervals would mostly give zero; variable child_domain_nvn eliminated

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