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

Last change on this file since 1792 was 1763, checked in by hellstea, 8 years ago

last commit documented

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