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

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

last commit documented

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