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

Last change on this file since 1929 was 1929, checked in by suehring, 8 years ago

several bugfixes in particle model and serial mode

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