source: palm/trunk/SOURCE/poismg_mod.f90 @ 4431

Last change on this file since 4431 was 4429, checked in by raasch, 4 years ago

serial (non-MPI) test case added, several bugfixes for the serial mode

  • Property svn:keywords set to Id
File size: 85.6 KB
Line 
1!> @file poismg.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
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-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: poismg_mod.f90 4429 2020-02-27 15:24:30Z gronemeier $
27! statement added to avoid compile error due to unused dummy argument
28! bugfix: cpp-directives added for serial mode
29!
30! 4360 2020-01-07 11:25:50Z suehring
31! Corrected "Former revisions" section
32!
33! 3725 2019-02-07 10:11:02Z raasch
34! unused subroutine removed
35!
36! 3655 2019-01-07 16:51:22Z knoop
37! unnecessary check eliminated
38!
39! Following optimisations have been made:
40! - vectorisation (for Intel-CPUs) of the red-black algorithm by resorting
41!   array elements with even and odd indices
42! - explicit boundary conditions for building walls removed (solver is
43!   running through the buildings
44! - reduced data transfer in case of ghost point exchange, because only
45!   "red" or "black" data points need to be exchanged. This is not applied
46!   for coarser grid levels, since for then the transfer time is latency bound
47!
48!
49! Description:
50! ------------
51!> Solves the Poisson equation for the perturbation pressure with a multigrid
52!> V- or W-Cycle scheme.
53!>
54!> This multigrid method was originally developed for PALM by Joerg Uhlenbrock,
55!> September 2000 - July 2001. It has been optimised for speed by Klaus
56!> Ketelsen in November 2014.
57!>
58!> @attention Loop unrolling and cache optimization in SOR-Red/Black method
59!>            still does not give the expected speedup!
60!>
61!> @todo Further work required.
62!------------------------------------------------------------------------------!
63 MODULE poismg_mod
64 
65    USE control_parameters,                                                    &
66        ONLY:  bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r,                 &
67               bc_dirichlet_s, bc_radiation_l, bc_radiation_n, bc_radiation_r, &
68               bc_radiation_s, grid_level, nesting_offline
69
70    USE cpulog,                                                                &
71        ONLY:  cpu_log, log_point_s
72
73    USE kinds
74
75    USE pegrid
76
77    PRIVATE
78
79    INTEGER, SAVE                             ::  ind_even_odd    !< border index between even and odd k index
80    INTEGER, DIMENSION(:), SAVE, ALLOCATABLE  ::  even_odd_level  !< stores ind_even_odd for all MG levels
81
82    REAL(wp), DIMENSION(:,:), SAVE, ALLOCATABLE ::  f1_mg_b, f2_mg_b, f3_mg_b  !< blocked version of f1_mg ...
83
84    INTERFACE poismg
85       MODULE PROCEDURE poismg
86    END INTERFACE poismg
87
88    INTERFACE sort_k_to_even_odd_blocks
89       MODULE PROCEDURE sort_k_to_even_odd_blocks
90!       MODULE PROCEDURE sort_k_to_even_odd_blocks_int
91       MODULE PROCEDURE sort_k_to_even_odd_blocks_1d
92    END INTERFACE sort_k_to_even_odd_blocks
93
94    PUBLIC poismg
95
96 CONTAINS
97
98!------------------------------------------------------------------------------!
99! Description:
100! ------------
101!> Solves the Poisson equation for the perturbation pressure with a multigrid
102!> V- or W-Cycle scheme.
103!------------------------------------------------------------------------------!
104    SUBROUTINE poismg( r )
105
106       USE arrays_3d,                                                          &
107           ONLY:  d, p_loc
108
109       USE control_parameters,                                                 &
110           ONLY:  bc_lr_cyc, bc_ns_cyc, gathered_size, grid_level,             &
111                  grid_level_count, ibc_p_t,                                   &
112                  maximum_grid_level, message_string, mgcycles, mg_cycles,     &
113                  mg_switch_to_pe0_level, residual_limit, subdomain_size
114
115       USE cpulog,                                                             &
116           ONLY:  cpu_log, log_point_s
117
118       USE indices,                                                            &
119           ONLY:  nxl, nxlg, nxl_mg, nxr, nxrg, nxr_mg, nys, nysg, nys_mg, nyn,&
120                  nyng, nyn_mg, nzb, nzt, nzt_mg
121
122       IMPLICIT NONE
123
124       REAL(wp) ::  maxerror          !<
125       REAL(wp) ::  maximum_mgcycles  !<
126       REAL(wp) ::  residual_norm     !<
127
128       REAL(wp), DIMENSION(nzb:nzt+1,nys-1:nyn+1,nxl-1:nxr+1) ::  r  !<
129
130       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  p3  !<
131
132
133       CALL cpu_log( log_point_s(29), 'poismg', 'start' )
134!
135!--    Initialize arrays and variables used in this subroutine
136
137!--    If the number of grid points of the gathered grid, which is collected
138!--    on PE0, is larger than the number of grid points of an PE, than array
139!--    p3 will be enlarged.
140       IF ( gathered_size > subdomain_size )  THEN
141          ALLOCATE( p3(nzb:nzt_mg(mg_switch_to_pe0_level)+1,nys_mg(            &
142                    mg_switch_to_pe0_level)-1:nyn_mg(mg_switch_to_pe0_level)+1,&
143                    nxl_mg(mg_switch_to_pe0_level)-1:nxr_mg(                   &
144                    mg_switch_to_pe0_level)+1) )
145       ELSE
146          ALLOCATE ( p3(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
147       ENDIF
148
149       p3 = 0.0_wp
150
151 
152!
153!--    Ghost boundaries have to be added to divergence array.
154!--    Exchange routine needs to know the grid level!
155       grid_level = maximum_grid_level
156       CALL exchange_horiz( d, 1)
157!
158!--    Set bottom and top boundary conditions
159       d(nzb,:,:) = d(nzb+1,:,:)
160       IF ( ibc_p_t == 1 )  d(nzt+1,:,: ) = d(nzt,:,:)
161!
162!--    Set lateral boundary conditions in non-cyclic case
163       IF ( .NOT. bc_lr_cyc )  THEN
164          IF ( bc_dirichlet_l  .OR.  bc_radiation_l )                          &
165             d(:,:,nxl-1) = d(:,:,nxl)
166          IF ( bc_dirichlet_r  .OR.  bc_radiation_r )                          &
167             d(:,:,nxr+1) = d(:,:,nxr)
168       ENDIF
169       IF ( .NOT. bc_ns_cyc )  THEN
170          IF ( bc_dirichlet_n  .OR.  bc_radiation_n )                          &
171             d(:,nyn+1,:) = d(:,nyn,:)
172          IF ( bc_dirichlet_s  .OR.  bc_radiation_s )                          &
173             d(:,nys-1,:) = d(:,nys,:)
174       ENDIF
175!
176!--    Initiation of the multigrid scheme. Does n cycles until the
177!--    residual is smaller than the given limit. The accuracy of the solution
178!--    of the poisson equation will increase with the number of cycles.
179!--    If the number of cycles is preset by the user, this number will be
180!--    carried out regardless of the accuracy.
181       grid_level_count =  0
182       mgcycles         =  0
183       IF ( mg_cycles == -1 )  THEN
184          maximum_mgcycles = 0
185          residual_norm    = 1.0_wp
186       ELSE
187          maximum_mgcycles = mg_cycles
188          residual_norm    = 0.0_wp
189       ENDIF
190
191!
192!--    Initial settings for sorting k-dimension from sequential order (alternate
193!--    even/odd) into blocks of even and odd or vice versa
194       CALL init_even_odd_blocks
195
196!
197!--    Sort input arrays in even/odd blocks along k-dimension
198       CALL sort_k_to_even_odd_blocks( d, grid_level )
199       CALL sort_k_to_even_odd_blocks( p_loc, grid_level )
200
201!
202!--    The complete multigrid cycles are running in block mode, i.e. over
203!--    seperate data blocks of even and odd indices
204       DO WHILE ( residual_norm > residual_limit  .OR. &
205                  mgcycles < maximum_mgcycles )
206 
207          CALL next_mg_level( d, p_loc, p3, r)
208
209!
210!--       Calculate the residual if the user has not preset the number of
211!--       cycles to be performed
212          IF ( maximum_mgcycles == 0 )  THEN
213             CALL resid( d, p_loc, r )
214             maxerror = SUM( r(nzb+1:nzt,nys:nyn,nxl:nxr)**2 )
215
216#if defined( __parallel )
217             IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
218                CALL MPI_ALLREDUCE( maxerror, residual_norm, 1, MPI_REAL,      &
219                                    MPI_SUM, comm2d, ierr)
220#else
221                residual_norm = maxerror
222#endif
223             residual_norm = SQRT( residual_norm )
224          ENDIF
225
226          mgcycles = mgcycles + 1
227
228!
229!--       If the user has not limited the number of cycles, stop the run in case
230!--       of insufficient convergence
231          IF ( mgcycles > 1000  .AND.  mg_cycles == -1 )  THEN
232             message_string = 'no sufficient convergence within 1000 cycles'
233             CALL message( 'poismg', 'PA0283', 1, 2, 0, 6, 0 )
234          ENDIF
235
236       ENDDO
237
238       DEALLOCATE( p3 )
239!
240!--    Result has to be sorted back from even/odd blocks to sequential order
241       CALL sort_k_to_sequential( p_loc )
242!
243!--    Unset the grid level. Variable is used to determine the MPI datatypes for
244!--    ghost point exchange
245       grid_level = 0
246
247       CALL cpu_log( log_point_s(29), 'poismg', 'stop' )
248
249    END SUBROUTINE poismg
250
251
252!------------------------------------------------------------------------------!
253! Description:
254! ------------
255!> Computes the residual of the perturbation pressure.
256!------------------------------------------------------------------------------!
257    SUBROUTINE resid( f_mg, p_mg, r )
258
259
260       USE arrays_3d,                                                          &
261           ONLY:  rho_air_mg
262
263       USE control_parameters,                                                 &
264           ONLY:  bc_lr_cyc, bc_ns_cyc, ibc_p_b, ibc_p_t
265       USE grid_variables,                                                     &
266           ONLY:  ddx2_mg, ddy2_mg
267
268       USE indices,                                                            &
269           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
270
271       IMPLICIT NONE
272
273       INTEGER(iwp) ::  i        !< index variable along x
274       INTEGER(iwp) ::  j        !< index variable along y
275       INTEGER(iwp) ::  k        !< index variable along z
276       INTEGER(iwp) ::  l        !< index indicating grid level
277       INTEGER(iwp) ::  km1      !< index variable along z dimension (k-1)
278       INTEGER(iwp) ::  kp1      !< index variable along z dimension (k+1)
279
280       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
281                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
282                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f_mg  !< velocity divergence
283       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
284                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
285                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  p_mg  !< perturbation pressure
286       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
287                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
288                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  r     !< residuum of perturbation pressure
289
290!
291!--    Calculate the residual
292       l = grid_level
293
294       CALL cpu_log( log_point_s(53), 'resid', 'start' )
295       !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
296       !$OMP DO
297       DO  i = nxl_mg(l), nxr_mg(l)
298          DO  j = nys_mg(l), nyn_mg(l)
299                !DIR$ IVDEP
300             DO k = ind_even_odd+1, nzt_mg(l)
301                km1 = k-ind_even_odd-1
302                kp1 = k-ind_even_odd
303                r(k,j,i) = f_mg(k,j,i)                                         &
304                      - rho_air_mg(k,l) * ddx2_mg(l) *                         &
305                      ( p_mg(k,j,i+1) +  p_mg(k,j,i-1)  )                      &
306                      - rho_air_mg(k,l) * ddy2_mg(l) *                         &
307                      ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )                       &
308                      - f2_mg_b(k,l) * p_mg(kp1,j,i)                           &
309                      - f3_mg_b(k,l) * p_mg(km1,j,i)                           &
310                      + f1_mg_b(k,l) * p_mg(k,j,i)
311             ENDDO
312             !DIR$ IVDEP
313             DO k = nzb+1, ind_even_odd
314                km1 = k+ind_even_odd
315                kp1 = k+ind_even_odd+1
316                r(k,j,i) = f_mg(k,j,i)                                         &
317                      - rho_air_mg(k,l) * ddx2_mg(l) *                         &
318                      ( p_mg(k,j,i+1) +  p_mg(k,j,i-1)  )                      &
319                      - rho_air_mg(k,l) * ddy2_mg(l) *                         &
320                      ( p_mg(k,j+1,i) + p_mg(k,j-1,i)  )                       &
321                      - f2_mg_b(k,l) * p_mg(kp1,j,i)                           &
322                      - f3_mg_b(k,l) * p_mg(km1,j,i)                           &
323                      + f1_mg_b(k,l) * p_mg(k,j,i)
324             ENDDO
325          ENDDO
326       ENDDO
327       !$OMP END PARALLEL
328!
329!--    Horizontal boundary conditions
330       CALL exchange_horiz( r, 1)
331
332       IF ( .NOT. bc_lr_cyc )  THEN
333          IF ( bc_dirichlet_l  .OR.  bc_radiation_l )  THEN
334             r(:,:,nxl_mg(l)-1) = r(:,:,nxl_mg(l))
335          ENDIF
336          IF ( bc_dirichlet_r  .OR.  bc_radiation_r )  THEN
337             r(:,:,nxr_mg(l)+1) = r(:,:,nxr_mg(l))
338          ENDIF
339       ENDIF
340
341       IF ( .NOT. bc_ns_cyc )  THEN
342          IF ( bc_dirichlet_n  .OR.  bc_radiation_n )  THEN
343             r(:,nyn_mg(l)+1,:) = r(:,nyn_mg(l),:)
344          ENDIF
345          IF ( bc_dirichlet_s  .OR.  bc_radiation_s )  THEN
346             r(:,nys_mg(l)-1,:) = r(:,nys_mg(l),:)
347          ENDIF
348       ENDIF
349
350!
351!--    Boundary conditions at bottom and top of the domain. Points may be within
352!--    buildings, but that doesn't matter.
353       IF ( ibc_p_b == 1 )  THEN
354!
355!--       equivalent to r(nzb,:,: ) = r(nzb+1,:,:)
356          r(nzb,:,: ) = r(ind_even_odd+1,:,:)
357       ELSE
358          r(nzb,:,: ) = 0.0_wp
359       ENDIF
360
361       IF ( ibc_p_t == 1 )  THEN
362!
363!--       equivalent to r(nzt_mg(l)+1,:,: ) = r(nzt_mg(l),:,:)
364          r(nzt_mg(l)+1,:,: ) = r(ind_even_odd,:,:)
365       ELSE
366          r(nzt_mg(l)+1,:,: ) = 0.0_wp
367       ENDIF
368
369       CALL cpu_log( log_point_s(53), 'resid', 'stop' )
370
371    END SUBROUTINE resid
372
373
374!------------------------------------------------------------------------------!
375! Description:
376! ------------
377!> Interpolates the residual on the next coarser grid with "full weighting"
378!> scheme
379!------------------------------------------------------------------------------!
380    SUBROUTINE restrict( f_mg, r )
381
382
383       USE control_parameters,                                                 &
384           ONLY:  bc_lr_cyc, bc_ns_cyc, ibc_p_b, ibc_p_t
385
386       USE indices,                                                            &
387           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
388
389       IMPLICIT NONE
390
391       INTEGER(iwp) ::  i    !< index variable along x on finer grid
392       INTEGER(iwp) ::  ic   !< index variable along x on coarser grid
393       INTEGER(iwp) ::  j    !< index variable along y on finer grid
394       INTEGER(iwp) ::  jc   !< index variable along y on coarser grid
395       INTEGER(iwp) ::  k    !< index variable along z on finer grid
396       INTEGER(iwp) ::  kc   !< index variable along z on coarser grid
397       INTEGER(iwp) ::  l    !< index indicating finer grid level
398       INTEGER(iwp) ::  km1  !< index variable along z dimension (k-1 on finer level)
399       INTEGER(iwp) ::  kp1  !< index variable along z dimension (k+1 on finer level)
400
401
402       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
403                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
404                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
405                                         f_mg  !< Residual on coarser grid level
406
407       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level+1)+1,                         &
408                           nys_mg(grid_level+1)-1:nyn_mg(grid_level+1)+1,      &
409                           nxl_mg(grid_level+1)-1:nxr_mg(grid_level+1)+1) ::   &
410                                         r !< Residual on finer grid level
411
412!
413!--    Interpolate the residual
414       l = grid_level
415
416       CALL cpu_log( log_point_s(54), 'restrict', 'start' )
417!
418!--    No wall treatment
419       !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,kc,km1,kp1)
420       !$OMP DO SCHEDULE( STATIC )
421       DO  ic = nxl_mg(l), nxr_mg(l)
422          i = 2*ic
423          DO  jc = nys_mg(l), nyn_mg(l)
424!
425!--          Calculation for the first point along k
426             j  = 2*jc
427!
428!--          Calculation for the other points along k
429             !DIR$ IVDEP
430             DO k = ind_even_odd+1, nzt_mg(l+1)    ! Fine grid at this point
431                km1 = k-ind_even_odd-1
432                kp1 = k-ind_even_odd
433                kc  = k-ind_even_odd               ! Coarse grid index
434
435                f_mg(kc,jc,ic) = 1.0_wp / 64.0_wp * (                      &
436                               8.0_wp * r(k,j,i)                            &
437                             + 4.0_wp * ( r(k,j,i-1)     + r(k,j,i+1)     + &
438                                          r(k,j+1,i)     + r(k,j-1,i)     ) &
439                             + 2.0_wp * ( r(k,j-1,i-1)   + r(k,j+1,i-1)   + &
440                                          r(k,j-1,i+1)   + r(k,j+1,i+1)   ) &
441                             + 4.0_wp * r(km1,j,i)                          &
442                             + 2.0_wp * ( r(km1,j,i-1)   + r(km1,j,i+1)   + &
443                                          r(km1,j+1,i)   + r(km1,j-1,i)   ) &
444                             +          ( r(km1,j-1,i-1) + r(km1,j+1,i-1) + &
445                                          r(km1,j-1,i+1) + r(km1,j+1,i+1) ) &
446                             + 4.0_wp * r(kp1,j,i)                          &
447                             + 2.0_wp * ( r(kp1,j,i-1)   + r(kp1,j,i+1)   + &
448                                          r(kp1,j+1,i)   + r(kp1,j-1,i)   ) &
449                             +          ( r(kp1,j-1,i-1) + r(kp1,j+1,i-1) + &
450                                          r(kp1,j-1,i+1) + r(kp1,j+1,i+1) ) &
451                                        )
452             ENDDO
453          ENDDO
454       ENDDO
455       !$OMP ENDDO
456       !$OMP END PARALLEL
457
458!
459!--    Ghost point exchange
460       CALL exchange_horiz( f_mg, 1)
461!
462!--    Horizontal boundary conditions
463       IF ( .NOT. bc_lr_cyc )  THEN
464          IF ( bc_dirichlet_l  .OR.  bc_radiation_l )  THEN
465             f_mg(:,:,nxl_mg(l)-1) = f_mg(:,:,nxl_mg(l))
466          ENDIF
467          IF ( bc_dirichlet_r  .OR.  bc_radiation_r )  THEN
468             f_mg(:,:,nxr_mg(l)+1) = f_mg(:,:,nxr_mg(l))
469          ENDIF
470       ENDIF
471
472       IF ( .NOT. bc_ns_cyc )  THEN
473          IF ( bc_dirichlet_n  .OR.  bc_radiation_n )  THEN
474             f_mg(:,nyn_mg(l)+1,:) = f_mg(:,nyn_mg(l),:)
475          ENDIF
476          IF ( bc_dirichlet_s  .OR.  bc_radiation_s )  THEN
477             f_mg(:,nys_mg(l)-1,:) = f_mg(:,nys_mg(l),:)
478          ENDIF
479       ENDIF
480
481!
482!--    Boundary conditions at bottom and top of the domain.
483!--    These points are not handled by the above loop. Points may be within
484!--    buildings, but that doesn't matter. Remark: f_mg is ordered sequentielly
485!--    after interpolation on coarse grid (is ordered in odd-even blocks further
486!--    below).
487       IF ( ibc_p_b == 1 )  THEN
488          f_mg(nzb,:,: ) = f_mg(nzb+1,:,:)
489       ELSE
490          f_mg(nzb,:,: ) = 0.0_wp
491       ENDIF
492
493       IF ( ibc_p_t == 1 )  THEN
494          f_mg(nzt_mg(l)+1,:,: ) = f_mg(nzt_mg(l),:,:)
495       ELSE
496          f_mg(nzt_mg(l)+1,:,: ) = 0.0_wp
497       ENDIF
498
499       CALL cpu_log( log_point_s(54), 'restrict', 'stop' )
500!
501!--    Since residual is in sequential order after interpolation, an additional
502!--    sorting in odd-even blocks along z dimension is required at this point.
503       CALL sort_k_to_even_odd_blocks( f_mg , l)
504
505    END SUBROUTINE restrict
506
507
508!------------------------------------------------------------------------------!
509! Description:
510! ------------
511!> Interpolates the correction of the perturbation pressure
512!> to the next finer grid.
513!------------------------------------------------------------------------------!
514    SUBROUTINE prolong( p, temp )
515
516
517       USE control_parameters,                                                 &
518           ONLY:  bc_lr_cyc, bc_ns_cyc, ibc_p_b, ibc_p_t
519       USE indices,                                                            &
520           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
521
522       IMPLICIT NONE
523
524       INTEGER(iwp) ::  i   !< index variable along x on coarser grid level
525       INTEGER(iwp) ::  j   !< index variable along y on coarser grid level
526       INTEGER(iwp) ::  k   !< index variable along z on coarser grid level
527       INTEGER(iwp) ::  l   !< index indicating finer grid level
528       INTEGER(iwp) ::  kp1 !< index variable along z
529       INTEGER(iwp) ::  ke  !< index for prolog even
530       INTEGER(iwp) ::  ko  !< index for prolog odd
531
532       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
533                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
534                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1 ) ::  &
535                               p     !< perturbation pressure on coarser grid level
536
537       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
538                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
539                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
540                               temp  !< perturbation pressure on finer grid level
541
542
543       CALL cpu_log( log_point_s(55), 'prolong', 'start' )
544
545!
546!--    First, store elements of the coarser grid on the next finer grid
547       l = grid_level
548       ind_even_odd = even_odd_level(grid_level-1)
549
550       !$OMP PARALLEL PRIVATE (i,j,k,kp1,ke,ko)
551       !$OMP DO
552       DO  i = nxl_mg(l-1), nxr_mg(l-1)
553          DO  j = nys_mg(l-1), nyn_mg(l-1)
554
555             !DIR$ IVDEP
556             DO k = ind_even_odd+1, nzt_mg(l-1)
557                kp1 = k - ind_even_odd
558                ke  = 2 * ( k-ind_even_odd - 1 ) + 1
559                ko  = 2 * k - 1
560!
561!--             Points of the coarse grid are directly stored on the next finer
562!--             grid
563                temp(ko,2*j,2*i)   = p(k,j,i)
564!
565!--             Points between two coarse-grid points
566                temp(ko,2*j,2*i+1) = 0.5_wp * ( p(k,j,i) + p(k,j,i+1) )
567                temp(ko,2*j+1,2*i) = 0.5_wp * ( p(k,j,i) + p(k,j+1,i) )
568                temp(ke,2*j,2*i)   = 0.5_wp * ( p(k,j,i) + p(kp1,j,i) )
569!
570!--             Points in the center of the planes stretched by four points
571!--             of the coarse grid cube
572                temp(ko,2*j+1,2*i+1) = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
573                                                   p(k,j+1,i) + p(k,j+1,i+1) )
574                temp(ke,2*j,2*i+1)   = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
575                                                   p(kp1,j,i) + p(kp1,j,i+1) )
576                temp(ke,2*j+1,2*i)   = 0.25_wp * ( p(k,j,i)   + p(k,j+1,i) +   &
577                                                   p(kp1,j,i) + p(kp1,j+1,i) )
578!
579!--             Points in the middle of coarse grid cube
580                temp(ke,2*j+1,2*i+1) = 0.125_wp *                              &
581                                               ( p(k,j,i)     + p(k,j,i+1)   + &
582                                                 p(k,j+1,i)   + p(k,j+1,i+1) + &
583                                                 p(kp1,j,i)   + p(kp1,j,i+1) + &
584                                                 p(kp1,j+1,i) + p(kp1,j+1,i+1) )
585
586             ENDDO
587
588             !DIR$ IVDEP
589             DO k = nzb+1, ind_even_odd
590                kp1 = k + ind_even_odd + 1
591                ke  = 2 * k
592                ko  = 2 * ( k + ind_even_odd )
593!
594!--             Points of the coarse grid are directly stored on the next finer
595!--             grid
596                temp(ko,2*j,2*i)   = p(k,j,i)
597!
598!--             Points between two coarse-grid points
599                temp(ko,2*j,2*i+1) = 0.5_wp * ( p(k,j,i) + p(k,j,i+1) )
600                temp(ko,2*j+1,2*i) = 0.5_wp * ( p(k,j,i) + p(k,j+1,i) )
601                temp(ke,2*j,2*i)   = 0.5_wp * ( p(k,j,i) + p(kp1,j,i) )
602!
603!--             Points in the center of the planes stretched by four points
604!--             of the coarse grid cube
605                temp(ko,2*j+1,2*i+1) = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
606                                                   p(k,j+1,i) + p(k,j+1,i+1) )
607                temp(ke,2*j,2*i+1)   = 0.25_wp * ( p(k,j,i)   + p(k,j,i+1) +   &
608                                                   p(kp1,j,i) + p(kp1,j,i+1) )
609                temp(ke,2*j+1,2*i)   = 0.25_wp * ( p(k,j,i)   + p(k,j+1,i) +   &
610                                                   p(kp1,j,i) + p(kp1,j+1,i) )
611!
612!--             Points in the middle of coarse grid cube
613                temp(ke,2*j+1,2*i+1) = 0.125_wp *                              &
614                                               ( p(k,j,i)     + p(k,j,i+1)   + &
615                                                 p(k,j+1,i)   + p(k,j+1,i+1) + &
616                                                 p(kp1,j,i)   + p(kp1,j,i+1) + &
617                                                 p(kp1,j+1,i) + p(kp1,j+1,i+1) )
618
619             ENDDO
620
621          ENDDO
622       ENDDO
623       !$OMP END PARALLEL
624
625       ind_even_odd = even_odd_level(grid_level)
626!
627!--    Horizontal boundary conditions
628       CALL exchange_horiz( temp, 1)
629
630       IF ( .NOT. bc_lr_cyc )  THEN
631          IF ( bc_dirichlet_l  .OR.  bc_radiation_l )  THEN
632             temp(:,:,nxl_mg(l)-1) = temp(:,:,nxl_mg(l))
633          ENDIF
634          IF ( bc_dirichlet_r  .OR.  bc_radiation_r )  THEN
635             temp(:,:,nxr_mg(l)+1) = temp(:,:,nxr_mg(l))
636          ENDIF
637       ENDIF
638
639       IF ( .NOT. bc_ns_cyc )  THEN
640          IF ( bc_dirichlet_n  .OR.  bc_radiation_n )  THEN
641             temp(:,nyn_mg(l)+1,:) = temp(:,nyn_mg(l),:)
642          ENDIF
643          IF ( bc_dirichlet_s  .OR.  bc_radiation_s )  THEN
644             temp(:,nys_mg(l)-1,:) = temp(:,nys_mg(l),:)
645          ENDIF
646       ENDIF
647
648!
649!--    Bottom and top boundary conditions
650       IF ( ibc_p_b == 1 )  THEN
651!
652!--       equivalent to temp(nzb,:,: ) = temp(nzb+1,:,:)
653          temp(nzb,:,: ) = temp(ind_even_odd+1,:,:)
654       ELSE
655          temp(nzb,:,: ) = 0.0_wp
656       ENDIF
657
658       IF ( ibc_p_t == 1 )  THEN
659!
660!--       equivalent to temp(nzt_mg(l)+1,:,: ) = temp(nzt_mg(l),:,:)
661          temp(nzt_mg(l)+1,:,: ) = temp(ind_even_odd,:,:)
662       ELSE
663          temp(nzt_mg(l)+1,:,: ) = 0.0_wp
664       ENDIF
665
666       CALL cpu_log( log_point_s(55), 'prolong', 'stop' )
667
668    END SUBROUTINE prolong
669
670
671!------------------------------------------------------------------------------!
672! Description:
673! ------------
674!> Relaxation method for the multigrid scheme. A Gauss-Seidel iteration with
675!> 3D-Red-Black decomposition (GS-RB) is used.
676!------------------------------------------------------------------------------!
677    SUBROUTINE redblack( f_mg, p_mg )
678
679
680       USE arrays_3d,                                                          &
681           ONLY:  rho_air_mg
682
683       USE control_parameters,                                                 &
684           ONLY:  bc_lr_cyc, bc_ns_cyc, ibc_p_b, ibc_p_t, ngsrb 
685
686       USE grid_variables,                                                     &
687           ONLY:  ddx2_mg, ddy2_mg
688
689       USE indices,                                                            &
690           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
691
692       IMPLICIT NONE
693
694       INTEGER(iwp) :: color    !< grid point color, either red or black
695       INTEGER(iwp) :: i        !< index variable along x
696       INTEGER(iwp) :: ic       !< index variable along x
697       INTEGER(iwp) :: j        !< index variable along y
698       INTEGER(iwp) :: jc       !< index variable along y
699       INTEGER(iwp) :: jj       !< index variable along y
700       INTEGER(iwp) :: k        !< index variable along z
701       INTEGER(iwp) :: l        !< grid level
702       INTEGER(iwp) :: n        !< loop variable Gauß-Seidel iterations
703       INTEGER(iwp) :: km1      !< index variable (k-1)
704       INTEGER(iwp) :: kp1      !< index variable (k+1)
705
706       LOGICAL      :: unroll   !< flag indicating whether loop unrolling is possible
707
708       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
709                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
710                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
711                                      f_mg  !< residual of perturbation pressure
712       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
713                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
714                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
715                                      p_mg  !< perturbation pressure
716
717       l = grid_level
718
719       unroll = ( MOD( nyn_mg(l)-nys_mg(l)+1, 4 ) == 0  .AND.                  &
720                  MOD( nxr_mg(l)-nxl_mg(l)+1, 2 ) == 0 )
721
722       DO  n = 1, ngsrb
723       
724          DO  color = 1, 2
725
726             IF ( .NOT. unroll )  THEN
727
728                CALL cpu_log( log_point_s(36), 'redblack_no_unroll_f', 'start' )
729!
730!--             Without unrolling of loops, no cache optimization
731                !$OMP PARALLEL PRIVATE (i,j,k,km1,kp1)
732                !$OMP DO
733                DO  i = nxl_mg(l), nxr_mg(l), 2
734                   DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
735                      !DIR$ IVDEP
736                      DO  k = ind_even_odd+1, nzt_mg(l)
737                         km1 = k-ind_even_odd-1
738                         kp1 = k-ind_even_odd
739                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
740                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
741                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
742                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
743                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
744                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
745                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
746                               - f_mg(k,j,i)                   )
747                      ENDDO
748                   ENDDO
749                ENDDO
750   
751                !$OMP DO
752                DO  i = nxl_mg(l)+1, nxr_mg(l), 2
753                   DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
754                       !DIR$ IVDEP
755                       DO  k = ind_even_odd+1, nzt_mg(l)
756                         km1 = k-ind_even_odd-1
757                         kp1 = k-ind_even_odd
758                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
759                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
760                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
761                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
762                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
763                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
764                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
765                               - f_mg(k,j,i)                   )
766                      ENDDO
767                   ENDDO
768                ENDDO
769 
770                !$OMP DO
771                DO  i = nxl_mg(l), nxr_mg(l), 2
772                   DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
773                      !DIR$ IVDEP
774                      DO  k = nzb+1, ind_even_odd
775                         km1 = k+ind_even_odd
776                         kp1 = k+ind_even_odd+1
777                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
778                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
779                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
780                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
781                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
782                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
783                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
784                               - f_mg(k,j,i)                   )
785                      ENDDO
786                   ENDDO
787                ENDDO
788
789                !$OMP DO
790                DO  i = nxl_mg(l)+1, nxr_mg(l), 2
791                   DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
792                      !DIR$ IVDEP
793                      DO  k = nzb+1, ind_even_odd
794                         km1 = k+ind_even_odd
795                         kp1 = k+ind_even_odd+1
796                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
797                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
798                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
799                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
800                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
801                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
802                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
803                               - f_mg(k,j,i)                   )
804                      ENDDO
805                   ENDDO
806                ENDDO
807                !$OMP END PARALLEL
808
809                CALL cpu_log( log_point_s(36), 'redblack_no_unroll_f', 'stop' )
810
811             ELSE
812!
813!--              Loop unrolling along y, only one i loop for better cache use
814                CALL cpu_log( log_point_s(38), 'redblack_unroll_f', 'start' )
815
816                !$OMP PARALLEL PRIVATE (i,j,k,ic,jc,km1,kp1,jj)
817                !$OMP DO
818                DO  ic = nxl_mg(l), nxr_mg(l), 2
819                   DO  jc = nys_mg(l), nyn_mg(l), 4
820                      i  = ic
821                      jj = jc+2-color
822                      !DIR$ IVDEP
823                      DO  k = ind_even_odd+1, nzt_mg(l)
824                         km1 = k-ind_even_odd-1
825                         kp1 = k-ind_even_odd
826                         j   = jj
827                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
828                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
829                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
830                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
831                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
832                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
833                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
834                               - f_mg(k,j,i)                   )
835                         j = jj+2
836                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
837                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
838                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
839                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
840                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
841                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
842                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
843                               - f_mg(k,j,i)                   )
844                      ENDDO
845
846                      i  = ic+1
847                      jj = jc+color-1
848                      !DIR$ IVDEP
849                      DO  k = ind_even_odd+1, nzt_mg(l)
850                         km1 = k-ind_even_odd-1
851                         kp1 = k-ind_even_odd
852                         j   = jj
853                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
854                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
855                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
856                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
857                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
858                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
859                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
860                               - f_mg(k,j,i)                   )
861                         j = jj+2
862                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
863                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
864                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
865                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
866                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
867                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
868                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
869                               - f_mg(k,j,i)                   )
870                      ENDDO
871
872                      i  = ic
873                      jj = jc+color-1
874                      !DIR$ IVDEP
875                      DO  k = nzb+1, ind_even_odd
876                         km1 = k+ind_even_odd
877                         kp1 = k+ind_even_odd+1
878                         j   = jj
879                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
880                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
881                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
882                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
883                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
884                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
885                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
886                               - f_mg(k,j,i)                   )
887                         j = jj+2
888                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
889                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
890                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
891                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
892                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
893                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
894                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
895                               - f_mg(k,j,i)                   )
896                      ENDDO
897
898                      i  = ic+1
899                      jj = jc+2-color
900                      !DIR$ IVDEP
901                      DO  k = nzb+1, ind_even_odd
902                         km1 = k+ind_even_odd
903                         kp1 = k+ind_even_odd+1
904                         j   = jj
905                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
906                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
907                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
908                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
909                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
910                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
911                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
912                               - f_mg(k,j,i)                   )
913                         j = jj+2
914                         p_mg(k,j,i) = 1.0_wp / f1_mg_b(k,l) * (               &
915                                 rho_air_mg(k,l) * ddx2_mg(l) *                &
916                               ( p_mg(k,j,i+1) + p_mg(k,j,i-1) )               &
917                               + rho_air_mg(k,l) * ddy2_mg(l) *                &
918                               ( p_mg(k,j+1,i) + p_mg(k,j-1,i) )               &
919                               + f2_mg_b(k,l) * p_mg(kp1,j,i)                  &
920                               + f3_mg_b(k,l) * p_mg(km1,j,i)                  &
921                               - f_mg(k,j,i)                   )
922                      ENDDO
923
924                   ENDDO
925                ENDDO
926                !$OMP END PARALLEL
927
928                CALL cpu_log( log_point_s(38), 'redblack_unroll_f', 'stop' )
929
930             ENDIF
931
932!
933!--          Horizontal boundary conditions
934             CALL special_exchange_horiz( p_mg, color )
935
936             IF ( .NOT. bc_lr_cyc )  THEN
937                IF ( bc_dirichlet_l  .OR.  bc_radiation_l )  THEN
938                   p_mg(:,:,nxl_mg(l)-1) = p_mg(:,:,nxl_mg(l))
939                ENDIF
940                IF ( bc_dirichlet_r  .OR.  bc_radiation_r )  THEN
941                   p_mg(:,:,nxr_mg(l)+1) = p_mg(:,:,nxr_mg(l))
942                ENDIF
943             ENDIF
944
945             IF ( .NOT. bc_ns_cyc )  THEN
946                IF ( bc_dirichlet_n  .OR.  bc_radiation_n )  THEN
947                   p_mg(:,nyn_mg(l)+1,:) = p_mg(:,nyn_mg(l),:)
948                ENDIF
949                IF ( bc_dirichlet_s  .OR.  bc_radiation_s )  THEN
950                   p_mg(:,nys_mg(l)-1,:) = p_mg(:,nys_mg(l),:)
951                ENDIF
952             ENDIF
953
954!
955!--          Bottom and top boundary conditions
956             IF ( ibc_p_b == 1 )  THEN
957!
958!--             equivalent to p_mg(nzb,:,: ) = p_mg(nzb+1,:,:)
959                p_mg(nzb,:,: ) = p_mg(ind_even_odd+1,:,:)
960             ELSE
961                p_mg(nzb,:,: ) = 0.0_wp
962             ENDIF
963
964             IF ( ibc_p_t == 1 )  THEN
965!
966!--             equivalent to p_mg(nzt_mg(l)+1,:,: ) = p_mg(nzt_mg(l),:,:)
967                p_mg(nzt_mg(l)+1,:,: ) = p_mg(ind_even_odd,:,:)
968             ELSE
969                p_mg(nzt_mg(l)+1,:,: ) = 0.0_wp
970             ENDIF
971
972          ENDDO
973
974       ENDDO
975
976    END SUBROUTINE redblack
977
978
979!------------------------------------------------------------------------------!
980! Description:
981! ------------
982!> Sort k-Dimension from sequential into blocks of even and odd.
983!> This is required to vectorize the red-black subroutine.
984!> Version for 3D-REAL arrays
985!------------------------------------------------------------------------------!
986    SUBROUTINE sort_k_to_even_odd_blocks( p_mg , glevel )
987
988
989       USE indices,                                                            &
990           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
991
992       IMPLICIT NONE
993
994       INTEGER(iwp), INTENT(IN) ::  glevel  !< grid level
995
996       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1,                               &
997                           nys_mg(glevel)-1:nyn_mg(glevel)+1,                  &
998                           nxl_mg(glevel)-1:nxr_mg(glevel)+1) ::               &
999                                      p_mg  !< array to be sorted
1000!
1001!--    Local variables
1002       INTEGER(iwp) :: i        !< index variable along x
1003       INTEGER(iwp) :: j        !< index variable along y
1004       INTEGER(iwp) :: k        !< index variable along z
1005       INTEGER(iwp) :: l        !< grid level
1006       INTEGER(iwp) :: ind      !< index variable along z
1007       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1) ::  tmp  !< odd-even sorted temporary array
1008
1009
1010       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'start' )
1011
1012       l = glevel
1013       ind_even_odd = even_odd_level(l)
1014
1015       !$OMP PARALLEL PRIVATE (i,j,k,ind,tmp)
1016       !$OMP DO
1017       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1018          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1019
1020!
1021!--          Sort the data with even k index
1022             ind = nzb-1
1023             DO  k = nzb, nzt_mg(l), 2
1024                ind = ind + 1
1025                tmp(ind) = p_mg(k,j,i)
1026             ENDDO
1027!
1028!--          Sort the data with odd k index
1029             DO  k = nzb+1, nzt_mg(l)+1, 2
1030                ind = ind + 1
1031                tmp(ind) = p_mg(k,j,i)
1032             ENDDO
1033
1034             p_mg(:,j,i) = tmp
1035
1036          ENDDO
1037       ENDDO
1038       !$OMP END PARALLEL
1039
1040       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'stop' )
1041
1042    END SUBROUTINE sort_k_to_even_odd_blocks
1043
1044
1045!------------------------------------------------------------------------------!
1046! Description:
1047! ------------
1048!> Sort k-Dimension from sequential into blocks of even and odd.
1049!> This is required to vectorize the red-black subroutine.
1050!> Version for 1D-REAL arrays
1051!------------------------------------------------------------------------------!
1052    SUBROUTINE sort_k_to_even_odd_blocks_1d( f_mg, f_mg_b, glevel )
1053
1054
1055       USE indices,                                                            &
1056           ONLY:  nzb, nzt_mg
1057
1058       IMPLICIT NONE
1059
1060       INTEGER(iwp), INTENT(IN) ::  glevel  !< grid level
1061
1062       REAL(wp), DIMENSION(nzb+1:nzt_mg(glevel)) ::  f_mg    !< 1D input array
1063       REAL(wp), DIMENSION(nzb:nzt_mg(glevel)+1) ::  f_mg_b  !< 1D output array
1064
1065!
1066!--    Local variables
1067       INTEGER(iwp) :: ind   !< index variable along z
1068       INTEGER(iwp) :: k     !< index variable along z
1069
1070
1071       ind = nzb - 1
1072!
1073!--    Sort the data with even k index
1074       DO  k = nzb, nzt_mg(glevel), 2
1075          ind = ind + 1
1076          IF ( k >= nzb+1  .AND.  k <= nzt_mg(glevel) )  THEN
1077             f_mg_b(ind) = f_mg(k)
1078          ENDIF
1079       ENDDO
1080!
1081!--    Sort the data with odd k index
1082       DO  k = nzb+1, nzt_mg(glevel)+1, 2
1083          ind = ind + 1
1084          IF( k >= nzb+1  .AND.  k <= nzt_mg(glevel) )  THEN
1085             f_mg_b(ind) = f_mg(k)
1086          ENDIF
1087       ENDDO
1088
1089    END SUBROUTINE sort_k_to_even_odd_blocks_1d
1090
1091
1092!------------------------------------------------------------------------------!
1093! Description:
1094! ------------
1095!> Sort k-Dimension from sequential into blocks of even and odd.
1096!> This is required to vectorize the red-black subroutine.
1097!> Version for 2D-INTEGER arrays
1098!------------------------------------------------------------------------------!
1099!    SUBROUTINE sort_k_to_even_odd_blocks_int( i_mg , glevel )
1100!
1101!
1102!       USE indices,                                                            &
1103!           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1104!
1105!       IMPLICIT NONE
1106!
1107!       INTEGER(iwp), INTENT(IN) ::  glevel  !< grid level
1108!
1109!       INTEGER(iwp), DIMENSION(nzb:nzt_mg(glevel)+1,                           &
1110!                               nys_mg(glevel)-1:nyn_mg(glevel)+1,              &
1111!                               nxl_mg(glevel)-1:nxr_mg(glevel)+1) ::           &
1112!                                    i_mg    !< array to be sorted
1113!
1114!--    Local variables
1115!       INTEGER(iwp) :: i        !< index variabel along x
1116!       INTEGER(iwp) :: j        !< index variable along y
1117!       INTEGER(iwp) :: k        !< index variable along z
1118!       INTEGER(iwp) :: l        !< grid level
1119!       INTEGER(iwp) :: ind      !< index variable along z
1120!       INTEGER(iwp),DIMENSION(nzb:nzt_mg(glevel)+1) ::  tmp  !< temporary odd-even sorted array
1121!
1122!
1123!       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'start' )
1124!
1125!       l = glevel
1126!       ind_even_odd = even_odd_level(l)
1127!
1128!       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1129!          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1130!
1131!
1132!--          Sort the data with even k index
1133!             ind = nzb-1
1134!             DO  k = nzb, nzt_mg(l), 2
1135!                ind = ind + 1
1136!                tmp(ind) = i_mg(k,j,i)
1137!             ENDDO
1138!
1139!--          Sort the data with odd k index
1140!             DO  k = nzb+1, nzt_mg(l)+1, 2
1141!                ind = ind + 1
1142!                tmp(ind) = i_mg(k,j,i)
1143!             ENDDO
1144!
1145!             i_mg(:,j,i) = tmp
1146!
1147!          ENDDO
1148!       ENDDO
1149!
1150!       CALL cpu_log( log_point_s(52), 'sort_k_to_even_odd', 'stop' )
1151!
1152!    END SUBROUTINE sort_k_to_even_odd_blocks_int
1153
1154
1155!------------------------------------------------------------------------------!
1156! Description:
1157! ------------
1158!> Sort k-dimension from blocks of even and odd into sequential
1159!------------------------------------------------------------------------------!
1160    SUBROUTINE sort_k_to_sequential( p_mg )
1161
1162
1163       USE control_parameters,                                                 &
1164           ONLY:  grid_level
1165
1166       USE indices,                                                            &
1167           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1168
1169       IMPLICIT NONE
1170
1171       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1172                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1173                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
1174                                     p_mg  !< array to be sorted
1175!
1176!--    Local variables
1177       INTEGER(iwp) :: i        !< index variable along x
1178       INTEGER(iwp) :: j        !< index variable along y
1179       INTEGER(iwp) :: k        !< index variable along z
1180       INTEGER(iwp) :: l        !< grid level
1181       INTEGER(iwp) :: ind      !< index variable along z
1182
1183       REAL(wp),DIMENSION(nzb:nzt_mg(grid_level)+1) ::  tmp
1184
1185
1186       l = grid_level
1187
1188       !$OMP PARALLEL PRIVATE (i,j,k,ind,tmp)
1189       !$OMP DO
1190       DO  i = nxl_mg(l)-1, nxr_mg(l)+1
1191          DO  j = nys_mg(l)-1, nyn_mg(l)+1
1192
1193             ind = nzb - 1
1194             tmp = p_mg(:,j,i)
1195             DO  k = nzb, nzt_mg(l), 2
1196                ind = ind + 1
1197                p_mg(k,j,i) = tmp(ind)
1198             ENDDO
1199
1200             DO  k = nzb+1, nzt_mg(l)+1, 2
1201                ind = ind + 1
1202                p_mg(k,j,i) = tmp(ind)
1203             ENDDO
1204          ENDDO
1205       ENDDO
1206       !$OMP END PARALLEL
1207
1208    END SUBROUTINE sort_k_to_sequential
1209
1210
1211!------------------------------------------------------------------------------!
1212! Description:
1213! ------------
1214!> Gather subdomain data from all PEs.
1215!------------------------------------------------------------------------------!
1216#if defined( __parallel )
1217    SUBROUTINE mg_gather( f2, f2_sub )
1218
1219       USE control_parameters,                                                 &
1220           ONLY:  grid_level
1221
1222       USE cpulog,                                                             &
1223           ONLY:  cpu_log, log_point_s
1224
1225       USE indices,                                                            &
1226           ONLY:  mg_loc_ind, nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1227
1228       IMPLICIT NONE
1229
1230       INTEGER(iwp) ::  i       !<
1231       INTEGER(iwp) ::  il      !<
1232       INTEGER(iwp) ::  ir      !<
1233       INTEGER(iwp) ::  j       !<
1234       INTEGER(iwp) ::  jn      !<
1235       INTEGER(iwp) ::  js      !<
1236       INTEGER(iwp) ::  k       !<
1237       INTEGER(iwp) ::  nwords  !<
1238
1239       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1240                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
1241                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f2    !<
1242       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1243                       nys_mg(grid_level)-1:nyn_mg(grid_level)+1,              &
1244                       nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::  f2_l  !<
1245
1246       REAL(wp), DIMENSION(nzb:mg_loc_ind(5,myid)+1,                           &
1247                           mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,          &
1248                           mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) ::  f2_sub  !<
1249
1250
1251       CALL cpu_log( log_point_s(34), 'mg_gather', 'start' )
1252
1253       f2_l = 0.0_wp
1254
1255!
1256!--    Store the local subdomain array on the total array
1257       js = mg_loc_ind(3,myid)
1258       IF ( south_border_pe )  js = js - 1
1259       jn = mg_loc_ind(4,myid)
1260       IF ( north_border_pe )  jn = jn + 1
1261       il = mg_loc_ind(1,myid)
1262       IF ( left_border_pe )   il = il - 1
1263       ir = mg_loc_ind(2,myid)
1264       IF ( right_border_pe )  ir = ir + 1
1265       DO  i = il, ir
1266          DO  j = js, jn
1267             DO  k = nzb, nzt_mg(grid_level)+1
1268                f2_l(k,j,i) = f2_sub(k,j,i)
1269             ENDDO
1270          ENDDO
1271       ENDDO
1272
1273!
1274!--    Find out the number of array elements of the total array
1275       nwords = SIZE( f2 )
1276
1277!
1278!--    Gather subdomain data from all PEs
1279       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1280       CALL MPI_ALLREDUCE( f2_l(nzb,nys_mg(grid_level)-1,nxl_mg(grid_level)-1), &
1281                           f2(nzb,nys_mg(grid_level)-1,nxl_mg(grid_level)-1),   &
1282                           nwords, MPI_REAL, MPI_SUM, comm2d, ierr )
1283
1284       CALL cpu_log( log_point_s(34), 'mg_gather', 'stop' )
1285   
1286    END SUBROUTINE mg_gather
1287#endif
1288
1289
1290!------------------------------------------------------------------------------!
1291! Description:
1292! ------------
1293!> @todo It might be possible to improve the speed of this routine by using
1294!>       non-blocking communication
1295!------------------------------------------------------------------------------!
1296#if defined( __parallel )
1297    SUBROUTINE mg_scatter( p2, p2_sub )
1298
1299       USE control_parameters,                                                 &
1300           ONLY:  grid_level
1301
1302       USE cpulog,                                                             &
1303           ONLY:  cpu_log, log_point_s
1304
1305       USE indices,                                                            &
1306           ONLY:  mg_loc_ind, nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1307
1308       IMPLICIT NONE
1309
1310       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1311                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1312                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  p2  !<
1313
1314       REAL(wp), DIMENSION(nzb:mg_loc_ind(5,myid)+1,                           &
1315                           mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,          &
1316                           mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) ::  p2_sub  !<
1317
1318
1319       CALL cpu_log( log_point_s(35), 'mg_scatter', 'start' )
1320
1321       p2_sub = p2(:,mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1, &
1322                     mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1)
1323
1324       CALL cpu_log( log_point_s(35), 'mg_scatter', 'stop' )
1325   
1326    END SUBROUTINE mg_scatter
1327#endif
1328
1329!------------------------------------------------------------------------------!
1330! Description:
1331! ------------
1332!> This is where the multigrid technique takes place. V- and W- Cycle are
1333!> implemented and steered by the parameter "gamma". Parameter "nue" determines
1334!> the convergence of the multigrid iterative solution. There are nue times
1335!> RB-GS iterations. It should be set to "1" or "2", considering the time effort
1336!> one would like to invest. Last choice shows a very good converging factor,
1337!> but leads to an increase in computing time.
1338!------------------------------------------------------------------------------!
1339    RECURSIVE SUBROUTINE next_mg_level( f_mg, p_mg, p3, r )
1340
1341       USE control_parameters,                                                 &
1342           ONLY:  bc_lr_dirrad, bc_lr_raddir, bc_ns_dirrad, bc_ns_raddir,      &
1343                  child_domain, gamma_mg, grid_level_count, maximum_grid_level,&
1344                  mg_switch_to_pe0_level, mg_switch_to_pe0, ngsrb
1345
1346       USE indices,                                                            &
1347           ONLY:  mg_loc_ind, nxl, nxl_mg, nxr, nxr_mg, nys, nys_mg, nyn,      &
1348                  nyn_mg, nzb, nzt, nzt_mg
1349
1350       IMPLICIT NONE
1351
1352       INTEGER(iwp) ::  i            !< index variable along x
1353       INTEGER(iwp) ::  j            !< index variable along y
1354       INTEGER(iwp) ::  k            !< index variable along z
1355       INTEGER(iwp) ::  nxl_mg_save  !<
1356       INTEGER(iwp) ::  nxr_mg_save  !<
1357       INTEGER(iwp) ::  nyn_mg_save  !<
1358       INTEGER(iwp) ::  nys_mg_save  !<
1359       INTEGER(iwp) ::  nzt_mg_save  !<
1360
1361       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1362                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1363                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: f_mg  !<
1364       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1365                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1366                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: p_mg  !<
1367       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1368                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1369                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: p3    !<
1370       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1371                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1372                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) :: r     !<
1373
1374       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1375                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1376                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  f2  !<
1377       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level-1)+1,                         &
1378                           nys_mg(grid_level-1)-1:nyn_mg(grid_level-1)+1,      &
1379                           nxl_mg(grid_level-1)-1:nxr_mg(grid_level-1)+1) ::  p2  !<
1380
1381       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  f2_sub  !<
1382
1383#if defined( __parallel )
1384       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  p2_sub  !<
1385#endif
1386
1387!
1388!--    Restriction to the coarsest grid
1389    10 IF ( grid_level == 1 )  THEN
1390
1391!
1392!--       Solution on the coarsest grid. Double the number of Gauss-Seidel
1393!--       iterations in order to get a more accurate solution.
1394          ngsrb = 2 * ngsrb
1395
1396          ind_even_odd = even_odd_level(grid_level)
1397
1398          CALL redblack( f_mg, p_mg )
1399
1400          ngsrb = ngsrb / 2
1401
1402
1403       ELSEIF ( grid_level /= 1 )  THEN
1404
1405          grid_level_count(grid_level) = grid_level_count(grid_level) + 1
1406
1407!
1408!--       Solution on the actual grid level
1409          ind_even_odd = even_odd_level(grid_level)
1410
1411          CALL redblack( f_mg, p_mg )
1412
1413!
1414!--       Determination of the actual residual
1415          CALL resid( f_mg, p_mg, r )
1416
1417!--       Restriction of the residual (finer grid values!) to the next coarser
1418!--       grid. Therefore, the grid level has to be decremented now. nxl..nzt have
1419!--       to be set to the coarse grid values, because these variables are needed
1420!--       for the exchange of ghost points in routine exchange_horiz
1421          grid_level = grid_level - 1
1422
1423          nxl = nxl_mg(grid_level)
1424          nys = nys_mg(grid_level)
1425          nxr = nxr_mg(grid_level)
1426          nyn = nyn_mg(grid_level)
1427          nzt = nzt_mg(grid_level)
1428
1429          IF ( grid_level == mg_switch_to_pe0_level )  THEN
1430
1431!
1432!--          From this level on, calculations are done on PE0 only.
1433!--          First, carry out restriction on the subdomain.
1434!--          Therefore, indices of the level have to be changed to subdomain
1435!--          values in between (otherwise, the restrict routine would expect
1436!--          the gathered array)
1437
1438             nxl_mg_save = nxl_mg(grid_level)
1439             nxr_mg_save = nxr_mg(grid_level)
1440             nys_mg_save = nys_mg(grid_level)
1441             nyn_mg_save = nyn_mg(grid_level)
1442             nzt_mg_save = nzt_mg(grid_level)
1443             nxl_mg(grid_level) = mg_loc_ind(1,myid)
1444             nxr_mg(grid_level) = mg_loc_ind(2,myid)
1445             nys_mg(grid_level) = mg_loc_ind(3,myid)
1446             nyn_mg(grid_level) = mg_loc_ind(4,myid)
1447             nzt_mg(grid_level) = mg_loc_ind(5,myid)
1448             nxl = mg_loc_ind(1,myid)
1449             nxr = mg_loc_ind(2,myid)
1450             nys = mg_loc_ind(3,myid)
1451             nyn = mg_loc_ind(4,myid)
1452             nzt = mg_loc_ind(5,myid)
1453
1454             ALLOCATE( f2_sub(nzb:nzt_mg(grid_level)+1,                    &
1455                              nys_mg(grid_level)-1:nyn_mg(grid_level)+1,   &
1456                              nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) )
1457
1458             CALL restrict( f2_sub, r )
1459
1460!
1461!--          Restore the correct indices of this level
1462             nxl_mg(grid_level) = nxl_mg_save
1463             nxr_mg(grid_level) = nxr_mg_save
1464             nys_mg(grid_level) = nys_mg_save
1465             nyn_mg(grid_level) = nyn_mg_save
1466             nzt_mg(grid_level) = nzt_mg_save
1467             nxl = nxl_mg(grid_level)
1468             nxr = nxr_mg(grid_level)
1469             nys = nys_mg(grid_level)
1470             nyn = nyn_mg(grid_level)
1471             nzt = nzt_mg(grid_level)
1472!
1473!--          Gather all arrays from the subdomains on PE0
1474#if defined( __parallel )
1475             CALL mg_gather( f2, f2_sub )
1476#endif
1477
1478!
1479!--          Set switch for routine exchange_horiz, that no ghostpoint exchange
1480!--          has to be carried out from now on
1481             mg_switch_to_pe0 = .TRUE.
1482
1483!
1484!--          In case of non-cyclic lateral boundary conditions, both in- and
1485!--          outflow conditions have to be used on all PEs after the switch,
1486!--          because then they have the total domain.
1487             IF ( bc_lr_dirrad )  THEN
1488                bc_dirichlet_l  = .TRUE.
1489                bc_dirichlet_r  = .FALSE.
1490                bc_radiation_l = .FALSE.
1491                bc_radiation_r = .TRUE.
1492             ELSEIF ( bc_lr_raddir )  THEN
1493                bc_dirichlet_l  = .FALSE.
1494                bc_dirichlet_r  = .TRUE.
1495                bc_radiation_l = .TRUE.
1496                bc_radiation_r = .FALSE.
1497             ELSEIF ( child_domain  .OR.  nesting_offline )  THEN
1498                bc_dirichlet_l = .TRUE.
1499                bc_dirichlet_r = .TRUE.
1500             ENDIF
1501
1502             IF ( bc_ns_dirrad )  THEN
1503                bc_dirichlet_n  = .TRUE.
1504                bc_dirichlet_s  = .FALSE.
1505                bc_radiation_n = .FALSE.
1506                bc_radiation_s = .TRUE.
1507             ELSEIF ( bc_ns_raddir )  THEN
1508                bc_dirichlet_n  = .FALSE.
1509                bc_dirichlet_s  = .TRUE.
1510                bc_radiation_n = .TRUE.
1511                bc_radiation_s = .FALSE.
1512             ELSEIF ( child_domain  .OR.  nesting_offline)  THEN
1513                bc_dirichlet_s = .TRUE.
1514                bc_dirichlet_n = .TRUE.
1515             ENDIF
1516
1517             DEALLOCATE( f2_sub )
1518
1519          ELSE
1520
1521             CALL restrict( f2, r )
1522
1523             ind_even_odd = even_odd_level(grid_level)  ! must be after restrict
1524
1525          ENDIF
1526
1527          p2 = 0.0_wp
1528
1529!
1530!--       Repeat the same procedure till the coarsest grid is reached
1531          CALL next_mg_level( f2, p2, p3, r )
1532
1533       ENDIF
1534
1535!
1536!--    Now follows the prolongation
1537       IF ( grid_level >= 2 )  THEN
1538
1539!
1540!--       Prolongation of the new residual. The values are transferred
1541!--       from the coarse to the next finer grid.
1542          IF ( grid_level == mg_switch_to_pe0_level+1 )  THEN
1543
1544#if defined( __parallel )
1545!
1546!--          At this level, the new residual first has to be scattered from
1547!--          PE0 to the other PEs
1548             ALLOCATE( p2_sub(nzb:mg_loc_ind(5,myid)+1,             &
1549                       mg_loc_ind(3,myid)-1:mg_loc_ind(4,myid)+1,   &
1550                       mg_loc_ind(1,myid)-1:mg_loc_ind(2,myid)+1) )
1551
1552             CALL mg_scatter( p2, p2_sub )
1553
1554!
1555!--          Therefore, indices of the previous level have to be changed to
1556!--          subdomain values in between (otherwise, the prolong routine would
1557!--          expect the gathered array)
1558             nxl_mg_save = nxl_mg(grid_level-1)
1559             nxr_mg_save = nxr_mg(grid_level-1)
1560             nys_mg_save = nys_mg(grid_level-1)
1561             nyn_mg_save = nyn_mg(grid_level-1)
1562             nzt_mg_save = nzt_mg(grid_level-1)
1563             nxl_mg(grid_level-1) = mg_loc_ind(1,myid)
1564             nxr_mg(grid_level-1) = mg_loc_ind(2,myid)
1565             nys_mg(grid_level-1) = mg_loc_ind(3,myid)
1566             nyn_mg(grid_level-1) = mg_loc_ind(4,myid)
1567             nzt_mg(grid_level-1) = mg_loc_ind(5,myid)
1568
1569!
1570!--          Set switch for routine exchange_horiz, that ghostpoint exchange
1571!--          has to be carried again out from now on
1572             mg_switch_to_pe0 = .FALSE.
1573
1574!
1575!--          For non-cyclic lateral boundary conditions and in case of nesting,
1576!--          restore the in-/outflow conditions.
1577             bc_dirichlet_l = .FALSE.;  bc_dirichlet_r = .FALSE.
1578             bc_dirichlet_n = .FALSE.;  bc_dirichlet_s = .FALSE.
1579             bc_radiation_l = .FALSE.;  bc_radiation_r = .FALSE.
1580             bc_radiation_n = .FALSE.;  bc_radiation_s = .FALSE.
1581
1582             IF ( pleft == MPI_PROC_NULL )  THEN
1583                IF ( bc_lr_dirrad  .OR.  child_domain  .OR.  nesting_offline ) &
1584                THEN
1585                   bc_dirichlet_l = .TRUE.
1586                ELSEIF ( bc_lr_raddir )  THEN
1587                   bc_radiation_l = .TRUE.
1588                ENDIF
1589             ENDIF
1590
1591             IF ( pright == MPI_PROC_NULL )  THEN
1592                IF ( bc_lr_dirrad )  THEN
1593                   bc_radiation_r = .TRUE.
1594                ELSEIF ( bc_lr_raddir  .OR.  child_domain  .OR.                &
1595                         nesting_offline )  THEN
1596                   bc_dirichlet_r = .TRUE.
1597                ENDIF
1598             ENDIF
1599
1600             IF ( psouth == MPI_PROC_NULL )  THEN
1601                IF ( bc_ns_dirrad )  THEN
1602                   bc_radiation_s = .TRUE.
1603                ELSEIF ( bc_ns_raddir  .OR.  child_domain  .OR.                &
1604                         nesting_offline )  THEN
1605                   bc_dirichlet_s = .TRUE.
1606                ENDIF
1607             ENDIF
1608
1609             IF ( pnorth == MPI_PROC_NULL )  THEN
1610                IF ( bc_ns_dirrad  .OR.  child_domain  .OR.  nesting_offline ) &
1611                THEN
1612                   bc_dirichlet_n = .TRUE.
1613                ELSEIF ( bc_ns_raddir )  THEN
1614                   bc_radiation_n = .TRUE.
1615                ENDIF
1616             ENDIF
1617
1618             CALL prolong( p2_sub, p3 )
1619
1620!
1621!--          Restore the correct indices of the previous level
1622             nxl_mg(grid_level-1) = nxl_mg_save
1623             nxr_mg(grid_level-1) = nxr_mg_save
1624             nys_mg(grid_level-1) = nys_mg_save
1625             nyn_mg(grid_level-1) = nyn_mg_save
1626             nzt_mg(grid_level-1) = nzt_mg_save
1627
1628             DEALLOCATE( p2_sub )
1629#endif
1630
1631          ELSE
1632
1633             CALL prolong( p2, p3 )
1634
1635          ENDIF
1636
1637!
1638!--       Computation of the new pressure correction. Therefore,
1639!--       values from prior grids are added up automatically stage by stage.
1640          DO  i = nxl_mg(grid_level)-1, nxr_mg(grid_level)+1
1641             DO  j = nys_mg(grid_level)-1, nyn_mg(grid_level)+1
1642                DO  k = nzb, nzt_mg(grid_level)+1
1643                   p_mg(k,j,i) = p_mg(k,j,i) + p3(k,j,i)
1644                ENDDO
1645             ENDDO
1646          ENDDO
1647
1648!
1649!--       Relaxation of the new solution
1650          CALL redblack( f_mg, p_mg )
1651
1652       ENDIF
1653
1654
1655!
1656!--    The following few lines serve the steering of the multigrid scheme
1657       IF ( grid_level == maximum_grid_level )  THEN
1658
1659          GOTO 20
1660
1661       ELSEIF ( grid_level /= maximum_grid_level  .AND.  grid_level /= 1  .AND. &
1662                grid_level_count(grid_level) /= gamma_mg )  THEN
1663
1664          GOTO 10
1665
1666       ENDIF
1667
1668!
1669!--    Reset counter for the next call of poismg
1670       grid_level_count(grid_level) = 0
1671
1672!
1673!--    Continue with the next finer level. nxl..nzt have to be
1674!--    set to the finer grid values, because these variables are needed for the
1675!--    exchange of ghost points in routine exchange_horiz
1676       grid_level = grid_level + 1
1677       ind_even_odd = even_odd_level(grid_level)
1678
1679       nxl = nxl_mg(grid_level)
1680       nxr = nxr_mg(grid_level)
1681       nys = nys_mg(grid_level)
1682       nyn = nyn_mg(grid_level)
1683       nzt = nzt_mg(grid_level)
1684
1685    20 CONTINUE
1686
1687    END SUBROUTINE next_mg_level
1688
1689
1690!------------------------------------------------------------------------------!
1691! Description:
1692! ------------
1693!> Initial settings for sorting k-dimension from sequential order (alternate
1694!> even/odd) into blocks of even and odd or vice versa
1695!------------------------------------------------------------------------------!
1696    SUBROUTINE init_even_odd_blocks
1697
1698
1699       USE arrays_3d,                                                          &
1700           ONLY:  f1_mg, f2_mg, f3_mg
1701
1702       USE control_parameters,                                                 &
1703           ONLY:  grid_level, maximum_grid_level
1704
1705       USE indices,                                                            &
1706           ONLY:  nzb, nzt, nzt_mg
1707
1708       USE indices,                                                            &
1709           ONLY:  nzb, nzt_mg
1710
1711       IMPLICIT NONE
1712!
1713!--    Local variables
1714       INTEGER(iwp) ::  i     !< 
1715       INTEGER(iwp) ::  l     !<
1716
1717       LOGICAL, SAVE ::  lfirst = .TRUE.
1718
1719
1720       IF ( .NOT. lfirst )  RETURN
1721
1722       ALLOCATE( even_odd_level(maximum_grid_level) )
1723
1724       ALLOCATE( f1_mg_b(nzb:nzt+1,maximum_grid_level),                        &
1725                 f2_mg_b(nzb:nzt+1,maximum_grid_level),                        &
1726                 f3_mg_b(nzb:nzt+1,maximum_grid_level) )
1727
1728!
1729!--    Set border index between the even and odd block
1730       DO  i = maximum_grid_level, 1, -1
1731          even_odd_level(i) = nzt_mg(i) / 2
1732       ENDDO
1733
1734!
1735!--    Sort grid coefficients used in red/black scheme and for calculating the
1736!--    residual to block (even/odd) structure
1737       DO  l = maximum_grid_level, 1 , -1
1738          CALL sort_k_to_even_odd_blocks( f1_mg(nzb+1:nzt_mg(grid_level),l),   &
1739                                          f1_mg_b(nzb:nzt_mg(grid_level)+1,l), &
1740                                          l )
1741          CALL sort_k_to_even_odd_blocks( f2_mg(nzb+1:nzt_mg(grid_level),l),   &
1742                                          f2_mg_b(nzb:nzt_mg(grid_level)+1,l), &
1743                                          l )
1744          CALL sort_k_to_even_odd_blocks( f3_mg(nzb+1:nzt_mg(grid_level),l),   &
1745                                          f3_mg_b(nzb:nzt_mg(grid_level)+1,l), &
1746                                          l )
1747       ENDDO
1748
1749       lfirst = .FALSE.
1750
1751     END SUBROUTINE init_even_odd_blocks
1752
1753
1754!------------------------------------------------------------------------------!
1755! Description:
1756! ------------
1757!> Special exchange_horiz subroutine for use in redblack. Transfers only
1758!> "red" or "black" data points.
1759!------------------------------------------------------------------------------!
1760     SUBROUTINE special_exchange_horiz( p_mg, color )
1761
1762
1763       USE control_parameters,                                                 &
1764           ONLY:  grid_level
1765
1766#if defined( __parallel )
1767       USE control_parameters,                                                 &
1768           ONLY:  mg_switch_to_pe0_level, synchronous_exchange
1769#endif
1770
1771       USE indices,                                                            &
1772           ONLY:  nxl_mg, nxr_mg, nys_mg, nyn_mg, nzb, nzt_mg
1773
1774#if defined( __parallel )
1775       USE indices,                                                            &
1776           ONLY:  nxl, nxr, nys, nyn, nzt
1777#endif
1778
1779       IMPLICIT NONE
1780
1781       REAL(wp), DIMENSION(nzb:nzt_mg(grid_level)+1,                           &
1782                           nys_mg(grid_level)-1:nyn_mg(grid_level)+1,          &
1783                           nxl_mg(grid_level)-1:nxr_mg(grid_level)+1) ::       &
1784                                    p_mg   !< treated array
1785
1786       INTEGER(iwp), INTENT(IN) ::  color  !< flag for grid point type (red or black)
1787
1788#if defined ( __parallel )
1789!
1790!--    Local variables
1791       INTEGER(iwp) ::  i        !< index variable along x
1792       INTEGER(iwp) ::  i1       !< index variable along x on coarse level
1793       INTEGER(iwp) ::  i2       !< index variable along x on coarse level
1794
1795       INTEGER(iwp) ::  j        !< index variable along y
1796       INTEGER(iwp) ::  j1       !< index variable along y on coarse level
1797       INTEGER(iwp) ::  j2       !< index variable along y on coarse level
1798       INTEGER(iwp) ::  k        !< index variable along z
1799       INTEGER(iwp) ::  l        !< short for grid level
1800       INTEGER(iwp) ::  jys      !< index for lower local PE boundary along y
1801       INTEGER(iwp) ::  jyn      !< index for upper local PE boundary along y
1802       INTEGER(iwp) ::  ixl      !< index for lower local PE boundary along x
1803       INTEGER(iwp) ::  ixr      !< index for upper local PE boundary along x
1804
1805       LOGICAL      ::  synchronous_exchange_save    !< dummy to reset synchronous_exchange to prescribed value
1806
1807       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  temp  !< temporary array on next coarser grid level
1808
1809       synchronous_exchange_save   = synchronous_exchange
1810       synchronous_exchange        = .FALSE.
1811
1812       l = grid_level
1813
1814       ind_even_odd = even_odd_level(grid_level)
1815
1816!
1817!--    Restricted transfer only on finer levels with enough data.
1818!--    Restricted transfer is not possible for levels smaller or equal to
1819!--    'switch to PE0 levels', since array bounds does not fit. Moreover,
1820!--    it is not possible for the coarsest grid level, since the dimensions
1821!--    of temp are not defined. For such cases, normal exchange_horiz is called.
1822       IF ( l > 1 .AND. l > mg_switch_to_pe0_level + 1 .AND.                   &
1823            ( ngp_xz(grid_level) >= 900 .OR. ngp_yz(grid_level) >= 900 ) )  THEN
1824
1825          jys = nys_mg(grid_level-1)
1826          jyn = nyn_mg(grid_level-1)
1827          ixl = nxl_mg(grid_level-1)
1828          ixr = nxr_mg(grid_level-1)
1829          ALLOCATE( temp(nzb:nzt_mg(l-1)+1,jys-1:jyn+1,ixl-1:ixr+1) )
1830!
1831!--       Handling the even k Values
1832!--       Collecting data for the north - south exchange
1833!--       Since only every second value has to be transfered, data are stored
1834!--       on the next coarser grid level, because the arrays on that level
1835!--       have just the required size
1836          i1 = nxl_mg(grid_level-1)
1837          i2 = nxl_mg(grid_level-1)
1838
1839          DO  i = nxl_mg(l), nxr_mg(l), 2
1840             DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
1841
1842                IF ( j == nys_mg(l) )  THEN
1843                   !DIR$ IVDEP
1844                   DO  k = ind_even_odd+1, nzt_mg(l)
1845                      temp(k-ind_even_odd,jys,i1) = p_mg(k,j,i)
1846                   ENDDO
1847                   i1 = i1 + 1
1848
1849                ENDIF
1850
1851                IF ( j == nyn_mg(l) )  THEN
1852                   !DIR$ IVDEP
1853                   DO  k = ind_even_odd+1, nzt_mg(l)
1854                      temp(k-ind_even_odd,jyn,i2) = p_mg(k,j,i)
1855                   ENDDO
1856                   i2 = i2 + 1
1857
1858                ENDIF
1859
1860             ENDDO
1861          ENDDO
1862
1863          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
1864             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
1865
1866                IF ( j == nys_mg(l) )  THEN
1867                   !DIR$ IVDEP
1868                   DO  k = ind_even_odd+1, nzt_mg(l)
1869                      temp(k-ind_even_odd,jys,i1) = p_mg(k,j,i)
1870                   ENDDO
1871                   i1 = i1 + 1
1872
1873                ENDIF
1874
1875                IF ( j == nyn_mg(l) )  THEN
1876                   !DIR$ IVDEP
1877                   DO  k = ind_even_odd+1, nzt_mg(l)
1878                      temp(k-ind_even_odd,jyn,i2) = p_mg(k,j,i)
1879                   ENDDO
1880                   i2 = i2 + 1
1881
1882                ENDIF
1883
1884             ENDDO
1885          ENDDO
1886
1887          grid_level = grid_level-1
1888
1889          nxl = nxl_mg(grid_level)
1890          nys = nys_mg(grid_level)
1891          nxr = nxr_mg(grid_level)
1892          nyn = nyn_mg(grid_level)
1893          nzt = nzt_mg(grid_level)
1894
1895          send_receive = 'ns'
1896          CALL exchange_horiz( temp, 1 )
1897
1898          grid_level = grid_level+1
1899
1900          i1 = nxl_mg(grid_level-1)
1901          i2 = nxl_mg(grid_level-1)
1902
1903          DO  i = nxl_mg(l), nxr_mg(l), 2
1904             DO  j = nys_mg(l) + 2 - color, nyn_mg(l), 2
1905
1906                IF ( j == nys_mg(l) )  THEN
1907                   !DIR$ IVDEP
1908                   DO  k = ind_even_odd+1, nzt_mg(l)
1909                      p_mg(k,nyn_mg(l)+1,i) = temp(k-ind_even_odd,jyn+1,i1)
1910                   ENDDO
1911                   i1 = i1 + 1
1912
1913                ENDIF
1914
1915                IF ( j == nyn_mg(l) )  THEN
1916                   !DIR$ IVDEP
1917                   DO  k = ind_even_odd+1, nzt_mg(l)
1918                      p_mg(k,nys_mg(l)-1,i) = temp(k-ind_even_odd,jys-1,i2)
1919                   ENDDO
1920                   i2 = i2 + 1
1921
1922                ENDIF
1923
1924             ENDDO
1925          ENDDO
1926
1927          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
1928             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
1929
1930                IF ( j == nys_mg(l) )  THEN
1931                   !DIR$ IVDEP
1932                   DO  k = ind_even_odd+1, nzt_mg(l)
1933                      p_mg(k,nyn_mg(l)+1,i) = temp(k-ind_even_odd,jyn+1,i1)
1934                   ENDDO
1935                   i1 = i1 + 1
1936
1937                ENDIF
1938
1939                IF ( j == nyn_mg(l) )  THEN
1940                   !DIR$ IVDEP
1941                   DO  k = ind_even_odd+1, nzt_mg(l)
1942                      p_mg(k,nys_mg(l)-1,i) = temp(k-ind_even_odd,jys-1,i2)
1943                   ENDDO
1944                   i2 = i2 + 1
1945
1946                ENDIF
1947
1948             ENDDO
1949          ENDDO
1950
1951!
1952!--       Collecting data for the left - right exchange
1953!--       Since only every second value has to be transfered, data are stored
1954!--       on the next coarser grid level, because the arrays on that level
1955!--       have just the required size
1956          j1 = nys_mg(grid_level-1)
1957          j2 = nys_mg(grid_level-1)
1958
1959          DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
1960             DO  i = nxl_mg(l), nxr_mg(l), 2
1961
1962                IF ( i == nxl_mg(l) )  THEN
1963                   !DIR$ IVDEP
1964                   DO  k = ind_even_odd+1, nzt_mg(l)
1965                      temp(k-ind_even_odd,j1,ixl) = p_mg(k,j,i)
1966                   ENDDO
1967                   j1 = j1 + 1
1968
1969                ENDIF
1970
1971                IF ( i == nxr_mg(l) )  THEN
1972                   !DIR$ IVDEP
1973                   DO  k = ind_even_odd+1, nzt_mg(l)
1974                      temp(k-ind_even_odd,j2,ixr) = p_mg(k,j,i)
1975                   ENDDO
1976                   j2 = j2 + 1
1977
1978                ENDIF
1979
1980             ENDDO
1981          ENDDO
1982
1983          DO j = nys_mg(l) + (color-1), nyn_mg(l), 2
1984             DO  i = nxl_mg(l)+1, nxr_mg(l), 2
1985
1986                IF ( i == nxl_mg(l) )  THEN
1987                   !DIR$ IVDEP
1988                   DO  k = ind_even_odd+1, nzt_mg(l)
1989                      temp(k-ind_even_odd,j1,ixl) = p_mg(k,j,i)
1990                   ENDDO
1991                   j1 = j1 + 1
1992
1993                ENDIF
1994
1995                IF ( i == nxr_mg(l) )  THEN
1996                   !DIR$ IVDEP
1997                   DO  k = ind_even_odd+1, nzt_mg(l)
1998                      temp(k-ind_even_odd,j2,ixr) = p_mg(k,j,i)
1999                   ENDDO
2000                   j2 = j2 + 1
2001
2002                ENDIF
2003
2004             ENDDO
2005          ENDDO
2006
2007          grid_level = grid_level-1
2008          send_receive = 'lr'
2009
2010          CALL exchange_horiz( temp, 1 )
2011
2012          grid_level = grid_level+1
2013
2014          j1 = nys_mg(grid_level-1)
2015          j2 = nys_mg(grid_level-1)
2016
2017          DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2018             DO  i = nxl_mg(l), nxr_mg(l), 2
2019
2020                IF ( i == nxl_mg(l) )  THEN
2021                   !DIR$ IVDEP
2022                   DO  k = ind_even_odd+1, nzt_mg(l)
2023                      p_mg(k,j,nxr_mg(l)+1)  = temp(k-ind_even_odd,j1,ixr+1)
2024                   ENDDO
2025                   j1 = j1 + 1
2026
2027                ENDIF
2028
2029                IF ( i == nxr_mg(l) )  THEN
2030                   !DIR$ IVDEP
2031                   DO  k = ind_even_odd+1, nzt_mg(l)
2032                      p_mg(k,j,nxl_mg(l)-1)  = temp(k-ind_even_odd,j2,ixl-1)
2033                   ENDDO
2034                   j2 = j2 + 1
2035
2036                ENDIF
2037
2038             ENDDO
2039          ENDDO
2040
2041          DO j = nys_mg(l) + (color-1), nyn_mg(l), 2
2042             DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2043
2044                IF ( i == nxl_mg(l) )  THEN
2045                   !DIR$ IVDEP
2046                   DO  k = ind_even_odd+1, nzt_mg(l)
2047                      p_mg(k,j,nxr_mg(l)+1)  = temp(k-ind_even_odd,j1,ixr+1)
2048                   ENDDO
2049                   j1 = j1 + 1
2050
2051                ENDIF
2052
2053                IF ( i == nxr_mg(l) )  THEN
2054                   !DIR$ IVDEP
2055                   DO  k = ind_even_odd+1, nzt_mg(l)
2056                      p_mg(k,j,nxl_mg(l)-1)  = temp(k-ind_even_odd,j2,ixl-1)
2057                   ENDDO
2058                   j2 = j2 + 1
2059
2060                ENDIF
2061
2062             ENDDO
2063          ENDDO
2064
2065!
2066!--       Now handling the even k values
2067!--       Collecting data for the north - south exchange
2068!--       Since only every second value has to be transfered, data are stored
2069!--       on the next coarser grid level, because the arrays on that level
2070!--       have just the required size
2071          i1 = nxl_mg(grid_level-1)
2072          i2 = nxl_mg(grid_level-1)
2073
2074          DO  i = nxl_mg(l), nxr_mg(l), 2
2075             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2076
2077                IF ( j == nys_mg(l) )  THEN
2078                   !DIR$ IVDEP
2079                   DO  k = nzb+1, ind_even_odd
2080                      temp(k,jys,i1) = p_mg(k,j,i)
2081                   ENDDO
2082                   i1 = i1 + 1
2083
2084                ENDIF
2085
2086                IF ( j == nyn_mg(l) )  THEN
2087                   !DIR$ IVDEP
2088                   DO  k = nzb+1, ind_even_odd
2089                      temp(k,jyn,i2) = p_mg(k,j,i)
2090                   ENDDO
2091                   i2 = i2 + 1
2092
2093                ENDIF
2094
2095             ENDDO
2096          ENDDO
2097
2098          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2099             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2100
2101                IF ( j == nys_mg(l) )  THEN
2102                   !DIR$ IVDEP
2103                   DO  k = nzb+1, ind_even_odd
2104                      temp(k,jys,i1) = p_mg(k,j,i)
2105                   ENDDO
2106                   i1 = i1 + 1
2107
2108                ENDIF
2109
2110                IF ( j == nyn_mg(l) )  THEN
2111                   !DIR$ IVDEP
2112                   DO  k = nzb+1, ind_even_odd
2113                      temp(k,jyn,i2) = p_mg(k,j,i)
2114                   ENDDO
2115                   i2 = i2 + 1
2116
2117                ENDIF
2118
2119             ENDDO
2120          ENDDO
2121
2122          grid_level = grid_level-1
2123
2124          send_receive = 'ns'
2125          CALL exchange_horiz( temp, 1 )
2126
2127          grid_level = grid_level+1
2128
2129          i1 = nxl_mg(grid_level-1)
2130          i2 = nxl_mg(grid_level-1)
2131
2132          DO  i = nxl_mg(l), nxr_mg(l), 2
2133             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2134
2135                IF ( j == nys_mg(l) )  THEN
2136                   !DIR$ IVDEP
2137                   DO  k = nzb+1, ind_even_odd
2138                      p_mg(k,nyn_mg(l)+1,i) = temp(k,jyn+1,i1)
2139                   ENDDO
2140                   i1 = i1 + 1
2141
2142                ENDIF
2143
2144                IF ( j == nyn_mg(l) )  THEN
2145                   !DIR$ IVDEP
2146                   DO  k = nzb+1, ind_even_odd
2147                      p_mg(k,nys_mg(l)-1,i) = temp(k,jys-1,i2)
2148                   ENDDO
2149                   i2 = i2 + 1
2150
2151                ENDIF
2152
2153             ENDDO
2154          ENDDO
2155
2156          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2157             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2158
2159                IF ( j == nys_mg(l) )  THEN
2160                   !DIR$ IVDEP
2161                   DO  k = nzb+1, ind_even_odd
2162                      p_mg(k,nyn_mg(l)+1,i) = temp(k,jyn+1,i1)
2163                   ENDDO
2164                   i1 = i1 + 1
2165
2166                ENDIF
2167
2168                IF ( j == nyn_mg(l) )  THEN
2169                   !DIR$ IVDEP
2170                   DO  k = nzb+1, ind_even_odd
2171                      p_mg(k,nys_mg(l)-1,i) = temp(k,jys-1,i2)
2172                   ENDDO
2173                   i2 = i2 + 1
2174
2175                ENDIF
2176
2177             ENDDO
2178          ENDDO
2179
2180          j1 = nys_mg(grid_level-1)
2181          j2 = nys_mg(grid_level-1)
2182
2183          DO  i = nxl_mg(l), nxr_mg(l), 2
2184             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2185
2186                IF ( i == nxl_mg(l) )  THEN
2187                   !DIR$ IVDEP
2188                   DO  k = nzb+1, ind_even_odd
2189                      temp(k,j1,ixl) = p_mg(k,j,i)
2190                   ENDDO
2191                   j1 = j1 + 1
2192
2193                ENDIF
2194
2195                IF ( i == nxr_mg(l) )  THEN
2196                   !DIR$ IVDEP
2197                   DO  k = nzb+1, ind_even_odd
2198                      temp(k,j2,ixr) = p_mg(k,j,i)
2199                   ENDDO
2200                   j2 = j2 + 1
2201
2202                ENDIF
2203
2204             ENDDO
2205          ENDDO
2206
2207          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2208             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2209
2210                IF ( i == nxl_mg(l) )  THEN
2211                   !DIR$ IVDEP
2212                   DO  k = nzb+1, ind_even_odd
2213                      temp(k,j1,ixl) = p_mg(k,j,i)
2214                   ENDDO
2215                   j1 = j1 + 1
2216
2217                ENDIF
2218
2219                IF ( i == nxr_mg(l) )  THEN
2220                   !DIR$ IVDEP
2221                   DO  k = nzb+1, ind_even_odd
2222                      temp(k,j2,ixr) = p_mg(k,j,i)
2223                   ENDDO
2224                   j2 = j2 + 1
2225
2226                ENDIF
2227
2228             ENDDO
2229          ENDDO
2230
2231          grid_level = grid_level-1
2232
2233          send_receive = 'lr'
2234          CALL exchange_horiz( temp, 1 )
2235
2236          grid_level = grid_level+1
2237
2238          nxl = nxl_mg(grid_level)
2239          nys = nys_mg(grid_level)
2240          nxr = nxr_mg(grid_level)
2241          nyn = nyn_mg(grid_level)
2242          nzt = nzt_mg(grid_level)
2243
2244          j1 = nys_mg(grid_level-1)
2245          j2 = nys_mg(grid_level-1)
2246
2247          DO  i = nxl_mg(l), nxr_mg(l), 2
2248             DO  j = nys_mg(l) + (color-1), nyn_mg(l), 2
2249
2250                IF ( i == nxl_mg(l) )  THEN
2251                   !DIR$ IVDEP
2252                   DO  k = nzb+1, ind_even_odd
2253                      p_mg(k,j,nxr_mg(l)+1)  = temp(k,j1,ixr+1)
2254                   ENDDO
2255                   j1 = j1 + 1
2256
2257                ENDIF
2258
2259                IF ( i == nxr_mg(l) )  THEN
2260                   !DIR$ IVDEP
2261                   DO  k = nzb+1, ind_even_odd
2262                      p_mg(k,j,nxl_mg(l)-1)  = temp(k,j2,ixl-1)
2263                   ENDDO
2264                   j2 = j2 + 1
2265
2266                ENDIF
2267
2268             ENDDO
2269          ENDDO
2270
2271          DO  i = nxl_mg(l)+1, nxr_mg(l), 2
2272             DO j = nys_mg(l) + 2 - color, nyn_mg(l), 2
2273
2274                IF ( i == nxl_mg(l) )  THEN
2275                   !DIR$ IVDEP       INTEGER(iwp) ::  ngp       !<
2276
2277                   DO  k = nzb+1, ind_even_odd
2278                      p_mg(k,j,nxr_mg(l)+1)  = temp(k,j1,ixr+1)
2279                   ENDDO
2280                   j1 = j1 + 1
2281
2282                ENDIF
2283
2284                IF ( i == nxr_mg(l) )  THEN
2285                   !DIR$ IVDEP
2286                   DO  k = nzb+1, ind_even_odd
2287                      p_mg(k,j,nxl_mg(l)-1)  = temp(k,j2,ixl-1)
2288                   ENDDO
2289                   j2 = j2 + 1
2290
2291                ENDIF
2292
2293             ENDDO
2294          ENDDO
2295
2296          DEALLOCATE( temp )
2297
2298       ELSE
2299
2300!
2301!--       Standard horizontal ghost boundary exchange for small coarse grid
2302!--       levels, where the transfer time is latency bound
2303          CALL exchange_horiz( p_mg, 1 )
2304
2305       ENDIF
2306
2307!
2308!--    Reset values to default PALM setup
2309       synchronous_exchange   = synchronous_exchange_save
2310       send_receive = 'al'
2311#else
2312
2313!
2314!--    Next line is to avoid compile error due to unused dummy argument
2315       IF ( color == 1234567 )  RETURN
2316!
2317!--    Standard horizontal ghost boundary exchange for small coarse grid
2318!--    levels, where the transfer time is latency bound
2319       CALL exchange_horiz( p_mg, 1 )
2320#endif
2321
2322    END SUBROUTINE special_exchange_horiz
2323
2324 END MODULE poismg_mod
Note: See TracBrowser for help on using the repository browser.