source: palm/trunk/SOURCE/poismg_fast.f90 @ 1575

Last change on this file since 1575 was 1575, checked in by raasch, 9 years ago

optimized multigrid method installed, new parameter seed_follows_topography for particle release, small adjustment in subjob for HLRN

  • Property svn:keywords set to Id
File size: 113.0 KB
Line 
1 MODULE poismg_mod
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
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!
17! Copyright 1997-2014 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Attention: Loop unrolling and cache optimization in SOR-Red/Black method
21!            still does not give the expected speedup! Further work required.
22!
23! Current revisions:
24! -----------------
25! Initial revision.
26! Routine re-written and optimised based on poismg.
27!
28! Following optimisations have been made:
29! - vectorisation (for Intel-CPUs) of the red-black algorithm by resorting
30!   array elements with even and odd indices
31! - explicit boundary conditions for building walls removed (solver is
32!   running through the buildings
33! - reduced data transfer in case of ghost point exchange, because only
34!   "red" or "black" data points need to be exchanged. This is not applied
35!   for coarser grid levels, since for then the transfer time is latency bound
36!
37! Former revisions:
38! -----------------
39! $Id: poismg_fast.f90 1575 2015-03-27 09:56:27Z raasch $
40!
41!
42!
43! Description:
44! ------------
45! Solves the Poisson equation for the perturbation pressure with a multigrid
46! V- or W-Cycle scheme.
47!
48! This multigrid method was originally developed for PALM by Joerg Uhlenbrock,
49! September 2000 - July 2001. It has been optimised for speed by Klaus
50! Ketelsen in November 2014.
51!------------------------------------------------------------------------------!
52
53    USE cpulog,                                                                &
54        ONLY:  cpu_log, log_point_s
55    USE kinds
56    USE MPI
57    USE pegrid
58
59    PRIVATE
60
61    INTEGER, SAVE                             ::  ind_even_odd    !: border index between even and odd k index
62    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  ::  even_odd_level  !: stores ind_even_odd for all MG levels
63
64    REAL(wp), DIMENSION(:,:), SAVE, ALLOCATABLE ::  f1_mg_b, f2_mg_b, f3_mg_b  !: blocked version of f1_mg ...
65
66    INTERFACE poismg_fast
67       MODULE PROCEDURE poismg_fast
68    END INTERFACE poismg_fast
69
70    INTERFACE sort_k_to_even_odd_blocks
71       MODULE PROCEDURE sort_k_to_even_odd_blocks
72       MODULE PROCEDURE sort_k_to_even_odd_blocks_int
73       MODULE PROCEDURE sort_k_to_even_odd_blocks_1d
74    END INTERFACE sort_k_to_even_odd_blocks
75
76    PUBLIC poismg_fast
77
78 CONTAINS
79
80    SUBROUTINE poismg_fast( r )
81
82       USE arrays_3d,                                                          &
83           ONLY:  d, p_loc
84
85       USE control_parameters,                                                 &
86           ONLY:  gathered_size, grid_level, grid_level_count,                 &
87                  maximum_grid_level, message_string, mgcycles, mg_cycles,     &
88                  mg_switch_to_pe0_level, residual_limit, subdomain_size
89
90       USE cpulog,                                                             &
91           ONLY:  cpu_log, log_point_s
92
93       USE indices,                                                            &
94           ONLY:  nxl, nxlg, nxl_mg, nxr, nxrg, nxr_mg, nys, nysg, nys_mg, nyn,&
95                  nyng, nyn_mg, nzb, nzt, nzt_mg
96
97       IMPLICIT NONE
98
99       REAL(wp) ::  maxerror          !:
100       REAL(wp) ::  maximum_mgcycles  !:
101       REAL(wp) ::  residual_norm     !:
102
103       REAL(wp), DIMENSION(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ::  r  !:
104
105       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  p3  !:
106
107
108       CALL cpu_log( log_point_s(29), 'poismg_fast', 'start' )
109!
110!--    Initialize arrays and variables used in this subroutine
111
112!--    If the number of grid points of the gathered grid, which is collected
113!--    on PE0, is larger than the number of grid points of an PE, than array
114!--    p3 will be enlarged.
115       IF ( gathered_size > subdomain_size )  THEN
116          ALLOCATE( p3(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg(            &
117                    mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,&
118                    nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg(                   &
119                    mg_switch_to_pe0_level)+1) )
120       ELSE
121          ALLOCATE ( p3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
122       ENDIF
123
124       p3 = 0.0_wp
125
126 
127!
128!--    Ghost boundaries have to be added to divergence array.
129!--    Exchange routine needs to know the grid level!
130       grid_level = maximum_grid_level
131       CALL exchange_horiz( d, 1)
132       d(nzb,:,:) = d(nzb+1,:,:)
133
134!
135!--    Initiation of the multigrid scheme. Does n cycles until the
136!--    residual is smaller than the given limit. The accuracy of the solution
137!--    of the poisson equation will increase with the number of cycles.
138!--    If the number of cycles is preset by the user, this number will be
139!--    carried out regardless of the accuracy.
140       grid_level_count =  0
141       mgcycles         =  0
142       IF ( mg_cycles == -1 )  THEN
143          maximum_mgcycles = 0
144          residual_norm    = 1.0_wp
145       ELSE
146          maximum_mgcycles = mg_cycles
147          residual_norm    = 0.0_wp
148       ENDIF
149
150!
151!--    Initial settings for sorting k-dimension from sequential order (alternate
152!--    even/odd) into blocks of even and odd or vice versa
153       CALL init_even_odd_blocks
154
155!
156!--    Sort input arrays in even/odd blocks along k-dimension
157       CALL sort_k_to_even_odd_blocks( d, grid_level )
158       CALL sort_k_to_even_odd_blocks( p_loc, grid_level )
159
160!
161!--    The complete multigrid cycles are running in block mode, i.e. over
162!--    seperate data blocks of even and odd indices
163       DO WHILE ( residual_norm > residual_limit  .OR. &
164                  mgcycles < maximum_mgcycles )
165 
166          CALL next_mg_level_fast( d, p_loc, p3, r)
167
168!
169!--       Calculate the residual if the user has not preset the number of
170!--       cycles to be performed
171          IF ( maximum_mgcycles == 0 )  THEN
172             CALL resid_fast( d, p_loc, r )
173             maxerror = SUM( r(nzb+1:nzt,nys:nyn,nxl:nxr)**2 )
174
175#if defined( __parallel )
176             IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
177                CALL MPI_ALLREDUCE( maxerror, residual_norm, 1, MPI_REAL,      &
178                                    MPI_SUM, comm2d, ierr)
179#else
180                residual_norm = maxerror
181#endif
182             residual_norm = SQRT( residual_norm )
183          ENDIF
184
185          mgcycles = mgcycles + 1
186
187!
188!--       If the user has not limited the number of cycles, stop the run in case
189!--       of insufficient convergence
190          IF ( mgcycles > 1000  .AND.  mg_cycles == -1 )  THEN
191             message_string = 'no sufficient convergence within 1000 cycles'
192             CALL message( 'poismg_fast', 'PA0283', 1, 2, 0, 6, 0 )
193          ENDIF
194
195       ENDDO
196
197       DEALLOCATE( p3 )
198!
199!--    Result has to be sorted back from even/odd blocks to sequential order
200       CALL sort_k_to_sequential( p_loc )
201!
202!--    Unset the grid level. Variable is used to determine the MPI datatypes for
203!--    ghost point exchange
204       grid_level = 0
205
206       CALL cpu_log( log_point_s(29), 'poismg_fast', 'stop' )
207
208    END SUBROUTINE poismg_fast
209
210
211
212    SUBROUTINE resid_fast( f_mg, p_mg, r )
213
214!------------------------------------------------------------------------------!
215! Description:
216! ------------
217! Computes the residual of the perturbation pressure.
218!------------------------------------------------------------------------------!
219
220       USE arrays_3d,                                                          &
221           ONLY:  f1_mg, f2_mg, f3_mg
222
223       USE control_parameters,                                                 &
224           ONLY:  bc_lr_cyc, bc_ns_cyc, grid_level, ibc_p_b, ibc_p_t, inflow_l,&
225                  inflow_n, inflow_r, inflow_s, masking_method, outflow_l,     &
226                  outflow_n, outflow_r, outflow_s, topography
227
228       USE grid_variables,                                                     &
229           ONLY:  ddx2_mg, ddy2_mg
230
231       USE indices,                                                            &
232           ONLY:  flags, wall_flags_1, wall_flags_2, wall_flags_3,             &
233                  wall_flags_4, wall_flags_5, wall_flags_6, wall_flags_7,      &
234                  wall_flags_8, wall_flags_9, wall_flags_10, nxl_mg, nxr_mg,   &
235                  nys_mg, nyn_mg, nzb, nzt_mg
236
237       IMPLICIT NONE
238
239       INTEGER(iwp) ::  i
240       INTEGER(iwp) ::  j
241       INTEGER(iwp) ::  k
242       INTEGER(iwp) ::  l
243       INTEGER(iwp) :: km1      !:
244       INTEGER(iwp) :: kp1      !:
245
246       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
247                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
248                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f_mg  !:
249       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
250                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
251                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  p_mg  !:
252       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
253                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
254                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  r     !:
255
256!
257!--    Calculate the residual
258       l = grid_level
259
260       CALL cpu_log( log_point_s(53), 'resid', 'start' )
261       IF ( topography == 'flat'  .OR.  masking_method )  THEN
262          !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
263          !$OMP DO
264          DO  i = nxl_mg(l), nxr_mg(l)
265             DO  j = nys_mg(l), nyn_mg(l)
266                !DIR$ IVDEP
267                DO k = ind_even_odd+1, nzt_mg(l)
268                   km1 = k-ind_even_odd-1
269                   kp1 = k-ind_even_odd
270                   r(k,j,i) = f_mg(k,j,i)                                 &
271                      - ddx2_mg(l) *                                      &
272                      ( p_mg(k,j,i+1) +  p_mg(k,j,i-1)  )                 &
273                      - ddy2_mg(l) *                                      &
274                      ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )                  &
275                      - f2_mg_b(k,l) * p_mg(kp1,j,i)                      &
276                      - f3_mg_b(k,l) * p_mg(km1,j,i)                      &
277                      + f1_mg_b(k,l) * p_mg(k,j,i)
278                ENDDO
279                !DIR$ IVDEP
280                DO k = nzb+1, ind_even_odd
281                   km1 = k+ind_even_odd
282                   kp1 = k+ind_even_odd+1
283                   r(k,j,i) = f_mg(k,j,i)                                 &
284                      - ddx2_mg(l) *                                      &
285                      ( p_mg(k,j,i+1) +  p_mg(k,j,i-1)  )                 &
286                      - ddy2_mg(l) *                                      &
287                      ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )                  &
288                      - f2_mg_b(k,l) * p_mg(kp1,j,i)                      &
289                      - f3_mg_b(k,l) * p_mg(km1,j,i)                      &
290                      + f1_mg_b(k,l) * p_mg(k,j,i)
291                ENDDO
292             ENDDO
293          ENDDO
294          !$OMP END PARALLEL
295       ELSE
296!
297!--       Choose flag array of this level
298          SELECT CASE ( l )
299             CASE ( 1 )
300                flags => wall_flags_1
301             CASE ( 2 )
302                flags => wall_flags_2
303             CASE ( 3 )
304                flags => wall_flags_3
305             CASE ( 4 )
306                flags => wall_flags_4
307             CASE ( 5 )
308                flags => wall_flags_5
309             CASE ( 6 )
310                flags => wall_flags_6
311             CASE ( 7 )
312                flags => wall_flags_7
313             CASE ( 8 )
314                flags => wall_flags_8
315             CASE ( 9 )
316                flags => wall_flags_9
317             CASE ( 10 )
318                flags => wall_flags_10
319          END SELECT
320
321          !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
322          !$OMP DO
323          DO  i = nxl_mg(l), nxr_mg(l)
324             DO  j = nys_mg(l), nyn_mg(l)
325
326                !DIR$ IVDEP
327                DO  k = ind_even_odd+1, nzt_mg(l)
328                   km1 = k-ind_even_odd-1
329                   kp1 = k-ind_even_odd
330
331                   r(k,j,i) = f_mg(k,j,i)                                         &
332                    - ddx2_mg(l) *                                                &
333                    ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5 ))   &
334                    + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4 )) ) &
335                    - ddy2_mg(l)                                                  &
336                    * ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3 )) &
337                    + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2 )) ) &
338                    - f2_mg(k,l) * p_mg(kp1,j,i)                                  &
339                    - f3_mg(k,l) *                                                &
340                      MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0 ))   &
341                    + f1_mg(k,l) * p_mg(k,j,i)
342                ENDDO
343
344                !DIR$ IVDEP
345                DO  k = nzb+1, ind_even_odd
346                   km1 = k+ind_even_odd
347                   kp1 = k+ind_even_odd+1
348
349                   r(k,j,i) = f_mg(k,j,i)                                         &
350                    - ddx2_mg(l) *                                                &
351                    ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5 ))   &
352                    + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4 )) ) &
353                    - ddy2_mg(l)                                                  &
354                    * ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3 )) &
355                    + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2 )) ) &
356                    - f2_mg(k,l) * p_mg(kp1,j,i)                                  &
357                    - f3_mg(k,l) *                                                &
358                      MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0 ))   &
359                    + f1_mg(k,l) * p_mg(k,j,i)
360                ENDDO
361!
362!--             The residual within topography should be zero
363                WHERE( BTEST(flags(nzb+1:nzt_mg(l),j,i), 6) )
364                   r(nzb+1:nzt_mg(l),j,i) = 0.0
365                END WHERE
366
367             ENDDO
368          ENDDO
369          !$OMP END PARALLEL
370
371       ENDIF
372
373!
374!--    Horizontal boundary conditions
375       CALL exchange_horiz( r, 1)
376
377       IF ( .NOT. bc_lr_cyc )  THEN
378          IF ( inflow_l .OR. outflow_l )  r(:,:,nxl_mg(l)-1) = r(:,:,nxl_mg(l))
379          IF ( inflow_r .OR. outflow_r )  r(:,:,nxr_mg(l)+1) = r(:,:,nxr_mg(l))
380       ENDIF
381
382       IF ( .NOT. bc_ns_cyc )  THEN
383          IF ( inflow_n .OR. outflow_n )  r(:,nyn_mg(l)+1,:) = r(:,nyn_mg(l),:)
384          IF ( inflow_s .OR. outflow_s )  r(:,nys_mg(l)-1,:) = r(:,nys_mg(l),:)
385       ENDIF
386
387!
388!--    Boundary conditions at bottom and top of the domain.outflow_l, outflow_n,
389!--    These points are not handled by the above loop. Points may be within
390!--    buildings, but that doesn't matter.
391       IF ( ibc_p_b == 1 )  THEN
392          r(nzb,:,: ) = r(nzb+1,:,:)
393       ELSE
394          r(nzb,:,: ) = 0.0_wp
395       ENDIF
396
397       IF ( ibc_p_t == 1 )  THEN
398          r(nzt_mg(l)+1,:,: ) = r(nzt_mg(l),:,:)
399       ELSE
400          r(nzt_mg(l)+1,:,: ) = 0.0_wp
401       ENDIF
402
403       CALL cpu_log( log_point_s(53), 'resid', 'stop' )
404
405    END SUBROUTINE resid_fast
406
407
408
409    SUBROUTINE restrict_fast( f_mg, r )
410
411!------------------------------------------------------------------------------!
412! Description:
413! ------------
414! Interpolates the residual on the next coarser grid with "full weighting"
415! scheme
416!------------------------------------------------------------------------------!
417
418       USE control_parameters,                                                 &
419           ONLY:  bc_lr_cyc, bc_ns_cyc, grid_level, ibc_p_b, ibc_p_t, inflow_l,&
420                  inflow_n, inflow_r, inflow_s, masking_method, outflow_l,     &
421                  outflow_n, outflow_r, outflow_s, topography
422
423       USE indices,                                                            &
424           ONLY:  flags, wall_flags_1, wall_flags_2, wall_flags_3,             &
425                  wall_flags_4, wall_flags_5, wall_flags_6, wall_flags_7,      &
426                  wall_flags_8, wall_flags_9, wall_flags_10, nxl_mg, nxr_mg,   &
427                  nys_mg, nyn_mg, nzb, nzt_mg
428
429       IMPLICIT NONE
430
431       INTEGER(iwp) ::  i    !:
432       INTEGER(iwp) ::  ic   !:
433       INTEGER(iwp) ::  j    !:
434       INTEGER(iwp) ::  jc   !:
435       INTEGER(iwp) ::  k    !:
436       INTEGER(iwp) ::  kc   !:
437       INTEGER(iwp) ::  l    !:
438       INTEGER(iwp) ::  km1  !:
439       INTEGER(iwp) ::  kp1  !:
440
441       REAL(wp) ::  rkjim    !:
442       REAL(wp) ::  rkjip    !:
443       REAL(wp) ::  rkjmi    !:
444       REAL(wp) ::  rkjmim   !:
445       REAL(wp) ::  rkjmip   !:
446       REAL(wp) ::  rkjpi    !:
447       REAL(wp) ::  rkjpim   !:
448       REAL(wp) ::  rkjpip   !:
449       REAL(wp) ::  rkmji    !:
450       REAL(wp) ::  rkmjim   !:
451       REAL(wp) ::  rkmjip   !:
452       REAL(wp) ::  rkmjmi   !:
453       REAL(wp) ::  rkmjmim  !:
454       REAL(wp) ::  rkmjmip  !:
455       REAL(wp) ::  rkmjpi   !:
456       REAL(wp) ::  rkmjpim  !:
457       REAL(wp) ::  rkmjpip  !:
458
459       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
460                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
461                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f_mg  !:
462
463       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level+1)+1,                         &
464                           nys_mg(grid_level+1)-1:nyn_mg(grid_level+1)+1,      &
465                           nxl_mg(grid_level+1)-1:nxr_mg(grid_level+1)+1) ::  r !:
466
467!
468!--    Interpolate the residual
469       l = grid_level
470
471       CALL cpu_log( log_point_s(54), 'restrict', 'start' )
472
473       IF ( topography == 'flat'  .OR.  masking_method )  THEN
474!
475!--       No wall treatment
476          !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,kc,km1,kp1)
477          DO  ic = nxl_mg(l), nxr_mg(l)
478             i = 2*ic
479             !$OMP DO SCHEDULE( STATIC )
480             DO  jc = nys_mg(l), nyn_mg(l)
481!
482!--             Calculation for the first point along k
483                j  = 2*jc
484                k   = ind_even_odd+1
485                kp1 = k-ind_even_odd
486                kc  = k-ind_even_odd
487                f_mg(kc,jc,ic) =  1.0_wp / 64.0_wp * (                         &
488                                  8.0_wp * r(k,j,i)                            &
489                               +  4.0_wp * ( r(k,j,i-1)     + r(k,j,i+1)     + &
490                                             r(k,j+1,i)     + r(k,j-1,i)     ) &
491                               +  2.0_wp * ( r(k,j-1,i-1)   + r(k,j+1,i-1)   + &
492                                             r(k,j-1,i+1)   + r(k,j+1,i+1)   ) &
493                               + 16.0_wp * r(k,j,i)                            &
494                               +  4.0_wp * r(kp1,j,i)                          &
495                               +  2.0_wp * ( r(kp1,j,i-1)   + r(kp1,j,i+1)   + &
496                                             r(kp1,j+1,i)   + r(kp1,j-1,i)   ) &
497                               +           ( r(kp1,j-1,i-1) + r(kp1,j+1,i-1) + &
498                                             r(kp1,j-1,i+1) + r(kp1,j+1,i+1) ) &
499                                                      )
500!
501!--             Calculation for the other points along k
502                !DIR$ IVDEP
503                DO k = ind_even_odd+2, nzt_mg(l+1)    ! Fine grid at this point
504                   km1 = k-ind_even_odd-1
505                   kp1 = k-ind_even_odd
506                   kc  = k-ind_even_odd               ! Coarse grid index
507
508                   f_mg(kc,jc,ic) = 1.0_wp / 64.0_wp * (                      &
509                                  8.0_wp * r(k,j,i)                            &
510                                + 4.0_wp * ( r(k,j,i-1)     + r(k,j,i+1)     + &
511                                             r(k,j+1,i)     + r(k,j-1,i)     ) &
512                                + 2.0_wp * ( r(k,j-1,i-1)   + r(k,j+1,i-1)   + &
513                                             r(k,j-1,i+1)   + r(k,j+1,i+1)   ) &
514                                + 4.0_wp * r(km1,j,i)                          &
515                                + 2.0_wp * ( r(km1,j,i-1)   + r(km1,j,i+1)   + &
516                                             r(km1,j+1,i)   + r(km1,j-1,i)   ) &
517                                +          ( r(km1,j-1,i-1) + r(km1,j+1,i-1) + &
518                                             r(km1,j-1,i+1) + r(km1,j+1,i+1) ) &
519                                + 4.0_wp * r(kp1,j,i)                          &
520                                + 2.0_wp * ( r(kp1,j,i-1)   + r(kp1,j,i+1)   + &
521                                             r(kp1,j+1,i)   + r(kp1,j-1,i)   ) &
522                                +          ( r(kp1,j-1,i-1) + r(kp1,j+1,i-1) + &
523                                             r(kp1,j-1,i+1) + r(kp1,j+1,i+1) ) &
524                                           )
525                ENDDO
526             ENDDO
527             !$OMP ENDDO nowait
528          ENDDO
529          !$OMP END PARALLEL
530
531       ELSE
532!
533!--       Choose flag array of the upper level
534          SELECT CASE ( l+1 )
535             CASE ( 1 )
536                flags => wall_flags_1
537             CASE ( 2 )
538                flags => wall_flags_2
539             CASE ( 3 )
540                flags => wall_flags_3
541             CASE ( 4 )
542                flags => wall_flags_4
543             CASE ( 5 )
544                flags => wall_flags_5
545             CASE ( 6 )
546                flags => wall_flags_6
547             CASE ( 7 )
548                flags => wall_flags_7
549             CASE ( 8 )
550                flags => wall_flags_8
551             CASE ( 9 )
552                flags => wall_flags_9
553             CASE ( 10 )
554                flags => wall_flags_10
555          END SELECT
556   
557          !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,kc,km1,kp1,rkjim,rkjip,rkjpi,    &
558          !$OMP                   rkjmim,rkjpim,rkjmip,rkjpip,rkmji,rkmjim,    &
559          !$OMP                   rkmjip,rkmjpi,rkmjmi,rkmjmim,rkmjpim,        &
560          !$OMP                   rkmjmip,rkmjpip)
561          !$OMP DO
562          DO  ic = nxl_mg(l), nxr_mg(l)
563             i = 2*ic
564             DO  jc = nys_mg(l), nyn_mg(l)
565                j = 2*jc
566                !DIR$ IVDEP
567                DO k = ind_even_odd+1, nzt_mg(l+1)      !fine grid at this point
568
569                    km1 = k-ind_even_odd-1
570                    kp1 = k-ind_even_odd
571                    kc  = k-ind_even_odd                !Coarse grid index
572!
573!--                 Use implicit Neumann BCs if the respective gridpoint is
574!--                 inside the building
575                    rkjim   = MERGE(r(k,j,i),r(k,j,i-1),BTEST( flags(k,j,i-1), 6))
576                    rkjip   = MERGE(r(k,j,i),r(k,j,i+1),BTEST( flags(k,j,i+1), 6))
577                    rkjpi   = MERGE(r(k,j,i),r(k,j+1,i),BTEST( flags(k,j+1,i), 6))
578                    rkjmi   = MERGE(r(k,j,i),r(k,j-1,i),BTEST( flags(k,j-1,i), 6))
579
580                    rkjmim  = MERGE(r(k,j,i),r(k,j-1,i-1),BTEST( flags(k,j-1,i-1), 6))
581                    rkjpim  = MERGE(r(k,j,i),r(k,j+1,i-1),BTEST( flags(k,j+1,i-1), 6))
582                    rkjmip  = MERGE(r(k,j,i),r(k,j-1,i+1),BTEST( flags(k,j-1,i+1), 6))
583                    rkjpip  = MERGE(r(k,j,i),r(k,j+1,i+1),BTEST( flags(k,j+1,i+1), 6))
584
585                    rkmji   = MERGE(r(k,j,i),r(km1,j,i)  ,BTEST( flags(km1,j,i)  , 6))
586                    rkmjim  = MERGE(r(k,j,i),r(km1,j,i-1),BTEST( flags(km1,j,i-1), 6))
587                    rkmjip  = MERGE(r(k,j,i),r(km1,j,i+1),BTEST( flags(km1,j,i+1), 6))
588                    rkmjpi  = MERGE(r(k,j,i),r(km1,j+1,i),BTEST( flags(km1,j+1,i), 6))
589                    rkmjmi  = MERGE(r(k,j,i),r(km1,j-1,i),BTEST( flags(km1,j-1,i), 6))
590
591                    rkmjmim = MERGE(r(k,j,i),r(km1,j-1,i-1),BTEST( flags(km1,j-1,i-1), 6))
592                    rkmjpim = MERGE(r(k,j,i),r(km1,j+1,i-1),BTEST( flags(km1,j+1,i-1), 6))
593                    rkmjmip = MERGE(r(k,j,i),r(km1,j-1,i+1),BTEST( flags(km1,j-1,i+1), 6))
594                    rkmjpip = MERGE(r(k,j,i),r(km1,j+1,i+1),BTEST( flags(km1,j+1,i+1), 6))
595
596                    f_mg(kc,jc,ic) = 1.0_wp / 64.0_wp * (                      &
597                           8.0_wp * r(k,j,i)                                   &
598                         + 4.0_wp * ( rkjim   + rkjip   +  rkjpi   + rkjmi   ) &
599                         + 2.0_wp * ( rkjmim  + rkjpim  +  rkjmip  + rkjpip  ) &
600                         + 4.0_wp * rkmji                                      &
601                         + 2.0_wp * ( rkmjim  + rkmjip  +  rkmjpi  + rkmjmi  ) &
602                         +          ( rkmjmim + rkmjpim +  rkmjmip + rkmjpip ) &
603                         + 4.0_wp * r(kp1,j,i)                                 &
604                         + 2.0_wp * ( r(kp1,j,i-1)   + r(kp1,j,i+1)   +        &
605                                      r(kp1,j+1,i)   + r(kp1,j-1,i)   )        &
606                         +          ( r(kp1,j-1,i-1) + r(kp1,j+1,i-1) +        &
607                                      r(kp1,j-1,i+1) + r(kp1,j+1,i+1) )        &
608                                                        )
609                ENDDO
610             ENDDO
611          ENDDO
612          !$OMP END PARALLEL
613
614       ENDIF
615
616!
617!--    Horizontal boundary conditions
618       CALL exchange_horiz( f_mg, 1)
619
620       IF ( .NOT. bc_lr_cyc )  THEN
621          IF (inflow_l .OR. outflow_l)  f_mg(:,:,nxl_mg(l)-1) = f_mg(:,:,nxl_mg(l))
622          IF (inflow_r .OR. outflow_r)  f_mg(:,:,nxr_mg(l)+1) = f_mg(:,:,nxr_mg(l))
623       ENDIF
624
625       IF ( .NOT. bc_ns_cyc )  THEN
626          IF (inflow_n .OR. outflow_n)  f_mg(:,nyn_mg(l)+1,:) = f_mg(:,nyn_mg(l),:)
627          IF (inflow_s .OR. outflow_s)  f_mg(:,nys_mg(l)-1,:) = f_mg(:,nys_mg(l),:)
628       ENDIF
629
630!
631!--    Boundary conditions at bottom and top of the domain.
632!--    These points are not handled by the above loop. Points may be within
633!--    buildings, but that doesn't matter.
634       IF ( ibc_p_b == 1 )  THEN
635          f_mg(nzb,:,: ) = f_mg(nzb+1,:,:)
636       ELSE
637          f_mg(nzb,:,: ) = 0.0_wp
638       ENDIF
639
640       IF ( ibc_p_t == 1 )  THEN
641          f_mg(nzt_mg(l)+1,:,: ) = f_mg(nzt_mg(l),:,:)
642       ELSE
643          f_mg(nzt_mg(l)+1,:,: ) = 0.0_wp
644       ENDIF
645
646       CALL cpu_log( log_point_s(54), 'restrict', 'stop' )
647!
648!--    Why do we need a sorting here?
649       CALL sort_k_to_even_odd_blocks( f_mg , l)
650
651    END SUBROUTINE restrict_fast
652
653
654
655    SUBROUTINE prolong_fast( p, temp )
656
657!------------------------------------------------------------------------------!
658! Description:
659! ------------
660! Interpolates the correction of the perturbation pressure
661! to the next finer grid.
662!------------------------------------------------------------------------------!
663
664       USE control_parameters,                                                 &
665           ONLY:  bc_lr_cyc, bc_ns_cyc, grid_level, ibc_p_b, ibc_p_t, inflow_l,&
666                  inflow_n, inflow_r, inflow_s, outflow_l, outflow_n,          &
667                  outflow_r, outflow_s
668
669       USE indices,                                                            &
670           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
671
672       IMPLICIT NONE
673
674       INTEGER(iwp) ::  i   !:
675       INTEGER(iwp) ::  j   !:
676       INTEGER(iwp) ::  k   !:
677       INTEGER(iwp) ::  l   !:
678       INTEGER(iwp) ::  kp1 !:
679       INTEGER(iwp) ::  ke  !: Index for prolog even
680       INTEGER(iwp) ::  ko  !: Index for prolog odd
681
682       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
683                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
684                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1 ) ::  p  !:
685
686       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
687                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
688                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  temp  !:
689
690
691       CALL cpu_log( log_point_s(55), 'prolong', 'start' )
692
693!
694!--    First, store elements of the coarser grid on the next finer grid
695       l = grid_level
696       ind_even_odd = even_odd_level(grid_level-1)
697
698       !$OMP PARALLEL PRIVATE (i,j,k,kp1,ke,ko)
699       !$OMP DO
700       DO  i = nxl_mg(l-1), nxr_mg(l-1)
701          DO  j = nys_mg(l-1), nyn_mg(l-1)
702
703             !DIR$ IVDEP
704             DO k = ind_even_odd+1, nzt_mg(l-1)
705                kp1 = k - ind_even_odd
706                ke  = 2 * ( k-ind_even_odd - 1 ) + 1
707                ko  = 2 * k - 1
708!
709!--             Points of the coarse grid are directly stored on the next finer
710!--             grid
711                temp(ko,2*j,2*i)   = p(k,j,i)
712!
713!--             Points between two coarse-grid points
714                temp(ko,2*j,2*i+1) = 0.5_wp * ( p(k,j,i) + p(k,j,i+1) )
715                temp(ko,2*j+1,2*i) = 0.5_wp * ( p(k,j,i) + p(k,j+1,i) )
716                temp(ke,2*j,2*i)   = 0.5_wp * ( p(k,j,i) + p(kp1,j,i) )
717!
718!--             Points in the center of the planes stretched by four points
719!--             of the coarse grid cube
720                temp(ko,2*j+1,2*i+1) = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
721                                                   p(k,j+1,i) + p(k,j+1,i+1) )
722                temp(ke,2*j,2*i+1)   = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
723                                                   p(kp1,j,i) + p(kp1,j,i+1) )
724                temp(ke,2*j+1,2*i)   = 0.25_wp * ( p(k,j,i)   + p(k,j+1,i) +   &
725                                                   p(kp1,j,i) + p(kp1,j+1,i) )
726!
727!--             Points in the middle of coarse grid cube
728                temp(ke,2*j+1,2*i+1) = 0.125_wp *                              &
729                                               ( p(k,j,i)     + p(k,j,i+1)   + &
730                                                 p(k,j+1,i)   + p(k,j+1,i+1) + &
731                                                 p(kp1,j,i)   + p(kp1,j,i+1) + &
732                                                 p(kp1,j+1,i) + p(kp1,j+1,i+1) )
733             ENDDO
734
735             !DIR$ IVDEP
736             DO k = nzb+1, ind_even_odd
737                kp1 = k + ind_even_odd + 1
738                ke  = 2 * k
739                ko  = 2 * ( k + ind_even_odd )
740!
741!--             Points of the coarse grid are directly stored on the next finer
742!--             grid
743                temp(ko,2*j,2*i)   = p(k,j,i)
744!
745!--             Points between two coarse-grid points
746                temp(ko,2*j,2*i+1) = 0.5_wp * ( p(k,j,i) + p(k,j,i+1) )
747                temp(ko,2*j+1,2*i) = 0.5_wp * ( p(k,j,i) + p(k,j+1,i) )
748                temp(ke,2*j,2*i)   = 0.5_wp * ( p(k,j,i) + p(kp1,j,i) )
749!
750!--             Points in the center of the planes stretched by four points
751!--             of the coarse grid cube
752                temp(ko,2*j+1,2*i+1) = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
753                                                   p(k,j+1,i) + p(k,j+1,i+1) )
754                temp(ke,2*j,2*i+1)   = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
755                                                   p(kp1,j,i) + p(kp1,j,i+1) )
756                temp(ke,2*j+1,2*i)   = 0.25_wp * ( p(k,j,i)   + p(k,j+1,i) +   &
757                                                   p(kp1,j,i) + p(kp1,j+1,i) )
758!
759!--             Points in the middle of coarse grid cube
760                temp(ke,2*j+1,2*i+1) = 0.125_wp *                              &
761                                               ( p(k,j,i)     + p(k,j,i+1)   + &
762                                                 p(k,j+1,i)   + p(k,j+1,i+1) + &
763                                                 p(kp1,j,i)   + p(kp1,j,i+1) + &
764                                                 p(kp1,j+1,i) + p(kp1,j+1,i+1) )
765             ENDDO
766
767          ENDDO
768       ENDDO
769       !$OMP END PARALLEL
770
771       ind_even_odd = even_odd_level(grid_level)
772!
773!--    Horizontal boundary conditions
774       CALL exchange_horiz( temp, 1)
775
776       IF ( .NOT. bc_lr_cyc )  THEN
777          IF (inflow_l .OR. outflow_l)  temp(:,:,nxl_mg(l)-1) = temp(:,:,nxl_mg(l))
778          IF (inflow_r .OR. outflow_r)  temp(:,:,nxr_mg(l)+1) = temp(:,:,nxr_mg(l))
779       ENDIF
780
781       IF ( .NOT. bc_ns_cyc )  THEN
782          IF (inflow_n .OR. outflow_n)  temp(:,nyn_mg(l)+1,:) = temp(:,nyn_mg(l),:)
783          IF (inflow_s .OR. outflow_s)  temp(:,nys_mg(l)-1,:) = temp(:,nys_mg(l),:)
784       ENDIF
785
786!
787!--    Bottom and top boundary conditions
788       IF ( ibc_p_b == 1 )  THEN
789          temp(nzb,:,: ) = temp(ind_even_odd+1,:,:)
790       ELSE
791          temp(nzb,:,: ) = 0.0_wp
792       ENDIF
793
794       IF ( ibc_p_t == 1 )  THEN
795          temp(nzt_mg(l)+1,:,: ) = temp(ind_even_odd,:,:)
796       ELSE
797          temp(nzt_mg(l)+1,:,: ) = 0.0_wp
798       ENDIF
799
800       CALL cpu_log( log_point_s(55), 'prolong', 'stop' )
801
802    END SUBROUTINE prolong_fast
803
804
805
806    SUBROUTINE redblack_fast( f_mg, p_mg )
807
808!------------------------------------------------------------------------------!
809! Description:
810! ------------
811! Relaxation method for the multigrid scheme. A Gauss-Seidel iteration with
812! 3D-Red-Black decomposition (GS-RB) is used.
813!------------------------------------------------------------------------------!
814
815       USE arrays_3d,                                                          &
816           ONLY:  f1_mg, f2_mg, f3_mg
817
818       USE control_parameters,                                                 &
819           ONLY:  bc_lr_cyc, bc_ns_cyc, grid_level, ibc_p_b, ibc_p_t, inflow_l,&
820                  inflow_n, inflow_r, inflow_s, masking_method, ngsrb,         &
821                  outflow_l, outflow_n, outflow_r, outflow_s, topography
822
823       USE grid_variables,                                                     &
824           ONLY:  ddx2_mg, ddy2_mg
825
826       USE indices,                                                            &
827           ONLY:  flags, wall_flags_1, wall_flags_2, wall_flags_3,             &
828                  wall_flags_4, wall_flags_5, wall_flags_6, wall_flags_7,      &
829                  wall_flags_8, wall_flags_9, wall_flags_10, nxl_mg, nxr_mg,   &
830                  nys_mg, nyn_mg, nzb, nzt_mg
831
832       IMPLICIT NONE
833
834       INTEGER(iwp) :: color    !:
835       INTEGER(iwp) :: i        !:
836       INTEGER(iwp) :: ic       !:
837       INTEGER(iwp) :: j        !:
838       INTEGER(iwp) :: jc       !:
839       INTEGER(iwp) :: jj       !:
840       INTEGER(iwp) :: k        !:
841       INTEGER(iwp) :: l        !:
842       INTEGER(iwp) :: n        !:
843       INTEGER(iwp) :: km1      !:
844       INTEGER(iwp) :: kp1      !:
845
846       LOGICAL :: unroll        !:
847
848       REAL(wp) ::  wall_left   !:
849       REAL(wp) ::  wall_north  !:
850       REAL(wp) ::  wall_right  !:
851       REAL(wp) ::  wall_south  !:
852       REAL(wp) ::  wall_total  !:
853       REAL(wp) ::  wall_top    !:
854
855       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
856                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
857                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f_mg  !:
858       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
859                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
860                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  p_mg  !:
861
862       l = grid_level
863
864!
865!--    Choose flag array of this level
866       IF ( topography /= 'flat'  .AND.  .NOT. masking_method )  THEN
867
868          SELECT CASE ( l )
869             CASE ( 1 )
870                flags => wall_flags_1
871             CASE ( 2 )
872                flags => wall_flags_2
873             CASE ( 3 )
874                flags => wall_flags_3
875             CASE ( 4 )
876                flags => wall_flags_4
877             CASE ( 5 )
878                flags => wall_flags_5
879             CASE ( 6 )
880                flags => wall_flags_6
881             CASE ( 7 )
882                flags => wall_flags_7
883             CASE ( 8 )
884                flags => wall_flags_8
885             CASE ( 9 )
886                flags => wall_flags_9
887             CASE ( 10 )
888                flags => wall_flags_10
889          END SELECT
890
891       ENDIF
892
893       unroll = ( MOD( nyn_mg(l)-nys_mg(l)+1, 4 ) == 0  .AND. &
894                  MOD( nxr_mg(l)-nxl_mg(l)+1, 2 ) == 0 )
895
896       DO  n = 1, ngsrb
897       
898          DO  color = 1, 2
899!
900!--          No wall treatment in case of flat topography
901             IF ( topography == 'flat'  .OR.  masking_method )  THEN
902
903                IF ( .NOT. unroll )  THEN
904
905                   CALL cpu_log( log_point_s(36), 'redblack_no_unroll_f', 'start' )
906!
907!--                Without unrolling of loops, no cache optimization
908                   !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
909                   !$OMP DO
910                   DO  i = nxl_mg(l), nxr_mg(l), 2
911                      DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
912                         !DIR$ IVDEP
913                         DO  k = ind_even_odd+1, nzt_mg(l)
914                            km1 = k-ind_even_odd-1
915                            kp1 = k-ind_even_odd
916                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
917                               ddx2_mg(l) *                                    &
918                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
919                               + ddy2_mg(l) *                                  &
920                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )              &
921                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
922                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
923                               - f_mg(k,j,i)               )
924                         ENDDO
925                      ENDDO
926                   ENDDO
927   
928                   !$OMP DO
929                   DO  i = nxl_mg(l)+1, nxr_mg(l), 2
930                      DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
931                          !DIR$ IVDEP
932                          DO  k = ind_even_odd+1, nzt_mg(l)
933                            km1 = k-ind_even_odd-1
934                            kp1 = k-ind_even_odd
935                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
936                               ddx2_mg(l) *                                    &
937                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
938                               + ddy2_mg(l) *                                  &
939                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
940                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
941                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
942                               - f_mg(k,j,i)               )
943                         ENDDO
944                      ENDDO
945                   ENDDO
946 
947                   !$OMP DO
948                   DO  i = nxl_mg(l), nxr_mg(l), 2
949                      DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
950                         !DIR$ IVDEP
951                         DO  k = nzb+1, ind_even_odd
952                            km1 = k+ind_even_odd
953                            kp1 = k+ind_even_odd+1
954                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
955                               ddx2_mg(l) *                                    &
956                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
957                               + ddy2_mg(l) *                                  &
958                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
959                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
960                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
961                               - f_mg(k,j,i)               )
962                         ENDDO
963                      ENDDO
964                   ENDDO
965
966                   !$OMP DO
967                   DO  i = nxl_mg(l)+1, nxr_mg(l), 2
968                      DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
969                         !DIR$ IVDEP
970                         DO  k = nzb+1, ind_even_odd
971                            km1 = k+ind_even_odd
972                            kp1 = k+ind_even_odd+1
973                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
974                               ddx2_mg(l) *                                    &
975                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
976                               + ddy2_mg(l) *                                  &
977                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
978                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
979                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
980                               - f_mg(k,j,i)               )
981                         ENDDO
982                      ENDDO
983                   ENDDO
984                   !$OMP END PARALLEL
985
986                   CALL cpu_log( log_point_s(36), 'redblack_no_unroll_f', 'stop' )
987
988                ELSE
989!
990!--                Loop unrolling along y, only one i loop for better cache use
991                   CALL cpu_log( log_point_s(38), 'redblack_unroll_f', 'start' )
992
993                   !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,km1,kp1,jj)
994                   !$OMP DO
995                   DO  ic = nxl_mg(l), nxr_mg(l), 2
996                      DO  jc = nys_mg(l), nyn_mg(l), 4
997                         i  = ic
998                         jj = jc+2-color
999                         !DIR$ IVDEP
1000                         DO  k = ind_even_odd+1, nzt_mg(l)
1001                            km1 = k-ind_even_odd-1
1002                            kp1 = k-ind_even_odd
1003                            j   = jj
1004                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1005                               ddx2_mg(l) *                                    &
1006                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
1007                               + ddy2_mg(l) *                                  &
1008                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
1009                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1010                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1011                               - f_mg(k,j,i)               )
1012                            j = jj+2
1013                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1014                               ddx2_mg(l) *                                    &
1015                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
1016                               + ddy2_mg(l) *                                  &
1017                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
1018                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1019                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1020                               - f_mg(k,j,i)               )
1021                         ENDDO
1022
1023                         i  = ic+1
1024                         jj = jc+color-1
1025                         !DIR$ IVDEP
1026                         DO  k = ind_even_odd+1, nzt_mg(l)
1027                            km1 = k-ind_even_odd-1
1028                            kp1 = k-ind_even_odd
1029                            j   = jj
1030                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1031                               ddx2_mg(l) *                                    &
1032                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
1033                               + ddy2_mg(l) *                                  &
1034                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
1035                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1036                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1037                               - f_mg(k,j,i)               )
1038                            j = jj+2
1039                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1040                               ddx2_mg(l) *                                    &
1041                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
1042                               + ddy2_mg(l) *                                  &
1043                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )              &
1044                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1045                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1046                               - f_mg(k,j,i)               )
1047                         ENDDO
1048
1049                         i  = ic
1050                         jj = jc+color-1
1051                         !DIR$ IVDEP
1052                         DO  k = nzb+1, ind_even_odd
1053                            km1 = k+ind_even_odd
1054                            kp1 = k+ind_even_odd+1
1055                            j   = jj
1056                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1057                               ddx2_mg(l) *                                    &
1058                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
1059                               + ddy2_mg(l) *                                  &
1060                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )              &
1061                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1062                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1063                               - f_mg(k,j,i)               )
1064                            j = jj+2
1065                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1066                               ddx2_mg(l) *                                    &
1067                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
1068                               + ddy2_mg(l) *                                  &
1069                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
1070                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1071                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1072                               - f_mg(k,j,i)               )
1073                         ENDDO
1074
1075                         i  = ic+1
1076                         jj = jc+2-color
1077                         !DIR$ IVDEP
1078                         DO  k = nzb+1, ind_even_odd
1079                            km1 = k+ind_even_odd
1080                            kp1 = k+ind_even_odd+1
1081                            j   = jj
1082                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1083                               ddx2_mg(l) *                                    &
1084                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1)  )              &
1085                               + ddy2_mg(l) *                                  &
1086                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
1087                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1088                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1089                               - f_mg(k,j,i)               )
1090                            j = jj+2
1091                            p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (            &
1092                               ddx2_mg(l) *                                    &
1093                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
1094                               + ddy2_mg(l) *                                  &
1095                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )              &
1096                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
1097                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
1098                               - f_mg(k,j,i)               )
1099                         ENDDO
1100
1101                      ENDDO
1102                   ENDDO
1103                   !$OMP END PARALLEL
1104
1105                   CALL cpu_log( log_point_s(38), 'redblack_unroll_f', 'stop' )
1106
1107                ENDIF
1108
1109             ELSE
1110
1111                IF ( .NOT. unroll )  THEN
1112
1113                   CALL cpu_log( log_point_s(36), 'redblack_no_unroll', 'start' )
1114!
1115!--                Without unrolling of loops, no cache optimization /
1116!--                vectorization. Therefore, the non-vectorized IBITS function
1117!--                is still used.
1118                   !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
1119                   !$OMP DO
1120                   DO  i = nxl_mg(l), nxr_mg(l), 2
1121                      DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
1122                         !DIR$ IVDEP
1123                         DO  k = ind_even_odd+1, nzt_mg(l)
1124                            km1 = k-ind_even_odd-1
1125                            kp1 = k-ind_even_odd
1126
1127                            p_mg(k,j,i) = 1.0_wp / f1_mg(k,l) * (              &
1128                               ddx2_mg(l) *                                    &
1129                               ( p_mg(k,j,i+1) + IBITS( flags(k,j,i), 5, 1 ) * &
1130                               ( p_mg(k,j,i) - p_mg(k,j,i+1) ) +               &
1131                               p_mg(k,j,i-1) + IBITS( flags(k,j,i), 4, 1 ) *   &
1132                               ( p_mg(k,j,i) - p_mg(k,j,i-1) ) )               &
1133                               + ddy2_mg(l) *                                  &
1134                               ( p_mg(k,j+1,i) + IBITS( flags(k,j,i), 3, 1 ) * &
1135                               ( p_mg(k,j,i) - p_mg(k,j+1,i) ) +               &
1136                               p_mg(k,j-1,i) + IBITS( flags(k,j,i), 2, 1 ) *   &
1137                               ( p_mg(k,j,i) - p_mg(k,j-1,i) ) )               &
1138                               + f2_mg(k,l) * p_mg(kp1,j,i)                    &
1139                               + f3_mg(k,l) *                                  &
1140                               ( p_mg(km1,j,i) + IBITS( flags(k,j,i), 0, 1 ) * &
1141                               ( p_mg(k,j,i) - p_mg(km1,j,i) ) )               &
1142                               - f_mg(k,j,i)               )
1143                         ENDDO
1144                      ENDDO
1145                   ENDDO
1146   
1147                   !$OMP DO
1148                   DO  i = nxl_mg(l)+1, nxr_mg(l), 2
1149                      DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
1150                          !DIR$ IVDEP
1151                          DO  k = ind_even_odd+1, nzt_mg(l)
1152                            km1 = k-ind_even_odd-1
1153                            kp1 = k-ind_even_odd
1154
1155                            p_mg(k,j,i) = 1.0_wp / f1_mg(k,l) * (              &
1156                               ddx2_mg(l) *                                    &
1157                               ( p_mg(k,j,i+1) + IBITS( flags(k,j,i), 5, 1 ) * &
1158                               ( p_mg(k,j,i) - p_mg(k,j,i+1) ) +               &
1159                               p_mg(k,j,i-1) + IBITS( flags(k,j,i), 4, 1 ) *   &
1160                               ( p_mg(k,j,i) - p_mg(k,j,i-1) ) )               &
1161                               + ddy2_mg(l) *                                  &
1162                               ( p_mg(k,j+1,i) + IBITS( flags(k,j,i), 3, 1 ) * &
1163                               ( p_mg(k,j,i) - p_mg(k,j+1,i) ) +               &
1164                               p_mg(k,j-1,i) + IBITS( flags(k,j,i), 2, 1 ) *   &
1165                               ( p_mg(k,j,i) - p_mg(k,j-1,i) ) )               &
1166                               + f2_mg(k,l) * p_mg(kp1,j,i)                    &
1167                               + f3_mg(k,l) *                                  &
1168                               ( p_mg(km1,j,i) + IBITS( flags(k,j,i), 0, 1 ) * &
1169                               ( p_mg(k,j,i) - p_mg(km1,j,i) ) )               &
1170                               - f_mg(k,j,i)               )
1171                         ENDDO
1172                      ENDDO
1173                   ENDDO
1174
1175                   !$OMP DO
1176                   DO  i = nxl_mg(l), nxr_mg(l), 2
1177                      DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
1178                         !DIR$ IVDEP
1179                         DO  k = nzb+1, ind_even_odd
1180                            km1 = k+ind_even_odd
1181                            kp1 = k+ind_even_odd+1
1182                            p_mg(k,j,i) = 1.0_wp / f1_mg(k,l) * (              &
1183                               ddx2_mg(l) *                                    &
1184                               ( p_mg(k,j,i+1) + IBITS( flags(k,j,i), 5, 1 ) * &
1185                               ( p_mg(k,j,i) - p_mg(k,j,i+1) ) +               &
1186                               p_mg(k,j,i-1) + IBITS( flags(k,j,i), 4, 1 ) *   &
1187                               ( p_mg(k,j,i) - p_mg(k,j,i-1) ) )               &
1188                               + ddy2_mg(l) *                                  &
1189                               ( p_mg(k,j+1,i) + IBITS( flags(k,j,i), 3, 1 ) * &
1190                               ( p_mg(k,j,i) - p_mg(k,j+1,i) ) +               &
1191                               p_mg(k,j-1,i) + IBITS( flags(k,j,i), 2, 1 ) *   &
1192                               ( p_mg(k,j,i) - p_mg(k,j-1,i) ) )               &
1193                               + f2_mg(k,l) * p_mg(kp1,j,i)                    &
1194                               + f3_mg(k,l) *                                  &
1195                               ( p_mg(km1,j,i) + IBITS( flags(k,j,i), 0, 1 ) * &
1196                               ( p_mg(k,j,i) - p_mg(km1,j,i) ) ) &
1197                               - f_mg(k,j,i)               )
1198                         ENDDO
1199                      ENDDO
1200                   ENDDO
1201
1202                   !$OMP DO
1203                   DO  i = nxl_mg(l)+1, nxr_mg(l), 2
1204                      DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
1205                         !DIR$ IVDEP
1206                         DO  k = nzb+1, ind_even_odd
1207                            km1 = k+ind_even_odd
1208                            kp1 = k+ind_even_odd+1
1209
1210                            p_mg(k,j,i) = 1.0_wp / f1_mg(k,l) * (              &
1211                               ddx2_mg(l) *                                    &
1212                               ( p_mg(k,j,i+1) + IBITS( flags(k,j,i), 5, 1 ) * &
1213                               ( p_mg(k,j,i) - p_mg(k,j,i+1) ) +               &
1214                               p_mg(k,j,i-1) + IBITS( flags(k,j,i), 4, 1 ) *   &
1215                               ( p_mg(k,j,i) - p_mg(k,j,i-1) ) )               &
1216                               + ddy2_mg(l) *                                  &
1217                               ( p_mg(k,j+1,i) + IBITS( flags(k,j,i), 3, 1 ) * &
1218                               ( p_mg(k,j,i) - p_mg(k,j+1,i) ) +               &
1219                               p_mg(k,j-1,i) + IBITS( flags(k,j,i), 2, 1 ) *   &
1220                               ( p_mg(k,j,i) - p_mg(k,j-1,i) ) )               &
1221                               + f2_mg(k,l) * p_mg(kp1,j,i)                    &
1222                               + f3_mg(k,l) *                                  &
1223                               ( p_mg(km1,j,i) + IBITS( flags(k,j,i), 0, 1 ) * &
1224                               ( p_mg(k,j,i) - p_mg(km1,j,i) ) )               &
1225                               - f_mg(k,j,i)               )
1226                         ENDDO
1227                      ENDDO
1228                   ENDDO
1229                   !$OMP END PARALLEL
1230
1231                   CALL cpu_log( log_point_s(36), 'redblack_no_unroll', 'stop' )
1232
1233                ELSE
1234!
1235!--                Loop unrolling along y, only one i loop for better cache use
1236                   CALL cpu_log( log_point_s(38), 'redblack_unroll', 'start' )
1237
1238                   !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,km1,kp1,jj)
1239                   !$OMP DO
1240                   DO  ic = nxl_mg(l), nxr_mg(l), 2
1241                      DO  jc = nys_mg(l), nyn_mg(l), 4
1242                         i  = ic
1243                         jj = jc+2-color
1244                         !DIR$ IVDEP
1245                         DO  k = ind_even_odd+1, nzt_mg(l)
1246
1247                            km1 = k-ind_even_odd-1
1248                            kp1 = k-ind_even_odd
1249
1250                            j = jj
1251                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1252                                ddx2_mg(l) *                                                 &
1253                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1254                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1255                                + ddy2_mg(l) *                                               &
1256                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1257                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1258                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1259                                + f3_mg(k,l) *                                               &
1260                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1261                                - f_mg(k,j,i)               )
1262                            j = jj+2
1263                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1264                                ddx2_mg(l) *                                                 &
1265                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1266                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1267                                + ddy2_mg(l) *                                               &
1268                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1269                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1270                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1271                                + f3_mg(k,l) *                                               &
1272                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1273                                - f_mg(k,j,i)               )
1274
1275                         ENDDO
1276
1277                         i  = ic+1
1278                         jj = jc+color-1
1279                         !DIR$ IVDEP
1280                         DO  k = ind_even_odd+1, nzt_mg(l)
1281
1282                            km1 = k-ind_even_odd-1
1283                            kp1 = k-ind_even_odd
1284
1285                            j =jj
1286                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1287                                ddx2_mg(l) *                                                 &
1288                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1289                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1290                                + ddy2_mg(l) *                                               &
1291                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1292                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1293                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1294                                + f3_mg(k,l) *                                               &
1295                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1296                                - f_mg(k,j,i)               )
1297
1298                            j = jj+2
1299                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1300                                ddx2_mg(l) *                                                 &
1301                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1302                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1303                                + ddy2_mg(l) *                                               &
1304                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1305                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1306                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1307                                + f3_mg(k,l) *                                               &
1308                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1309                                - f_mg(k,j,i)               )
1310
1311                         ENDDO
1312
1313                         i  = ic
1314                         jj = jc+color-1
1315                         !DIR$ IVDEP
1316                         DO  k = nzb+1, ind_even_odd
1317
1318                            km1 = k+ind_even_odd
1319                            kp1 = k+ind_even_odd+1
1320
1321                            j =jj
1322                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1323                                ddx2_mg(l) *                                                 &
1324                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1325                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1326                                + ddy2_mg(l) *                                               &
1327                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1328                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1329                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1330                                + f3_mg(k,l) *                                               &
1331                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1332                                - f_mg(k,j,i)               )
1333
1334                            j = jj+2
1335                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1336                                ddx2_mg(l) *                                                 &
1337                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1338                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1339                                + ddy2_mg(l) *                                               &
1340                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1341                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1342                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1343                                + f3_mg(k,l) *                                               &
1344                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1345                                - f_mg(k,j,i)               )
1346
1347                         ENDDO
1348
1349                         i  = ic+1
1350                         jj = jc+2-color
1351                         !DIR$ IVDEP
1352                         DO  k = nzb+1, ind_even_odd
1353
1354                            km1 = k+ind_even_odd
1355                            kp1 = k+ind_even_odd+1
1356
1357                            j =jj
1358                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1359                                ddx2_mg(l) *                                                 &
1360                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1361                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1362                                + ddy2_mg(l) *                                               &
1363                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1364                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1365                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1366                                + f3_mg(k,l) *                                               &
1367                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1368                                - f_mg(k,j,i)               )
1369
1370                            j = jj+2
1371                            p_mg(k,j,i) = 1.0 / f1_mg(k,l) * (                               &
1372                                ddx2_mg(l) *                                                 &
1373                                ( MERGE(p_mg(k,j,i),p_mg(k,j,i+1),BTEST( flags(k,j,i), 5))   &
1374                                + MERGE(p_mg(k,j,i),p_mg(k,j,i-1),BTEST( flags(k,j,i), 4)) ) &
1375                                + ddy2_mg(l) *                                               &
1376                                ( MERGE(p_mg(k,j,i),p_mg(k,j+1,i),BTEST( flags(k,j,i), 3))   &
1377                                + MERGE(p_mg(k,j,i),p_mg(k,j-1,i),BTEST( flags(k,j,i), 2)) ) &
1378                                + f2_mg(k,l) * p_mg(kp1,j,i)                                 &
1379                                + f3_mg(k,l) *                                               &
1380                                MERGE(p_mg(k,j,i),p_mg(km1,j,i),BTEST( flags(k,j,i), 0))     &
1381                                - f_mg(k,j,i)               )
1382
1383                         ENDDO
1384
1385                      ENDDO
1386                   ENDDO
1387                   !$OMP END PARALLEL
1388
1389                   CALL cpu_log( log_point_s(38), 'redblack_unroll', 'stop' )
1390
1391                ENDIF
1392
1393             ENDIF
1394
1395!
1396!--          Horizontal boundary conditions
1397             CALL special_exchange_horiz( p_mg, color )
1398
1399              IF ( .NOT. bc_lr_cyc )  THEN
1400                IF ( inflow_l .OR. outflow_l )  THEN
1401                   p_mg(:,:,nxl_mg(l)-1) = p_mg(:,:,nxl_mg(l))
1402                ENDIF
1403                IF ( inflow_r .OR. outflow_r )  THEN
1404                   p_mg(:,:,nxr_mg(l)+1) = p_mg(:,:,nxr_mg(l))
1405                ENDIF
1406             ENDIF
1407
1408             IF ( .NOT. bc_ns_cyc )  THEN
1409                IF ( inflow_n .OR. outflow_n )  THEN
1410                   p_mg(:,nyn_mg(l)+1,:) = p_mg(:,nyn_mg(l),:)
1411                ENDIF
1412                IF ( inflow_s .OR. outflow_s )  THEN
1413                   p_mg(:,nys_mg(l)-1,:) = p_mg(:,nys_mg(l),:)
1414                ENDIF
1415             ENDIF
1416
1417!
1418!--          Bottom and top boundary conditions
1419             IF ( ibc_p_b == 1 )  THEN
1420                p_mg(nzb,:,: ) = p_mg(ind_even_odd+1,:,:)
1421             ELSE
1422                p_mg(nzb,:,: ) = 0.0_wp
1423             ENDIF
1424
1425             IF ( ibc_p_t == 1 )  THEN
1426                p_mg(nzt_mg(l)+1,:,: ) = p_mg(ind_even_odd,:,:)
1427             ELSE
1428                p_mg(nzt_mg(l)+1,:,: ) = 0.0_wp
1429             ENDIF
1430
1431          ENDDO
1432       ENDDO
1433!
1434!--    Set pressure within topography and at the topography surfaces
1435       IF ( topography /= 'flat'  .AND.  .NOT. masking_method )  THEN
1436
1437          !$OMP PARALLEL PRIVATE (i,j,k,kp1,wall_left,wall_north,wall_right,   &
1438          !$OMP                   wall_south,wall_top,wall_total)
1439          !$OMP DO
1440          DO  i = nxl_mg(l), nxr_mg(l)
1441             DO  j = nys_mg(l), nyn_mg(l)
1442                !DIR$ IVDEP
1443                DO  k = ind_even_odd+1, nzt_mg(l)
1444                   km1 = k-ind_even_odd-1
1445                   kp1 = k-ind_even_odd
1446!
1447!--                First, set pressure inside topography to zero
1448                   p_mg(k,j,i) = p_mg(k,j,i) *                                 &
1449                                 ( 1.0_wp - IBITS( flags(k,j,i), 6, 1 ) )
1450!
1451!--                Second, determine if the gridpoint inside topography is
1452!--                adjacent to a wall and set its value to a value given by the
1453!--                average of those values obtained from Neumann boundary
1454!--                condition
1455                   wall_left  = IBITS( flags(k,j,i-1), 5, 1 )
1456                   wall_right = IBITS( flags(k,j,i+1), 4, 1 )
1457                   wall_south = IBITS( flags(k,j-1,i), 3, 1 )
1458                   wall_north = IBITS( flags(k,j+1,i), 2, 1 )
1459                   wall_top   = IBITS( flags(kp1,j,i), 0, 1 )
1460                   wall_total = wall_left + wall_right + wall_south +          &
1461                                wall_north + wall_top
1462
1463                   IF ( wall_total > 0.0_wp )  THEN
1464                      p_mg(k,j,i) = 1.0_wp / wall_total *                      &
1465                                    ( wall_left  * p_mg(k,j,i-1) +             &
1466                                      wall_right * p_mg(k,j,i+1) +             &
1467                                      wall_south * p_mg(k,j-1,i) +             &
1468                                      wall_north * p_mg(k,j+1,i) +             &
1469                                      wall_top   * p_mg(kp1,j,i) )
1470                   ENDIF
1471                ENDDO
1472
1473                !DIR$ IVDEP
1474                DO  k = nzb+1, ind_even_odd
1475                   km1 = k+ind_even_odd
1476                   kp1 = k+ind_even_odd+1
1477!
1478!--                First, set pressure inside topography to zero
1479                   p_mg(k,j,i) = p_mg(k,j,i) *                                 &
1480                                 ( 1.0_wp - IBITS( flags(k,j,i), 6, 1 ) )
1481!
1482!--                Second, determine if the gridpoint inside topography is
1483!--                adjacent to a wall and set its value to a value given by the
1484!--                average of those values obtained from Neumann boundary
1485!--                condition
1486                   wall_left  = IBITS( flags(k,j,i-1), 5, 1 )
1487                   wall_right = IBITS( flags(k,j,i+1), 4, 1 )
1488                   wall_south = IBITS( flags(k,j-1,i), 3, 1 )
1489                   wall_north = IBITS( flags(k,j+1,i), 2, 1 )
1490                   wall_top   = IBITS( flags(kp1,j,i), 0, 1 )
1491                   wall_total = wall_left + wall_right + wall_south +          &
1492                                wall_north + wall_top
1493
1494                   IF ( wall_total > 0.0_wp )  THEN
1495                      p_mg(k,j,i) = 1.0_wp / wall_total *                      &
1496                                    ( wall_left  * p_mg(k,j,i-1) +             &
1497                                      wall_right * p_mg(k,j,i+1) +             &
1498                                      wall_south * p_mg(k,j-1,i) +             &
1499                                      wall_north * p_mg(k,j+1,i) +             &
1500                                      wall_top   * p_mg(kp1,j,i) )
1501                   ENDIF
1502                ENDDO
1503             ENDDO
1504          ENDDO
1505          !$OMP END PARALLEL
1506
1507!
1508!--       One more time horizontal boundary conditions
1509          CALL exchange_horiz( p_mg, 1)
1510
1511       ENDIF
1512
1513    END SUBROUTINE redblack_fast
1514
1515
1516
1517    SUBROUTINE sort_k_to_even_odd_blocks( p_mg , glevel )
1518
1519!------------------------------------------------------------------------------!
1520! Description:
1521! ------------
1522! Sort k-Dimension from sequential into blocks of even and odd.
1523! This is required to vectorize the red-black subroutine.
1524! Version for 3D-REAL arrays
1525!------------------------------------------------------------------------------!
1526
1527       USE control_parameters,                                                 &
1528           ONLY:  grid_level
1529
1530       USE indices,                                                            &
1531           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1532
1533       IMPLICIT NONE
1534
1535       INTEGER(iwp), INTENT(IN) ::  glevel
1536
1537       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1,                               &
1538                           nys_mg(glevel)-1:nyn_mg(glevel)+1,                  &
1539                           nxl_mg(glevel)-1:nxr_mg(glevel)+1) ::  p_mg  !:
1540!
1541!--    Local variables
1542       INTEGER(iwp) :: i        !:
1543       INTEGER(iwp) :: j        !:
1544       INTEGER(iwp) :: k        !:
1545       INTEGER(iwp) :: l        !:
1546       INTEGER(iwp) :: ind      !:
1547       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1) ::  tmp  !:
1548
1549
1550       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'start' )
1551
1552       l = glevel
1553       ind_even_odd = even_odd_level(l)
1554
1555       !$OMP PARALLEL PRIVATE (i,j,k,ind,tmp)
1556       !$OMP DO
1557       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1558          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1559
1560!
1561!--          Sort the data with even k index
1562             ind = nzb-1
1563             DO  k = nzb, nzt_mg(l), 2
1564                ind = ind + 1
1565                tmp(ind) = p_mg(k,j,i)
1566             ENDDO
1567!
1568!--          Sort the data with odd k index
1569             DO  k = nzb+1, nzt_mg(l)+1, 2
1570                ind = ind + 1
1571                tmp(ind) = p_mg(k,j,i)
1572             ENDDO
1573
1574             p_mg(:,j,i) = tmp
1575
1576          ENDDO
1577       ENDDO
1578       !$OMP END PARALLEL
1579
1580       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'stop' )
1581
1582    END SUBROUTINE sort_k_to_even_odd_blocks
1583
1584
1585
1586    SUBROUTINE sort_k_to_even_odd_blocks_1d( f_mg, f_mg_b, glevel )
1587
1588!------------------------------------------------------------------------------!
1589! Description:
1590! ------------
1591! Sort k-Dimension from sequential into blocks of even and odd.
1592! This is required to vectorize the red-black subroutine.
1593! Version for 1D-REAL arrays
1594!------------------------------------------------------------------------------!
1595
1596       USE indices,                                                            &
1597           ONLY:  nzb, nzt_mg
1598
1599       IMPLICIT NONE
1600
1601       INTEGER(iwp), INTENT(IN) ::  glevel
1602
1603       REAL(wp), DIMENSION(nzb+1:nzt_mg(glevel)) ::  f_mg
1604       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1) ::  f_mg_b
1605
1606!
1607!--    Local variables
1608       INTEGER(iwp) :: ind   !:
1609       INTEGER(iwp) :: k     !:
1610
1611
1612       ind = nzb - 1
1613!
1614!--    Sort the data with even k index
1615       DO  k = nzb, nzt_mg(glevel), 2
1616          ind = ind + 1
1617          IF ( k >= nzb+1  .AND.  k <= nzt_mg(glevel) )  THEN
1618             f_mg_b(ind) = f_mg(k)
1619          ENDIF
1620       ENDDO
1621!
1622!--    Sort the data with odd k index
1623       DO  k = nzb+1, nzt_mg(glevel)+1, 2
1624          ind = ind + 1
1625          IF( k >= nzb+1  .AND.  k <= nzt_mg(glevel) )  THEN
1626             f_mg_b(ind) = f_mg(k)
1627          ENDIF
1628       ENDDO
1629
1630    END SUBROUTINE sort_k_to_even_odd_blocks_1d
1631
1632
1633
1634    SUBROUTINE sort_k_to_even_odd_blocks_int( i_mg , glevel )
1635
1636!------------------------------------------------------------------------------!
1637! Description:
1638! ------------
1639! Sort k-Dimension from sequential into blocks of even and odd.
1640! This is required to vectorize the red-black subroutine.
1641! Version for 2D-INTEGER arrays
1642!------------------------------------------------------------------------------!
1643
1644       USE control_parameters,                                                 &
1645           ONLY:  grid_level
1646
1647       USE indices,                                                            &
1648           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1649
1650       IMPLICIT NONE
1651
1652       INTEGER(iwp), INTENT(IN) ::  glevel
1653
1654       INTEGER(iwp), DIMENSION(nzb:nzt_mg(glevel)+1,                           &
1655                               nys_mg(glevel)-1:nyn_mg(glevel)+1,              &
1656                               nxl_mg(glevel)-1:nxr_mg(glevel)+1) ::  i_mg  !:
1657!
1658!--    Local variables
1659       INTEGER(iwp) :: i        !:
1660       INTEGER(iwp) :: j        !:
1661       INTEGER(iwp) :: k        !:
1662       INTEGER(iwp) :: l        !:
1663       INTEGER(iwp) :: ind      !:
1664       INTEGER(iwp),DIMENSION(nzb:nzt_mg(glevel)+1) ::  tmp
1665
1666
1667       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'start' )
1668
1669       l = glevel
1670       ind_even_odd = even_odd_level(l)
1671
1672       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1673          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1674
1675!
1676!--          Sort the data with even k index
1677             ind = nzb-1
1678             DO  k = nzb, nzt_mg(l), 2
1679                ind = ind + 1
1680                tmp(ind) = i_mg(k,j,i)
1681             ENDDO
1682!
1683!++          ATTENTION: Check reason for this error. Remove it or replace WRITE
1684!++                     by PALM message
1685             IF ( ind /= ind_even_odd )  THEN
1686                WRITE (0,*) 'ERROR ==> illegal ind_even_odd ',ind,ind_even_odd,l
1687                CALL MPI_ABORT(MPI_COMM_WORLD,i,j)
1688             ENDIF
1689!
1690!--          Sort the data with odd k index
1691             DO  k = nzb+1, nzt_mg(l)+1, 2
1692                ind = ind + 1
1693                tmp(ind) = i_mg(k,j,i)
1694             ENDDO
1695
1696             i_mg(:,j,i) = tmp
1697
1698          ENDDO
1699       ENDDO
1700
1701       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'stop' )
1702
1703    END SUBROUTINE sort_k_to_even_odd_blocks_int
1704
1705
1706
1707    SUBROUTINE sort_k_to_sequential( p_mg )
1708
1709!------------------------------------------------------------------------------!
1710! Description:
1711! ------------
1712! Sort k-dimension from blocks of even and odd into sequential
1713!------------------------------------------------------------------------------!
1714
1715       USE control_parameters,                                                 &
1716           ONLY:  grid_level
1717
1718       USE indices,                                                            &
1719           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1720
1721       IMPLICIT NONE
1722
1723       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1724                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1725                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  p_mg  !:
1726!
1727!--    Local variables
1728       INTEGER(iwp) :: i        !:
1729       INTEGER(iwp) :: j        !:
1730       INTEGER(iwp) :: k        !:
1731       INTEGER(iwp) :: l        !:
1732       INTEGER(iwp) :: ind      !:
1733
1734       REAL(wp),DIMENSION(nzb:nzt_mg(grid_level)+1) ::  tmp
1735
1736
1737       l = grid_level
1738
1739       !$OMP PARALLEL PRIVATE (i,j,k,ind,tmp)
1740       !$OMP DO
1741       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1742          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1743
1744             ind = nzb - 1
1745             tmp = p_mg(:,j,i)
1746             DO  k = nzb, nzt_mg(l), 2
1747                ind = ind + 1
1748                p_mg(k,j,i) = tmp(ind)
1749             ENDDO
1750
1751             DO  k = nzb+1, nzt_mg(l)+1, 2
1752                ind = ind + 1
1753                p_mg(k,j,i) = tmp(ind)
1754             ENDDO
1755          ENDDO
1756       ENDDO
1757       !$OMP END PARALLEL
1758
1759    END SUBROUTINE sort_k_to_sequential
1760
1761
1762
1763    SUBROUTINE mg_gather_fast( f2, f2_sub )
1764
1765       USE control_parameters,                                                 &
1766           ONLY:  grid_level
1767
1768       USE cpulog,                                                             &
1769           ONLY:  cpu_log, log_point_s
1770
1771       USE indices,                                                            &
1772           ONLY:  mg_loc_ind, nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1773
1774       IMPLICIT NONE
1775
1776       INTEGER(iwp) ::  i       !:
1777       INTEGER(iwp) ::  il      !:
1778       INTEGER(iwp) ::  ir      !:
1779       INTEGER(iwp) ::  j       !:
1780       INTEGER(iwp) ::  jn      !:
1781       INTEGER(iwp) ::  js      !:
1782       INTEGER(iwp) ::  k       !:
1783       INTEGER(iwp) ::  nwords  !:
1784
1785       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1786                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
1787                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f2    !:
1788       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1789                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
1790                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f2_l  !:
1791
1792       REAL(wp), DIMENSION(nzb:mg_loc_ind(5,myid)+1,                           &
1793                           mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,          &
1794                           mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) ::  f2_sub  !:
1795
1796
1797#if defined( __parallel )
1798       CALL cpu_log( log_point_s(34), 'mg_gather', 'start' )
1799
1800       f2_l = 0.0_wp
1801
1802!
1803!--    Store the local subdomain array on the total array
1804       js = mg_loc_ind(3,myid)
1805       IF ( south_border_pe )  js = js - 1
1806       jn = mg_loc_ind(4,myid)
1807       IF ( north_border_pe )  jn = jn + 1
1808       il = mg_loc_ind(1,myid)
1809       IF ( left_border_pe )   il = il - 1
1810       ir = mg_loc_ind(2,myid)
1811       IF ( right_border_pe )  ir = ir + 1
1812       DO  i = il, ir
1813          DO  j = js, jn
1814             DO  k = nzb, nzt_mg(grid_level)+1
1815                f2_l(k,j,i) = f2_sub(k,j,i)
1816             ENDDO
1817          ENDDO
1818       ENDDO
1819
1820!
1821!--    Find out the number of array elements of the total array
1822       nwords = SIZE( f2 )
1823
1824!
1825!--    Gather subdomain data from all PEs
1826       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1827       CALL MPI_ALLREDUCE( f2_l(nzb,nys_mg(grid_level)-1,nxl_mg(grid_level)-1), &
1828                           f2(nzb,nys_mg(grid_level)-1,nxl_mg(grid_level)-1),   &
1829                           nwords, MPI_REAL, MPI_SUM, comm2d, ierr )
1830
1831       CALL cpu_log( log_point_s(34), 'mg_gather', 'stop' )
1832#endif
1833   
1834    END SUBROUTINE mg_gather_fast
1835
1836
1837
1838    SUBROUTINE mg_scatter_fast( p2, p2_sub )
1839!
1840!--    TODO: It might be possible to improve the speed of this routine by using
1841!--    non-blocking communication
1842
1843       USE control_parameters,                                                 &
1844           ONLY:  grid_level
1845
1846       USE cpulog,                                                             &
1847           ONLY:  cpu_log, log_point_s
1848
1849       USE indices,                                                            &
1850           ONLY:  mg_loc_ind, nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1851
1852       IMPLICIT NONE
1853
1854       INTEGER(iwp) ::  nwords  !:
1855
1856       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1857                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1858                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  p2  !:
1859
1860       REAL(wp), DIMENSION(nzb:mg_loc_ind(5,myid)+1,                           &
1861                           mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,          &
1862                           mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) ::  p2_sub  !:
1863
1864!
1865!--    Find out the number of array elements of the subdomain array
1866       nwords = SIZE( p2_sub )
1867
1868#if defined( __parallel )
1869       CALL cpu_log( log_point_s(35), 'mg_scatter', 'start' )
1870
1871       p2_sub = p2(:,mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1, &
1872                     mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1)
1873
1874       CALL cpu_log( log_point_s(35), 'mg_scatter', 'stop' )
1875#endif
1876   
1877    END SUBROUTINE mg_scatter_fast
1878
1879
1880
1881    RECURSIVE SUBROUTINE next_mg_level_fast( f_mg, p_mg, p3, r )
1882
1883!------------------------------------------------------------------------------!
1884! Description:
1885! ------------
1886! This is where the multigrid technique takes place. V- and W- Cycle are
1887! implemented and steered by the parameter "gamma". Parameter "nue" determines
1888! the convergence of the multigrid iterative solution. There are nue times
1889! RB-GS iterations. It should be set to "1" or "2", considering the time effort
1890! one would like to invest. Last choice shows a very good converging factor,
1891! but leads to an increase in computing time.
1892!------------------------------------------------------------------------------!
1893
1894       USE control_parameters,                                                 &
1895           ONLY:  bc_lr_dirrad, bc_lr_raddir, bc_ns_dirrad, bc_ns_raddir,      &
1896                  gamma_mg, grid_level, grid_level_count, ibc_p_b, ibc_p_t,    &
1897                  inflow_l, inflow_n, inflow_r, inflow_s, maximum_grid_level,  &
1898                  mg_switch_to_pe0_level, mg_switch_to_pe0, ngsrb, outflow_l,  &
1899                  outflow_n, outflow_r, outflow_s
1900
1901       USE indices,                                                            &
1902           ONLY:  mg_loc_ind, nxl, nxl_mg, nxr, nxr_mg, nys, nys_mg, nyn,      &
1903                  nyn_mg, nzb, nzt, nzt_mg
1904
1905       IMPLICIT NONE
1906
1907       INTEGER(iwp) ::  i            !:
1908       INTEGER(iwp) ::  j            !:
1909       INTEGER(iwp) ::  k            !:
1910       INTEGER(iwp) ::  nxl_mg_save  !:
1911       INTEGER(iwp) ::  nxr_mg_save  !:
1912       INTEGER(iwp) ::  nyn_mg_save  !:
1913       INTEGER(iwp) ::  nys_mg_save  !:
1914       INTEGER(iwp) ::  nzt_mg_save  !:
1915
1916       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1917                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1918                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: f_mg  !:
1919       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1920                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1921                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: p_mg  !:
1922       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1923                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1924                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: p3    !:
1925       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1926                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1927                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: r     !:
1928
1929       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1930                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1931                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  f2  !:
1932       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1933                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1934                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  p2  !:
1935
1936       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  f2_sub  !:
1937       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  p2_sub  !:
1938
1939!
1940!--    Restriction to the coarsest grid
1941    10 IF ( grid_level == 1 )  THEN
1942
1943!
1944!--       Solution on the coarsest grid. Double the number of Gauss-Seidel
1945!--       iterations in order to get a more accurate solution.
1946          ngsrb = 2 * ngsrb
1947
1948          ind_even_odd = even_odd_level(grid_level)
1949
1950          CALL redblack_fast( f_mg, p_mg )
1951
1952          ngsrb = ngsrb / 2
1953
1954
1955       ELSEIF ( grid_level /= 1 )  THEN
1956
1957          grid_level_count(grid_level) = grid_level_count(grid_level) + 1
1958
1959!
1960!--       Solution on the actual grid level
1961          ind_even_odd = even_odd_level(grid_level)
1962
1963          CALL redblack_fast( f_mg, p_mg )
1964
1965!
1966!--       Determination of the actual residual
1967          CALL resid_fast( f_mg, p_mg, r )
1968
1969!--       Restriction of the residual (finer grid values!) to the next coarser
1970!--       grid. Therefore, the grid level has to be decremented now. nxl..nzt have
1971!--       to be set to the coarse grid values, because these variables are needed
1972!--       for the exchange of ghost points in routine exchange_horiz
1973          grid_level = grid_level - 1
1974
1975          nxl = nxl_mg(grid_level)
1976          nys = nys_mg(grid_level)
1977          nxr = nxr_mg(grid_level)
1978          nyn = nyn_mg(grid_level)
1979          nzt = nzt_mg(grid_level)
1980
1981          IF ( grid_level == mg_switch_to_pe0_level )  THEN
1982
1983!
1984!--          From this level on, calculations are done on PE0 only.
1985!--          First, carry out restriction on the subdomain.
1986!--          Therefore, indices of the level have to be changed to subdomain
1987!--          values in between (otherwise, the restrict routine would expect
1988!--          the gathered array)
1989
1990             nxl_mg_save = nxl_mg(grid_level)
1991             nxr_mg_save = nxr_mg(grid_level)
1992             nys_mg_save = nys_mg(grid_level)
1993             nyn_mg_save = nyn_mg(grid_level)
1994             nzt_mg_save = nzt_mg(grid_level)
1995             nxl_mg(grid_level) = mg_loc_ind(1,myid)
1996             nxr_mg(grid_level) = mg_loc_ind(2,myid)
1997             nys_mg(grid_level) = mg_loc_ind(3,myid)
1998             nyn_mg(grid_level) = mg_loc_ind(4,myid)
1999             nzt_mg(grid_level) = mg_loc_ind(5,myid)
2000             nxl = mg_loc_ind(1,myid)
2001             nxr = mg_loc_ind(2,myid)
2002             nys = mg_loc_ind(3,myid)
2003             nyn = mg_loc_ind(4,myid)
2004             nzt = mg_loc_ind(5,myid)
2005
2006             ALLOCATE( f2_sub(nzb:nzt_mg(grid_level)+1,                    &
2007                              nys_mg(grid_level)-1:nyn_mg(grid_level)+1,   &
2008                              nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) )
2009
2010             CALL restrict_fast( f2_sub, r )
2011
2012!
2013!--          Restore the correct indices of this level
2014             nxl_mg(grid_level) = nxl_mg_save
2015             nxr_mg(grid_level) = nxr_mg_save
2016             nys_mg(grid_level) = nys_mg_save
2017             nyn_mg(grid_level) = nyn_mg_save
2018             nzt_mg(grid_level) = nzt_mg_save
2019             nxl = nxl_mg(grid_level)
2020             nxr = nxr_mg(grid_level)
2021             nys = nys_mg(grid_level)
2022             nyn = nyn_mg(grid_level)
2023             nzt = nzt_mg(grid_level)
2024!
2025!--          Gather all arrays from the subdomains on PE0
2026             CALL mg_gather_fast( f2, f2_sub )
2027
2028!
2029!--          Set switch for routine exchange_horiz, that no ghostpoint exchange
2030!--          has to be carried out from now on
2031             mg_switch_to_pe0 = .TRUE.
2032
2033!
2034!--          In case of non-cyclic lateral boundary conditions, both in- and
2035!--          outflow conditions have to be used on all PEs after the switch,
2036!--          because then they have the total domain.
2037             IF ( bc_lr_dirrad )  THEN
2038                inflow_l  = .TRUE.
2039                inflow_r  = .FALSE.
2040                outflow_l = .FALSE.
2041                outflow_r = .TRUE.
2042             ELSEIF ( bc_lr_raddir )  THEN
2043                inflow_l  = .FALSE.
2044                inflow_r  = .TRUE.
2045                outflow_l = .TRUE.
2046                outflow_r = .FALSE.
2047             ENDIF
2048
2049             IF ( bc_ns_dirrad )  THEN
2050                inflow_n  = .TRUE.
2051                inflow_s  = .FALSE.
2052                outflow_n = .FALSE.
2053                outflow_s = .TRUE.
2054             ELSEIF ( bc_ns_raddir )  THEN
2055                inflow_n  = .FALSE.
2056                inflow_s  = .TRUE.
2057                outflow_n = .TRUE.
2058                outflow_s = .FALSE.
2059             ENDIF
2060
2061             DEALLOCATE( f2_sub )
2062
2063          ELSE
2064
2065             CALL restrict_fast( f2, r )
2066
2067             ind_even_odd = even_odd_level(grid_level)  ! must be after restrict
2068
2069          ENDIF
2070
2071          p2 = 0.0_wp
2072
2073!
2074!--       Repeat the same procedure till the coarsest grid is reached
2075          CALL next_mg_level_fast( f2, p2, p3, r )
2076
2077       ENDIF
2078
2079!
2080!--    Now follows the prolongation
2081       IF ( grid_level >= 2 )  THEN
2082
2083!
2084!--       Prolongation of the new residual. The values are transferred
2085!--       from the coarse to the next finer grid.
2086          IF ( grid_level == mg_switch_to_pe0_level+1 )  THEN
2087
2088#if defined( __parallel )
2089!
2090!--          At this level, the new residual first has to be scattered from
2091!--          PE0 to the other PEs
2092             ALLOCATE( p2_sub(nzb:mg_loc_ind(5,myid)+1,             &
2093                       mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,   &
2094                       mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) )
2095
2096             CALL mg_scatter_fast( p2, p2_sub )
2097
2098!
2099!--          Therefore, indices of the previous level have to be changed to
2100!--          subdomain values in between (otherwise, the prolong routine would
2101!--          expect the gathered array)
2102             nxl_mg_save = nxl_mg(grid_level-1)
2103             nxr_mg_save = nxr_mg(grid_level-1)
2104             nys_mg_save = nys_mg(grid_level-1)
2105             nyn_mg_save = nyn_mg(grid_level-1)
2106             nzt_mg_save = nzt_mg(grid_level-1)
2107             nxl_mg(grid_level-1) = mg_loc_ind(1,myid)
2108             nxr_mg(grid_level-1) = mg_loc_ind(2,myid)
2109             nys_mg(grid_level-1) = mg_loc_ind(3,myid)
2110             nyn_mg(grid_level-1) = mg_loc_ind(4,myid)
2111             nzt_mg(grid_level-1) = mg_loc_ind(5,myid)
2112
2113!
2114!--          Set switch for routine exchange_horiz, that ghostpoint exchange
2115!--          has to be carried again out from now on
2116             mg_switch_to_pe0 = .FALSE.
2117
2118!
2119!--          For non-cyclic lateral boundary conditions, restore the
2120!--          in-/outflow conditions
2121             inflow_l  = .FALSE.;  inflow_r  = .FALSE.
2122             inflow_n  = .FALSE.;  inflow_s  = .FALSE.
2123             outflow_l = .FALSE.;  outflow_r = .FALSE.
2124             outflow_n = .FALSE.;  outflow_s = .FALSE.
2125
2126             IF ( pleft == MPI_PROC_NULL )  THEN
2127                IF ( bc_lr_dirrad )  THEN
2128                   inflow_l  = .TRUE.
2129                ELSEIF ( bc_lr_raddir )  THEN
2130                   outflow_l = .TRUE.
2131                ENDIF
2132             ENDIF
2133
2134             IF ( pright == MPI_PROC_NULL )  THEN
2135                IF ( bc_lr_dirrad )  THEN
2136                   outflow_r = .TRUE.
2137                ELSEIF ( bc_lr_raddir )  THEN
2138                   inflow_r  = .TRUE.
2139                ENDIF
2140             ENDIF
2141
2142             IF ( psouth == MPI_PROC_NULL )  THEN
2143                IF ( bc_ns_dirrad )  THEN
2144                   outflow_s = .TRUE.
2145                ELSEIF ( bc_ns_raddir )  THEN
2146                   inflow_s  = .TRUE.
2147                ENDIF
2148             ENDIF
2149
2150             IF ( pnorth == MPI_PROC_NULL )  THEN
2151                IF ( bc_ns_dirrad )  THEN
2152                   inflow_n  = .TRUE.
2153                ELSEIF ( bc_ns_raddir )  THEN
2154                   outflow_n = .TRUE.
2155                ENDIF
2156             ENDIF
2157
2158             CALL prolong_fast( p2_sub, p3 )
2159
2160!
2161!--          Restore the correct indices of the previous level
2162             nxl_mg(grid_level-1) = nxl_mg_save
2163             nxr_mg(grid_level-1) = nxr_mg_save
2164             nys_mg(grid_level-1) = nys_mg_save
2165             nyn_mg(grid_level-1) = nyn_mg_save
2166             nzt_mg(grid_level-1) = nzt_mg_save
2167
2168             DEALLOCATE( p2_sub )
2169#endif
2170
2171          ELSE
2172
2173             CALL prolong_fast( p2, p3 )
2174
2175          ENDIF
2176
2177!
2178!--       Computation of the new pressure correction. Therefore,
2179!--       values from prior grids are added up automatically stage by stage.
2180          DO  i = nxl_mg(grid_level)-1, nxr_mg(grid_level)+1
2181             DO  j = nys_mg(grid_level)-1, nyn_mg(grid_level)+1
2182                DO  k = nzb, nzt_mg(grid_level)+1
2183                   p_mg(k,j,i) = p_mg(k,j,i) + p3(k,j,i)
2184                ENDDO
2185             ENDDO
2186          ENDDO
2187
2188!
2189!--       Relaxation of the new solution
2190          CALL redblack_fast( f_mg, p_mg )
2191
2192       ENDIF
2193
2194
2195!
2196!--    The following few lines serve the steering of the multigrid scheme
2197       IF ( grid_level == maximum_grid_level )  THEN
2198
2199          GOTO 20
2200
2201       ELSEIF ( grid_level /= maximum_grid_level  .AND.  grid_level /= 1  .AND. &
2202                grid_level_count(grid_level) /= gamma_mg )  THEN
2203
2204          GOTO 10
2205
2206       ENDIF
2207
2208!
2209!--    Reset counter for the next call of poismg_fast
2210       grid_level_count(grid_level) = 0
2211
2212!
2213!--    Continue with the next finer level. nxl..nzt have to be
2214!--    set to the finer grid values, because these variables are needed for the
2215!--    exchange of ghost points in routine exchange_horiz
2216       grid_level = grid_level + 1
2217       ind_even_odd = even_odd_level(grid_level)
2218
2219       nxl = nxl_mg(grid_level)
2220       nxr = nxr_mg(grid_level)
2221       nys = nys_mg(grid_level)
2222       nyn = nyn_mg(grid_level)
2223       nzt = nzt_mg(grid_level)
2224
2225    20 CONTINUE
2226
2227    END SUBROUTINE next_mg_level_fast
2228
2229
2230
2231    SUBROUTINE init_even_odd_blocks
2232
2233!------------------------------------------------------------------------------!
2234! Description:
2235! ------------
2236! Initial settings for sorting k-dimension from sequential order (alternate
2237! even/odd) into blocks of even and odd or vice versa
2238!------------------------------------------------------------------------------!
2239
2240       USE arrays_3d,                                                          &
2241           ONLY:  f1_mg, f2_mg, f3_mg
2242
2243       USE control_parameters,                                                 &
2244           ONLY:  grid_level, masking_method, maximum_grid_level, topography
2245
2246       USE indices,                                                            &
2247           ONLY:  nzb, nzt, nzt_mg
2248
2249       USE indices,                                                            &
2250           ONLY:  flags, wall_flags_1, wall_flags_2, wall_flags_3,             &
2251                  wall_flags_4, wall_flags_5, wall_flags_6, wall_flags_7,      &
2252                  wall_flags_8, wall_flags_9, wall_flags_10, nxl_mg, nxr_mg,   &
2253                  nys_mg, nyn_mg, nzb, nzt_mg
2254
2255       IMPLICIT NONE
2256!
2257!--    Local variables
2258       INTEGER(iwp) ::  i     !:
2259       INTEGER(iwp) ::  l     !:
2260
2261       LOGICAL, SAVE ::  lfirst = .TRUE.
2262
2263
2264       IF ( .NOT. lfirst )  RETURN
2265
2266       ALLOCATE( even_odd_level(maximum_grid_level) )
2267
2268       ALLOCATE( f1_mg_b(nzb:nzt+1,maximum_grid_level),                        &
2269                 f2_mg_b(nzb:nzt+1,maximum_grid_level),                        &
2270                 f3_mg_b(nzb:nzt+1,maximum_grid_level) )
2271
2272!
2273!--    Set border index between the even and odd block
2274       DO  i = maximum_grid_level, 1, -1
2275          even_odd_level(i) = nzt_mg(i) / 2
2276       ENDDO
2277
2278!
2279!--    Sort grid coefficients used in red/black scheme and for calculating the
2280!--    residual to block (even/odd) structure
2281       DO  l = maximum_grid_level, 1 , -1
2282          CALL sort_k_to_even_odd_blocks( f1_mg(nzb+1:nzt_mg(grid_level),l),   &
2283                                          f1_mg_b(nzb:nzt_mg(grid_level)+1,l), &
2284                                          l )
2285          CALL sort_k_to_even_odd_blocks( f2_mg(nzb+1:nzt_mg(grid_level),l),   &
2286                                          f2_mg_b(nzb:nzt_mg(grid_level)+1,l), &
2287                                          l )
2288          CALL sort_k_to_even_odd_blocks( f3_mg(nzb+1:nzt_mg(grid_level),l),   &
2289                                          f3_mg_b(nzb:nzt_mg(grid_level)+1,l), &
2290                                          l )
2291       ENDDO
2292
2293!
2294!--    Sort flags array for all levels to block (even/odd) structure
2295       IF ( topography /= 'flat'  .AND.  .NOT. masking_method )  THEN
2296
2297          DO  l = maximum_grid_level, 1 , -1
2298!
2299!--          Assign the flag level to be calculated
2300             SELECT CASE ( l )
2301                CASE ( 1 )
2302                   flags => wall_flags_1
2303                CASE ( 2 )
2304                   flags => wall_flags_2
2305                CASE ( 3 )
2306                   flags => wall_flags_3
2307                CASE ( 4 )
2308                   flags => wall_flags_4
2309                CASE ( 5 )
2310                   flags => wall_flags_5
2311                CASE ( 6 )
2312                   flags => wall_flags_6
2313                CASE ( 7 )
2314                   flags => wall_flags_7
2315                CASE ( 8 )
2316                   flags => wall_flags_8
2317                CASE ( 9 )
2318                   flags => wall_flags_9
2319                CASE ( 10 )
2320                   flags => wall_flags_10
2321             END SELECT
2322
2323             CALL sort_k_to_even_odd_blocks( flags, l )
2324
2325          ENDDO
2326
2327       ENDIF
2328
2329       lfirst = .FALSE.
2330
2331     END SUBROUTINE init_even_odd_blocks
2332
2333
2334
2335     SUBROUTINE special_exchange_horiz ( p_mg, color )
2336
2337!------------------------------------------------------------------------------!
2338! Description:
2339! ------------
2340! Special exchange_horiz subroutine for use in redblack. Transfers only
2341! "red" or "black" data points.
2342!------------------------------------------------------------------------------!
2343
2344       USE control_parameters,                                                 &
2345           ONLY:  bc_lr_cyc, bc_ns_cyc, grid_level, ibc_p_b, ibc_p_t,          &
2346                  inflow_l, inflow_n, inflow_r, inflow_s, maximum_grid_level,  &
2347                  outflow_l, outflow_n, outflow_r, outflow_s,                  &
2348                  synchronous_exchange
2349
2350       USE indices,                                                            &
2351           ONLY:  mg_loc_ind, nxl, nxl_mg, nxr, nxr_mg, nys, nys_mg, nyn,      &
2352                  nyn_mg, nzb, nzt, nzt_mg
2353
2354       USE pegrid,                                                             &
2355           ONLY:  ngp_xz, ngp_yz
2356
2357       IMPLICIT NONE
2358
2359       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
2360                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
2361                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  p_mg  !:
2362
2363       INTEGER(iwp), intent(IN) ::  color  !:
2364!
2365!--    Local variables
2366       INTEGER(iwp) ::  i,i1,i2  !:
2367       INTEGER(iwp) ::  j,j1,j2  !:
2368       INTEGER(iwp) ::  k        !:
2369       INTEGER(iwp) ::  l        !:
2370       INTEGER(iwp) ::  jys      !:
2371       INTEGER(iwp) ::  jyn      !:
2372       INTEGER(iwp) ::  ixl      !:
2373       INTEGER(iwp) ::  ixr      !:
2374       logical      ::  sendrecv_in_background_save  !:
2375       logical      ::  synchronous_exchange_save    !:
2376
2377       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
2378                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
2379                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  temp  !:
2380
2381
2382       sendrecv_in_background_save = sendrecv_in_background
2383       sendrecv_in_background      = .FALSE.
2384       synchronous_exchange_save   = synchronous_exchange
2385       synchronous_exchange        = .FALSE.
2386
2387       l = grid_level
2388
2389       ind_even_odd = even_odd_level(grid_level)
2390
2391       jys = nys_mg(grid_level-1)
2392       jyn = nyn_mg(grid_level-1)
2393       ixl = nxl_mg(grid_level-1)
2394       ixr = nxr_mg(grid_level-1)
2395
2396!
2397!--    Restricted transfer only on finer levels with enough data
2398       IF ( ngp_xz(grid_level) >= 900 .OR. ngp_yz(grid_level) >= 900 )  THEN
2399!
2400!--       Handling the even k Values
2401!--       Collecting data for the north - south exchange
2402!--       Since only every second value has to be transfered, data are stored
2403!--       on the next coarser grid level, because the arrays on that level
2404!--       have just the required size
2405          i1 = nxl_mg(grid_level-1)
2406          i2 = nxl_mg(grid_level-1)
2407
2408          DO  i = nxl_mg(l), nxr_mg(l), 2
2409             DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2410
2411                IF ( j == nys_mg(l) )  THEN
2412                   !DIR$ IVDEP
2413                   DO  k = ind_even_odd+1, nzt_mg(l)
2414                      temp(k-ind_even_odd,jys,i1) = p_mg(k,j,i)
2415                   ENDDO
2416                   i1 = i1 + 1
2417
2418                ENDIF
2419
2420                IF ( j == nyn_mg(l) )  THEN
2421                   !DIR$ IVDEP
2422                   DO  k = ind_even_odd+1, nzt_mg(l)
2423                      temp(k-ind_even_odd,jyn,i2) = p_mg(k,j,i)
2424                   ENDDO
2425                   i2 = i2 + 1
2426
2427                ENDIF
2428
2429             ENDDO
2430          ENDDO
2431
2432          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2433             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2434
2435                IF ( j == nys_mg(l) )  THEN
2436                   !DIR$ IVDEP
2437                   DO  k = ind_even_odd+1, nzt_mg(l)
2438                      temp(k-ind_even_odd,jys,i1) = p_mg(k,j,i)
2439                   ENDDO
2440                   i1 = i1 + 1
2441
2442                ENDIF
2443
2444                IF ( j == nyn_mg(l) )  THEN
2445                   !DIR$ IVDEP
2446                   DO  k = ind_even_odd+1, nzt_mg(l)
2447                      temp(k-ind_even_odd,jyn,i2) = p_mg(k,j,i)
2448                   ENDDO
2449                   i2 = i2 + 1
2450
2451                ENDIF
2452
2453             ENDDO
2454          ENDDO
2455
2456          grid_level = grid_level-1
2457
2458          nxl = nxl_mg(grid_level)
2459          nys = nys_mg(grid_level)
2460          nxr = nxr_mg(grid_level)
2461          nyn = nyn_mg(grid_level)
2462          nzt = nzt_mg(grid_level)
2463
2464          send_receive = 'ns'
2465          CALL exchange_horiz( temp, 1 )
2466
2467          grid_level = grid_level+1
2468
2469          i1 = nxl_mg(grid_level-1)
2470          i2 = nxl_mg(grid_level-1)
2471
2472          DO  i = nxl_mg(l), nxr_mg(l), 2
2473             DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2474
2475                IF ( j == nys_mg(l) )  THEN
2476                   !DIR$ IVDEP
2477                   DO  k = ind_even_odd+1, nzt_mg(l)
2478                      p_mg(k,nyn_mg(l)+1,i) = temp(k-ind_even_odd,jyn+1,i1)
2479                   ENDDO
2480                   i1 = i1 + 1
2481
2482                ENDIF
2483
2484                IF ( j == nyn_mg(l) )  THEN
2485                   !DIR$ IVDEP
2486                   DO  k = ind_even_odd+1, nzt_mg(l)
2487                      p_mg(k,nys_mg(l)-1,i) = temp(k-ind_even_odd,jys-1,i2)
2488                   ENDDO
2489                   i2 = i2 + 1
2490
2491                ENDIF
2492
2493             ENDDO
2494          ENDDO
2495
2496          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2497             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2498
2499                IF ( j == nys_mg(l) )  THEN
2500                   !DIR$ IVDEP
2501                   DO  k = ind_even_odd+1, nzt_mg(l)
2502                      p_mg(k,nyn_mg(l)+1,i) = temp(k-ind_even_odd,jyn+1,i1)
2503                   ENDDO
2504                   i1 = i1 + 1
2505
2506                ENDIF
2507
2508                IF ( j == nyn_mg(l) )  THEN
2509                   !DIR$ IVDEP
2510                   DO  k = ind_even_odd+1, nzt_mg(l)
2511                      p_mg(k,nys_mg(l)-1,i) = temp(k-ind_even_odd,jys-1,i2)
2512                   ENDDO
2513                   i2 = i2 + 1
2514
2515                ENDIF
2516
2517             ENDDO
2518          ENDDO
2519
2520!
2521!--       Collecting data for the left - right exchange
2522!--       Since only every second value has to be transfered, data are stored
2523!--       on the next coarser grid level, because the arrays on that level
2524!--       have just the required size
2525          j1 = nys_mg(grid_level-1)
2526          j2 = nys_mg(grid_level-1)
2527
2528          DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2529             DO  i = nxl_mg(l), nxr_mg(l), 2
2530
2531                IF ( i == nxl_mg(l) )  THEN
2532                   !DIR$ IVDEP
2533                   DO  k = ind_even_odd+1, nzt_mg(l)
2534                      temp(k-ind_even_odd,j1,ixl) = p_mg(k,j,i)
2535                   ENDDO
2536                   j1 = j1 + 1
2537
2538                ENDIF
2539
2540                IF ( i == nxr_mg(l) )  THEN
2541                   !DIR$ IVDEP
2542                   DO  k = ind_even_odd+1, nzt_mg(l)
2543                      temp(k-ind_even_odd,j2,ixr) = p_mg(k,j,i)
2544                   ENDDO
2545                   j2 = j2 + 1
2546
2547                ENDIF
2548
2549             ENDDO
2550          ENDDO
2551
2552          DO j = nys_mg(l) + (color-1), nyn_mg(l), 2
2553             DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2554
2555                IF ( i == nxl_mg(l) )  THEN
2556                   !DIR$ IVDEP
2557                   DO  k = ind_even_odd+1, nzt_mg(l)
2558                      temp(k-ind_even_odd,j1,ixl) = p_mg(k,j,i)
2559                   ENDDO
2560                   j1 = j1 + 1
2561
2562                ENDIF
2563
2564                IF ( i == nxr_mg(l) )  THEN
2565                   !DIR$ IVDEP
2566                   DO  k = ind_even_odd+1, nzt_mg(l)
2567                      temp(k-ind_even_odd,j2,ixr) = p_mg(k,j,i)
2568                   ENDDO
2569                   j2 = j2 + 1
2570
2571                ENDIF
2572
2573             ENDDO
2574          ENDDO
2575
2576          grid_level = grid_level-1
2577          send_receive = 'lr'
2578
2579          CALL exchange_horiz( temp, 1 )
2580
2581          grid_level = grid_level+1
2582
2583          j1 = nys_mg(grid_level-1)
2584          j2 = nys_mg(grid_level-1)
2585
2586          DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2587             DO  i = nxl_mg(l), nxr_mg(l), 2
2588
2589                IF ( i == nxl_mg(l) )  THEN
2590                   !DIR$ IVDEP
2591                   DO  k = ind_even_odd+1, nzt_mg(l)
2592                      p_mg(k,j,nxr_mg(l)+1)  = temp(k-ind_even_odd,j1,ixr+1)
2593                   ENDDO
2594                   j1 = j1 + 1
2595
2596                ENDIF
2597
2598                IF ( i == nxr_mg(l) )  THEN
2599                   !DIR$ IVDEP
2600                   DO  k = ind_even_odd+1, nzt_mg(l)
2601                      p_mg(k,j,nxl_mg(l)-1)  = temp(k-ind_even_odd,j2,ixl-1)
2602                   ENDDO
2603                   j2 = j2 + 1
2604
2605                ENDIF
2606
2607             ENDDO
2608          ENDDO
2609
2610          DO j = nys_mg(l) + (color-1), nyn_mg(l), 2
2611             DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2612
2613                IF ( i == nxl_mg(l) )  THEN
2614                   !DIR$ IVDEP
2615                   DO  k = ind_even_odd+1, nzt_mg(l)
2616                      p_mg(k,j,nxr_mg(l)+1)  = temp(k-ind_even_odd,j1,ixr+1)
2617                   ENDDO
2618                   j1 = j1 + 1
2619
2620                ENDIF
2621
2622                IF ( i == nxr_mg(l) )  THEN
2623                   !DIR$ IVDEP
2624                   DO  k = ind_even_odd+1, nzt_mg(l)
2625                      p_mg(k,j,nxl_mg(l)-1)  = temp(k-ind_even_odd,j2,ixl-1)
2626                   ENDDO
2627                   j2 = j2 + 1
2628
2629                ENDIF
2630
2631             ENDDO
2632          ENDDO
2633
2634!
2635!--       Now handling the even k values
2636!--       Collecting data for the north - south exchange
2637!--       Since only every second value has to be transfered, data are stored
2638!--       on the next coarser grid level, because the arrays on that level
2639!--       have just the required size
2640          i1 = nxl_mg(grid_level-1)
2641          i2 = nxl_mg(grid_level-1)
2642
2643          DO  i = nxl_mg(l), nxr_mg(l), 2
2644             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2645
2646                IF ( j == nys_mg(l) )  THEN
2647                   !DIR$ IVDEP
2648                   DO  k = nzb+1, ind_even_odd
2649                      temp(k,jys,i1) = p_mg(k,j,i)
2650                   ENDDO
2651                   i1 = i1 + 1
2652
2653                ENDIF
2654
2655                IF ( j == nyn_mg(l) )  THEN
2656                   !DIR$ IVDEP
2657                   DO  k = nzb+1, ind_even_odd
2658                      temp(k,jyn,i2) = p_mg(k,j,i)
2659                   ENDDO
2660                   i2 = i2 + 1
2661
2662                ENDIF
2663
2664             ENDDO
2665          ENDDO
2666
2667          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2668             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2669
2670                IF ( j == nys_mg(l) )  THEN
2671                   !DIR$ IVDEP
2672                   DO  k = nzb+1, ind_even_odd
2673                      temp(k,jys,i1) = p_mg(k,j,i)
2674                   ENDDO
2675                   i1 = i1 + 1
2676
2677                ENDIF
2678
2679                IF ( j == nyn_mg(l) )  THEN
2680                   !DIR$ IVDEP
2681                   DO  k = nzb+1, ind_even_odd
2682                      temp(k,jyn,i2) = p_mg(k,j,i)
2683                   ENDDO
2684                   i2 = i2 + 1
2685
2686                ENDIF
2687
2688             ENDDO
2689          ENDDO
2690
2691          grid_level = grid_level-1
2692
2693          send_receive = 'ns'
2694          CALL exchange_horiz( temp, 1 )
2695
2696          grid_level = grid_level+1
2697
2698          i1 = nxl_mg(grid_level-1)
2699          i2 = nxl_mg(grid_level-1)
2700
2701          DO  i = nxl_mg(l), nxr_mg(l), 2
2702             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2703
2704                IF ( j == nys_mg(l) )  THEN
2705                   !DIR$ IVDEP
2706                   DO  k = nzb+1, ind_even_odd
2707                      p_mg(k,nyn_mg(l)+1,i) = temp(k,jyn+1,i1)
2708                   ENDDO
2709                   i1 = i1 + 1
2710
2711                ENDIF
2712
2713                IF ( j == nyn_mg(l) )  THEN
2714                   !DIR$ IVDEP
2715                   DO  k = nzb+1, ind_even_odd
2716                      p_mg(k,nys_mg(l)-1,i) = temp(k,jys-1,i2)
2717                   ENDDO
2718                   i2 = i2 + 1
2719
2720                ENDIF
2721
2722             ENDDO
2723          ENDDO
2724
2725          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2726             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2727
2728                IF ( j == nys_mg(l) )  THEN
2729                   !DIR$ IVDEP
2730                   DO  k = nzb+1, ind_even_odd
2731                      p_mg(k,nyn_mg(l)+1,i) = temp(k,jyn+1,i1)
2732                   ENDDO
2733                   i1 = i1 + 1
2734
2735                ENDIF
2736
2737                IF ( j == nyn_mg(l) )  THEN
2738                   !DIR$ IVDEP
2739                   DO  k = nzb+1, ind_even_odd
2740                      p_mg(k,nys_mg(l)-1,i) = temp(k,jys-1,i2)
2741                   ENDDO
2742                   i2 = i2 + 1
2743
2744                ENDIF
2745
2746             ENDDO
2747          ENDDO
2748
2749          j1 = nys_mg(grid_level-1)
2750          j2 = nys_mg(grid_level-1)
2751
2752          DO  i = nxl_mg(l), nxr_mg(l), 2
2753             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2754
2755                IF ( i == nxl_mg(l) )  THEN
2756                   !DIR$ IVDEP
2757                   DO  k = nzb+1, ind_even_odd
2758                      temp(k,j1,ixl) = p_mg(k,j,i)
2759                   ENDDO
2760                   j1 = j1 + 1
2761
2762                ENDIF
2763
2764                IF ( i == nxr_mg(l) )  THEN
2765                   !DIR$ IVDEP
2766                   DO  k = nzb+1, ind_even_odd
2767                      temp(k,j2,ixr) = p_mg(k,j,i)
2768                   ENDDO
2769                   j2 = j2 + 1
2770
2771                ENDIF
2772
2773             ENDDO
2774          ENDDO
2775
2776          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2777             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2778
2779                IF ( i == nxl_mg(l) )  THEN
2780                   !DIR$ IVDEP
2781                   DO  k = nzb+1, ind_even_odd
2782                      temp(k,j1,ixl) = p_mg(k,j,i)
2783                   ENDDO
2784                   j1 = j1 + 1
2785
2786                ENDIF
2787
2788                IF ( i == nxr_mg(l) )  THEN
2789                   !DIR$ IVDEP
2790                   DO  k = nzb+1, ind_even_odd
2791                      temp(k,j2,ixr) = p_mg(k,j,i)
2792                   ENDDO
2793                   j2 = j2 + 1
2794
2795                ENDIF
2796
2797             ENDDO
2798          ENDDO
2799
2800          grid_level = grid_level-1
2801
2802          send_receive = 'lr'
2803          CALL exchange_horiz( temp, 1 )
2804
2805          grid_level = grid_level+1
2806
2807          nxl = nxl_mg(grid_level)
2808          nys = nys_mg(grid_level)
2809          nxr = nxr_mg(grid_level)
2810          nyn = nyn_mg(grid_level)
2811          nzt = nzt_mg(grid_level)
2812
2813          j1 = nys_mg(grid_level-1)
2814          j2 = nys_mg(grid_level-1)
2815
2816          DO  i = nxl_mg(l), nxr_mg(l), 2
2817             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2818
2819                IF ( i == nxl_mg(l) )  THEN
2820                   !DIR$ IVDEP
2821                   DO  k = nzb+1, ind_even_odd
2822                      p_mg(k,j,nxr_mg(l)+1)  = temp(k,j1,ixr+1)
2823                   ENDDO
2824                   j1 = j1 + 1
2825
2826                ENDIF
2827
2828                IF ( i == nxr_mg(l) )  THEN
2829                   !DIR$ IVDEP
2830                   DO  k = nzb+1, ind_even_odd
2831                      p_mg(k,j,nxl_mg(l)-1)  = temp(k,j2,ixl-1)
2832                   ENDDO
2833                   j2 = j2 + 1
2834
2835                ENDIF
2836
2837             ENDDO
2838          ENDDO
2839
2840          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2841             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2842
2843                IF ( i == nxl_mg(l) )  THEN
2844                   !DIR$ IVDEP
2845                   DO  k = nzb+1, ind_even_odd
2846                      p_mg(k,j,nxr_mg(l)+1)  = temp(k,j1,ixr+1)
2847                   ENDDO
2848                   j1 = j1 + 1
2849
2850                ENDIF
2851
2852                IF ( i == nxr_mg(l) )  THEN
2853                   !DIR$ IVDEP
2854                   DO  k = nzb+1, ind_even_odd
2855                      p_mg(k,j,nxl_mg(l)-1)  = temp(k,j2,ixl-1)
2856                   ENDDO
2857                   j2 = j2 + 1
2858
2859                ENDIF
2860
2861             ENDDO
2862          ENDDO
2863
2864       ELSE
2865!
2866!--       Standard horizontal ghost boundary exchange for small coarse grid
2867!--       levels, where the transfer time is latency bound
2868          CALL exchange_horiz( p_mg, 1 )
2869
2870       ENDIF
2871
2872!
2873!--    Reset values to default PALM setup
2874       sendrecv_in_background = sendrecv_in_background_save
2875       synchronous_exchange   = synchronous_exchange_save
2876       send_receive = 'al'
2877
2878    END SUBROUTINE special_exchange_horiz
2879
2880 END MODULE poismg_mod
Note: See TracBrowser for help on using the repository browser.