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

Last change on this file since 3529 was 3524, checked in by raasch, 6 years ago

unused variables removed, missing working precision added, missing preprocessor directives added, bugfix concerning allocation of t_surf_wall_v in nopointer case, declaration statements rearranged to avoid compile time errors, mpi_abort arguments replaced to avoid compile errors

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