source: palm/trunk/SOURCE/init_grid.f90 @ 1779

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

pmc array management changed from linked list to sequential loop; further small changes and cosmetics for the pmc

  • Property svn:keywords set to Id
File size: 69.2 KB
Line 
1!> @file init_grid.f90
2!-------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2015 Leibniz Universitaet Hannover
17!-------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21! coupling_char is trimmed at every place it occurs, because it can have
22! different length now
23!
24! Former revisions:
25! -----------------
26! $Id: init_grid.f90 1779 2016-03-03 08:01:28Z raasch $
27!
28! 1762 2016-02-25 12:31:13Z hellstea
29! Introduction of nested domain feature
30!
31! 1743 2016-01-13 10:23:51Z raasch
32! Bugfix for calculation of nzb_s_outer and nzb_u_outer at north boundary of
33! total domain
34!
35! 1691 2015-10-26 16:17:44Z maronga
36! Renamed prandtl_layer to constant_flux_layer.
37!
38! 1682 2015-10-07 23:56:08Z knoop
39! Code annotations made doxygen readable
40!
41! 1677 2015-10-02 13:25:23Z boeske
42! Bugfix: Ghost points are included in wall_flags_0 and wall_flags_00
43!
44! 1675 2015-10-02 08:28:59Z gronemeier
45! Bugfix: Definition of topography grid levels
46!
47! 1660 2015-09-21 08:15:16Z gronemeier
48! Bugfix: Definition of topography grid levels if vertical grid stretching
49!         starts below the maximum topography height.
50!
51! 1580 2015-04-10 13:43:49Z suehring
52! Bugfix: setting flags for 5th order scheme near buildings
53!
54! 1575 2015-03-27 09:56:27Z raasch
55! adjustments for psolver-queries
56!
57! 1557 2015-03-05 16:43:04Z suehring
58! Adjustment for monotoinic limiter
59!
60! 1418 2014-06-06 13:05:08Z fricke
61! Bugfix: Change if-condition for stretched grid in the ocean, with the old
62!          condition and a negative value for dz_stretch_level the condition
63!          was always true for the whole model domain
64!
65! 1409 2014-05-23 12:11:32Z suehring
66! Bugfix: set wall_flags_0 at inflow and outflow boundary also for i <= nxlu
67! j <= nysv
68!
69! 1353 2014-04-08 15:21:23Z heinze
70! REAL constants provided with KIND-attribute
71!
72! 1322 2014-03-20 16:38:49Z raasch
73! REAL constants defined as wp-kind
74!
75! 1320 2014-03-20 08:40:49Z raasch
76! ONLY-attribute added to USE-statements,
77! kind-parameters added to all INTEGER and REAL declaration statements,
78! kinds are defined in new module kinds,
79! revision history before 2012 removed,
80! comment fields (!:) to be used for variable explanations added to
81! all variable declaration statements
82!
83! 1221 2013-09-10 08:59:13Z raasch
84! wall_flags_00 introduced to hold bits 32-63,
85! additional 3D-flag arrays for replacing the 2D-index array nzb_s_inner in
86! loops optimized for openACC (pres + flow_statistics)
87!
88! 1092 2013-02-02 11:24:22Z raasch
89! unused variables removed
90!
91! 1069 2012-11-28 16:18:43Z maronga
92! bugfix: added coupling_char to TOPOGRAPHY_DATA to allow topography in the
93!         ocean model in case of coupled runs
94!
95! 1036 2012-10-22 13:43:42Z raasch
96! code put under GPL (PALM 3.9)
97!
98! 1015 2012-09-27 09:23:24Z raasch
99! lower index for calculating wall_flags_0 set to nzb_w_inner instead of
100! nzb_w_inner+1
101!
102! 996 2012-09-07 10:41:47Z raasch
103! little reformatting
104!
105! 978 2012-08-09 08:28:32Z fricke
106! Bugfix: nzb_max is set to nzt at non-cyclic lateral boundaries
107! Bugfix: Set wall_flags_0 for inflow boundary
108!
109! 927 2012-06-06 19:15:04Z raasch
110! Wall flags are not set for multigrid method in case of masking method
111!
112! 864 2012-03-27 15:10:33Z gryschka
113! In case of ocean and Dirichlet bottom bc for u and v dzu_mg and ddzu_pres
114! were not correctly defined for k=1.
115!
116! 861 2012-03-26 14:18:34Z suehring
117! Set wall_flags_0. The array is needed for degradation in ws-scheme near walls,
118! inflow and outflow boundaries as well as near the bottom and the top of the
119! model domain.!
120! Initialization of nzb_s_inner and nzb_w_inner.
121! gls has to be at least nbgp to do not exceed the array bounds of nzb_local
122! while setting wall_flags_0
123!
124! 843 2012-02-29 15:16:21Z gryschka
125! In case of ocean and dirichlet bc for u and v at the bottom
126! the first u-level ist defined at same height as the first w-level
127!
128! 818 2012-02-08 16:11:23Z maronga
129! Bugfix: topo_height is only required if topography is used. It is thus now
130! allocated in the topography branch
131!
132! 809 2012-01-30 13:32:58Z maronga
133! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
134!
135! 807 2012-01-25 11:53:51Z maronga
136! New cpp directive "__check" implemented which is used by check_namelist_files
137!
138! Revision 1.1  1997/08/11 06:17:45  raasch
139! Initial revision (Testversion)
140!
141!
142! Description:
143! ------------
144!> Creating grid depending constants
145!------------------------------------------------------------------------------!
146 SUBROUTINE init_grid
147 
148
149    USE arrays_3d,                                                             &
150        ONLY:  dd2zu, ddzu, ddzu_pres, ddzw, dzu, dzu_mg, dzw, dzw_mg, f1_mg,  &
151               f2_mg, f3_mg, l_grid, l_wall, zu, zw
152       
153    USE control_parameters,                                                    &
154        ONLY:  bc_lr, bc_ns, building_height, building_length_x,               &
155               building_length_y, building_wall_left, building_wall_south,     &
156               canyon_height, canyon_wall_left, canyon_wall_south,             &
157               canyon_width_x, canyon_width_y, constant_flux_layer,            &
158               coupling_char, dp_level_ind_b, dz, dz_max, dz_stretch_factor,   &
159               dz_stretch_level, dz_stretch_level_index, ibc_uv_b, io_blocks,  &
160               io_group, inflow_l, inflow_n, inflow_r, inflow_s,               &
161               masking_method, maximum_grid_level, message_string,             &
162               momentum_advec, nest_domain, nest_bound_l, nest_bound_n,        &
163               nest_bound_r, nest_bound_s, ocean, outflow_l, outflow_n,        &
164               outflow_r, outflow_s, psolver, scalar_advec, topography,        &
165               topography_grid_convention, use_surface_fluxes, use_top_fluxes, &
166               wall_adjustment_factor 
167       
168    USE grid_variables,                                                        &
169        ONLY:  ddx, ddx2, ddx2_mg, ddy, ddy2, ddy2_mg, dx, dx2, dy, dy2, fwxm, &
170               fwxp, fwym, fwyp, fxm, fxp, fym, fyp, wall_e_x, wall_e_y,       &
171               wall_u, wall_v, wall_w_x, wall_w_y, zu_s_inner, zw_w_inner
172       
173    USE indices,                                                               &
174        ONLY:  flags, nbgp, nx, nxl, nxlg, nxlu, nxl_mg, nxr, nxrg, nxr_mg,    &
175               ny, nyn, nyng, nyn_mg, nys, nysv, nys_mg, nysg, nz, nzb,        &
176               nzb_2d, nzb_diff, nzb_diff_s_inner, nzb_diff_s_outer,           &
177               nzb_diff_u, nzb_diff_v, nzb_max, nzb_s_inner, nzb_s_outer,      &
178               nzb_u_inner, nzb_u_outer, nzb_v_inner, nzb_v_outer,             &
179               nzb_w_inner, nzb_w_outer, nzt, nzt_diff, nzt_mg, rflags_invers, &
180               rflags_s_inner, wall_flags_0, wall_flags_00, wall_flags_1,      &
181               wall_flags_10, wall_flags_2, wall_flags_3,  wall_flags_4,       &
182               wall_flags_5, wall_flags_6, wall_flags_7, wall_flags_8,         &
183               wall_flags_9
184   
185    USE kinds
186   
187    USE pegrid
188
189    IMPLICIT NONE
190
191    INTEGER(iwp) ::  bh      !<
192    INTEGER(iwp) ::  blx     !<
193    INTEGER(iwp) ::  bly     !<
194    INTEGER(iwp) ::  bxl     !<
195    INTEGER(iwp) ::  bxr     !<
196    INTEGER(iwp) ::  byn     !<
197    INTEGER(iwp) ::  bys     !<
198    INTEGER(iwp) ::  ch      !<
199    INTEGER(iwp) ::  cwx     !<
200    INTEGER(iwp) ::  cwy     !<
201    INTEGER(iwp) ::  cxl     !<
202    INTEGER(iwp) ::  cxr     !<
203    INTEGER(iwp) ::  cyn     !<
204    INTEGER(iwp) ::  cys     !<
205    INTEGER(iwp) ::  gls     !<
206    INTEGER(iwp) ::  i       !<
207    INTEGER(iwp) ::  ii      !<
208    INTEGER(iwp) ::  inc     !<
209    INTEGER(iwp) ::  j       !<
210    INTEGER(iwp) ::  k       !<
211    INTEGER(iwp) ::  l       !<
212    INTEGER(iwp) ::  nxl_l   !<
213    INTEGER(iwp) ::  nxr_l   !<
214    INTEGER(iwp) ::  nyn_l   !<
215    INTEGER(iwp) ::  nys_l   !<
216    INTEGER(iwp) ::  nzb_si  !<
217    INTEGER(iwp) ::  nzt_l   !<
218    INTEGER(iwp) ::  vi      !<
219
220    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::  vertical_influence  !<
221
222    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nl  !<
223    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nr  !<
224    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sl  !<
225    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sr  !<
226    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_l     !<
227    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_n     !<
228    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_r     !<
229    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_s     !<
230    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_local  !<
231    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_tmp    !<
232
233    LOGICAL :: flag_set = .FALSE.  !<
234
235    REAL(wp) ::  dx_l          !<
236    REAL(wp) ::  dy_l          !<
237    REAL(wp) ::  dz_stretched  !<
238
239    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  topo_height  !<
240
241   
242!
243!-- Calculation of horizontal array bounds including ghost layers
244    nxlg = nxl - nbgp
245    nxrg = nxr + nbgp
246    nysg = nys - nbgp
247    nyng = nyn + nbgp
248
249!
250!-- Allocate grid arrays
251    ALLOCATE( ddzu(1:nzt+1), ddzw(1:nzt+1), dd2zu(1:nzt), dzu(1:nzt+1),        &
252              dzw(1:nzt+1), l_grid(1:nzt), zu(nzb:nzt+1), zw(nzb:nzt+1) )
253
254!
255!-- Compute height of u-levels from constant grid length and dz stretch factors
256    IF ( dz == -1.0_wp )  THEN
257       message_string = 'missing dz'
258       CALL message( 'init_grid', 'PA0200', 1, 2, 0, 6, 0 ) 
259    ELSEIF ( dz <= 0.0_wp )  THEN
260       WRITE( message_string, * ) 'dz=',dz,' <= 0.0'
261       CALL message( 'init_grid', 'PA0201', 1, 2, 0, 6, 0 )
262    ENDIF
263
264!
265!-- Define the vertical grid levels
266    IF ( .NOT. ocean )  THEN
267!
268!--    Grid for atmosphere with surface at z=0 (k=0, w-grid).
269!--    The second u-level (k=1) corresponds to the top of the
270!--    Prandtl-layer.
271
272       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2 ) THEN
273          zu(0) = 0.0_wp
274      !    zu(0) = - dz * 0.5_wp
275       ELSE
276          zu(0) = - dz * 0.5_wp
277       ENDIF
278       zu(1) =   dz * 0.5_wp
279
280       dz_stretch_level_index = nzt+1
281       dz_stretched = dz
282       DO  k = 2, nzt+1
283          IF ( dz_stretch_level <= zu(k-1)  .AND.  dz_stretched < dz_max )  THEN
284             dz_stretched = dz_stretched * dz_stretch_factor
285             dz_stretched = MIN( dz_stretched, dz_max )
286             IF ( dz_stretch_level_index == nzt+1 ) dz_stretch_level_index = k-1
287          ENDIF
288          zu(k) = zu(k-1) + dz_stretched
289       ENDDO
290
291!
292!--    Compute the w-levels. They are always staggered half-way between the
293!--    corresponding u-levels. In case of dirichlet bc for u and v at the
294!--    ground the first u- and w-level (k=0) are defined at same height (z=0).
295!--    The top w-level is extrapolated linearly.
296       zw(0) = 0.0_wp
297       DO  k = 1, nzt
298          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
299       ENDDO
300       zw(nzt+1) = zw(nzt) + 2.0_wp * ( zu(nzt+1) - zw(nzt) )
301
302    ELSE
303!
304!--    Grid for ocean with free water surface is at k=nzt (w-grid).
305!--    In case of neumann bc at the ground the first first u-level (k=0) lies
306!--    below the first w-level (k=0). In case of dirichlet bc the first u- and
307!--    w-level are defined at same height, but staggered from the second level.
308!--    The second u-level (k=1) corresponds to the top of the Prandtl-layer.
309       zu(nzt+1) =   dz * 0.5_wp
310       zu(nzt)   = - dz * 0.5_wp
311
312       dz_stretch_level_index = 0
313       dz_stretched = dz
314       DO  k = nzt-1, 0, -1
315!
316!--       The default value of dz_stretch_level is positive, thus the first
317!--       condition is always true. Hence, the second condition is necessary.
318          IF ( dz_stretch_level >= zu(k+1)  .AND.  dz_stretch_level <= 0.0  &
319               .AND.  dz_stretched < dz_max )  THEN
320             dz_stretched = dz_stretched * dz_stretch_factor
321             dz_stretched = MIN( dz_stretched, dz_max )
322             IF ( dz_stretch_level_index == 0 ) dz_stretch_level_index = k+1
323          ENDIF
324          zu(k) = zu(k+1) - dz_stretched
325       ENDDO
326
327!
328!--    Compute the w-levels. They are always staggered half-way between the
329!--    corresponding u-levels, except in case of dirichlet bc for u and v
330!--    at the ground. In this case the first u- and w-level are defined at
331!--    same height. The top w-level (nzt+1) is not used but set for
332!--    consistency, since w and all scalar variables are defined up tp nzt+1.
333       zw(nzt+1) = dz
334       zw(nzt)   = 0.0_wp
335       DO  k = 0, nzt
336          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
337       ENDDO
338
339!
340!--    In case of dirichlet bc for u and v the first u- and w-level are defined
341!--    at same height.
342       IF ( ibc_uv_b == 0 ) THEN
343          zu(0) = zw(0)
344       ENDIF
345
346    ENDIF
347
348!
349!-- Compute grid lengths.
350    DO  k = 1, nzt+1
351       dzu(k)  = zu(k) - zu(k-1)
352       ddzu(k) = 1.0_wp / dzu(k)
353       dzw(k)  = zw(k) - zw(k-1)
354       ddzw(k) = 1.0_wp / dzw(k)
355    ENDDO
356
357    DO  k = 1, nzt
358       dd2zu(k) = 1.0_wp / ( dzu(k) + dzu(k+1) )
359    ENDDO
360   
361!   
362!-- The FFT- SOR-pressure solvers assume grid spacings of a staggered grid
363!-- everywhere. For the actual grid, the grid spacing at the lowest level
364!-- is only dz/2, but should be dz. Therefore, an additional array
365!-- containing with appropriate grid information is created for these
366!-- solvers.
367    IF ( psolver(1:9) /= 'multigrid' )  THEN
368       ALLOCATE( ddzu_pres(1:nzt+1) )
369       ddzu_pres = ddzu
370       ddzu_pres(1) = ddzu_pres(2)  ! change for lowest level
371    ENDIF   
372
373!
374!-- In case of multigrid method, compute grid lengths and grid factors for the
375!-- grid levels
376    IF ( psolver(1:9) == 'multigrid' )  THEN
377
378       ALLOCATE( ddx2_mg(maximum_grid_level), ddy2_mg(maximum_grid_level), &
379                 dzu_mg(nzb+1:nzt+1,maximum_grid_level),                   &
380                 dzw_mg(nzb+1:nzt+1,maximum_grid_level),                   &
381                 f1_mg(nzb+1:nzt,maximum_grid_level),                      &
382                 f2_mg(nzb+1:nzt,maximum_grid_level),                      &
383                 f3_mg(nzb+1:nzt,maximum_grid_level) )
384
385       dzu_mg(:,maximum_grid_level) = dzu
386!       
387!--    Next line to ensure an equally spaced grid.
388       dzu_mg(1,maximum_grid_level) = dzu(2)
389
390       dzw_mg(:,maximum_grid_level) = dzw
391       nzt_l = nzt
392       DO  l = maximum_grid_level-1, 1, -1
393           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
394           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
395           nzt_l = nzt_l / 2
396           DO  k = 2, nzt_l+1
397              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
398              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
399           ENDDO
400       ENDDO
401
402       nzt_l = nzt
403       dx_l  = dx
404       dy_l  = dy
405       DO  l = maximum_grid_level, 1, -1
406          ddx2_mg(l) = 1.0_wp / dx_l**2
407          ddy2_mg(l) = 1.0_wp / dy_l**2
408          DO  k = nzb+1, nzt_l
409             f2_mg(k,l) = 1.0_wp / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
410             f3_mg(k,l) = 1.0_wp / ( dzu_mg(k,l)   * dzw_mg(k,l) )
411             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) + &
412                          f2_mg(k,l) + f3_mg(k,l)
413          ENDDO
414          nzt_l = nzt_l / 2
415          dx_l  = dx_l * 2.0_wp
416          dy_l  = dy_l * 2.0_wp
417       ENDDO
418
419    ENDIF
420
421!
422!-- Compute the reciprocal values of the horizontal grid lengths.
423    ddx = 1.0_wp / dx
424    ddy = 1.0_wp / dy
425    dx2 = dx * dx
426    dy2 = dy * dy
427    ddx2 = 1.0_wp / dx2
428    ddy2 = 1.0_wp / dy2
429
430!
431!-- Compute the grid-dependent mixing length.
432    DO  k = 1, nzt
433       l_grid(k)  = ( dx * dy * dzw(k) )**0.33333333333333_wp
434    ENDDO
435
436!
437!-- Allocate outer and inner index arrays for topography and set
438!-- defaults.
439!-- nzb_local has to contain additional layers of ghost points for calculating
440!-- the flag arrays needed for the multigrid method
441    gls = 2**( maximum_grid_level )
442    IF ( gls < nbgp )  gls = nbgp
443
444    ALLOCATE( corner_nl(nys:nyn,nxl:nxr), corner_nr(nys:nyn,nxl:nxr),       &
445              corner_sl(nys:nyn,nxl:nxr), corner_sr(nys:nyn,nxl:nxr),       &
446              nzb_local(-gls:ny+gls,-gls:nx+gls),                                   &
447              nzb_tmp(-nbgp:ny+nbgp,-nbgp:nx+nbgp),                         &
448              wall_l(nys:nyn,nxl:nxr), wall_n(nys:nyn,nxl:nxr),             &
449              wall_r(nys:nyn,nxl:nxr), wall_s(nys:nyn,nxl:nxr) )
450    ALLOCATE( fwxm(nysg:nyng,nxlg:nxrg), fwxp(nysg:nyng,nxlg:nxrg),         &
451              fwym(nysg:nyng,nxlg:nxrg), fwyp(nysg:nyng,nxlg:nxrg),         &
452              fxm(nysg:nyng,nxlg:nxrg), fxp(nysg:nyng,nxlg:nxrg),           &
453              fym(nysg:nyng,nxlg:nxrg), fyp(nysg:nyng,nxlg:nxrg),           &
454              nzb_s_inner(nysg:nyng,nxlg:nxrg),                             &
455              nzb_s_outer(nysg:nyng,nxlg:nxrg),                             &
456              nzb_u_inner(nysg:nyng,nxlg:nxrg),                             &
457              nzb_u_outer(nysg:nyng,nxlg:nxrg),                             &
458              nzb_v_inner(nysg:nyng,nxlg:nxrg),                             &
459              nzb_v_outer(nysg:nyng,nxlg:nxrg),                             &
460              nzb_w_inner(nysg:nyng,nxlg:nxrg),                             &
461              nzb_w_outer(nysg:nyng,nxlg:nxrg),                             &
462              nzb_diff_s_inner(nysg:nyng,nxlg:nxrg),                        &
463              nzb_diff_s_outer(nysg:nyng,nxlg:nxrg),                        &
464              nzb_diff_u(nysg:nyng,nxlg:nxrg),                              &
465              nzb_diff_v(nysg:nyng,nxlg:nxrg),                              &
466              nzb_2d(nysg:nyng,nxlg:nxrg),                                  &
467              rflags_s_inner(nzb:nzt+2,nysg:nyng,nxlg:nxrg),                &
468              rflags_invers(nysg:nyng,nxlg:nxrg,nzb:nzt+2),                 &
469              wall_e_x(nysg:nyng,nxlg:nxrg),                                &
470              wall_e_y(nysg:nyng,nxlg:nxrg),                                &
471              wall_u(nysg:nyng,nxlg:nxrg),                                  &
472              wall_v(nysg:nyng,nxlg:nxrg),                                  &
473              wall_w_x(nysg:nyng,nxlg:nxrg),                                &
474              wall_w_y(nysg:nyng,nxlg:nxrg) )
475
476
477
478    ALLOCATE( l_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
479
480
481    nzb_s_inner = nzb;  nzb_s_outer = nzb
482    nzb_u_inner = nzb;  nzb_u_outer = nzb
483    nzb_v_inner = nzb;  nzb_v_outer = nzb
484    nzb_w_inner = nzb;  nzb_w_outer = nzb
485
486    rflags_s_inner = 1.0_wp
487    rflags_invers  = 1.0_wp
488
489!
490!-- Define vertical gridpoint from (or to) which on the usual finite difference
491!-- form (which does not use surface fluxes) is applied
492    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
493       nzb_diff = nzb + 2
494    ELSE
495       nzb_diff = nzb + 1
496    ENDIF
497    IF ( use_top_fluxes )  THEN
498       nzt_diff = nzt - 1
499    ELSE
500       nzt_diff = nzt
501    ENDIF
502
503    nzb_diff_s_inner = nzb_diff;  nzb_diff_s_outer = nzb_diff
504    nzb_diff_u = nzb_diff;  nzb_diff_v = nzb_diff
505
506    wall_e_x = 0.0_wp;  wall_e_y = 0.0_wp;  wall_u = 0.0_wp;  wall_v = 0.0_wp
507    wall_w_x = 0.0_wp;  wall_w_y = 0.0_wp
508    fwxp = 1.0_wp;  fwxm = 1.0_wp;  fwyp = 1.0_wp;  fwym = 1.0_wp
509    fxp  = 1.0_wp;  fxm  = 1.0_wp;  fyp  = 1.0_wp;  fym  = 1.0_wp
510
511!
512!-- Initialize near-wall mixing length l_wall only in the vertical direction
513!-- for the moment,
514!-- multiplication with wall_adjustment_factor near the end of this routine
515    l_wall(nzb,:,:)   = l_grid(1)
516    DO  k = nzb+1, nzt
517       l_wall(k,:,:)  = l_grid(k)
518    ENDDO
519    l_wall(nzt+1,:,:) = l_grid(nzt)
520
521    ALLOCATE ( vertical_influence(nzb:nzt) )
522    DO  k = 1, nzt
523       vertical_influence(k) = MIN ( INT( l_grid(k) / &
524                     ( wall_adjustment_factor * dzw(k) ) + 0.5_wp ), nzt - k )
525    ENDDO
526
527    DO  k = 1, MAXVAL( nzb_s_inner )
528       IF ( l_grid(k) > 1.5_wp * dx * wall_adjustment_factor .OR.  &
529            l_grid(k) > 1.5_wp * dy * wall_adjustment_factor )  THEN
530          WRITE( message_string, * ) 'grid anisotropy exceeds ', &
531                                     'threshold given by only local', &
532                                     ' &horizontal reduction of near_wall ', &
533                                     'mixing length l_wall', &
534                                     ' &starting from height level k = ', k, '.'
535          CALL message( 'init_grid', 'PA0202', 0, 1, 0, 6, 0 )
536          EXIT
537       ENDIF
538    ENDDO
539    vertical_influence(0) = vertical_influence(1)
540
541    DO  i = nxlg, nxrg
542       DO  j = nysg, nyng
543          DO  k = nzb_s_inner(j,i) + 1, &
544                  nzb_s_inner(j,i) + vertical_influence(nzb_s_inner(j,i))
545             l_wall(k,j,i) = zu(k) - zw(nzb_s_inner(j,i))
546          ENDDO
547       ENDDO
548    ENDDO
549
550!
551!-- Set outer and inner index arrays for non-flat topography.
552!-- Here consistency checks concerning domain size and periodicity are
553!-- necessary.
554!-- Within this SELECT CASE structure only nzb_local is initialized
555!-- individually depending on the chosen topography type, all other index
556!-- arrays are initialized further below.
557    SELECT CASE ( TRIM( topography ) )
558
559       CASE ( 'flat' )
560!
561!--       nzb_local is required for the multigrid solver
562          nzb_local = 0
563
564       CASE ( 'single_building' )
565!
566!--       Single rectangular building, by default centered in the middle of the
567!--       total domain
568          blx = NINT( building_length_x / dx )
569          bly = NINT( building_length_y / dy )
570          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
571          IF ( ABS( zw(bh  ) - building_height ) == &
572               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
573
574          IF ( building_wall_left == 9999999.9_wp )  THEN
575             building_wall_left = ( nx + 1 - blx ) / 2 * dx
576          ENDIF
577          bxl = NINT( building_wall_left / dx )
578          bxr = bxl + blx
579
580          IF ( building_wall_south == 9999999.9_wp )  THEN
581             building_wall_south = ( ny + 1 - bly ) / 2 * dy
582          ENDIF
583          bys = NINT( building_wall_south / dy )
584          byn = bys + bly
585
586!
587!--       Building size has to meet some requirements
588          IF ( ( bxl < 1 ) .OR. ( bxr > nx-1 ) .OR. ( bxr < bxl+3 ) .OR.  &
589               ( bys < 1 ) .OR. ( byn > ny-1 ) .OR. ( byn < bys+3 ) )  THEN
590             WRITE( message_string, * ) 'inconsistent building parameters:',   &
591                                      '& bxl=', bxl, 'bxr=', bxr, 'bys=', bys, &
592                                      'byn=', byn, 'nx=', nx, 'ny=', ny
593             CALL message( 'init_grid', 'PA0203', 1, 2, 0, 6, 0 )
594          ENDIF
595
596!
597!--       Define the building.
598          nzb_local = 0
599          nzb_local(bys:byn,bxl:bxr) = bh
600
601       CASE ( 'single_street_canyon' )
602!
603!--       Single quasi-2D street canyon of infinite length in x or y direction.
604!--       The canyon is centered in the other direction by default.
605          IF ( canyon_width_x /= 9999999.9_wp )  THEN
606!
607!--          Street canyon in y direction
608             cwx = NINT( canyon_width_x / dx )
609             IF ( canyon_wall_left == 9999999.9_wp )  THEN
610                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
611             ENDIF
612             cxl = NINT( canyon_wall_left / dx )
613             cxr = cxl + cwx
614
615          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
616!
617!--          Street canyon in x direction
618             cwy = NINT( canyon_width_y / dy )
619             IF ( canyon_wall_south == 9999999.9_wp )  THEN
620                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
621             ENDIF
622             cys = NINT( canyon_wall_south / dy )
623             cyn = cys + cwy
624
625          ELSE
626             
627             message_string = 'no street canyon width given'
628             CALL message( 'init_grid', 'PA0204', 1, 2, 0, 6, 0 )
629 
630          ENDIF
631
632          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
633          IF ( ABS( zw(ch  ) - canyon_height ) == &
634               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
635
636          dp_level_ind_b = ch
637!
638!--       Street canyon size has to meet some requirements
639          IF ( canyon_width_x /= 9999999.9_wp )  THEN
640             IF ( ( cxl < 1 ) .OR. ( cxr > nx-1 ) .OR. ( cwx < 3 ) .OR.        &
641               ( ch < 3 ) )  THEN
642                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
643                                           '&cxl=', cxl, 'cxr=', cxr,          &
644                                           'cwx=', cwx,                        &
645                                           'ch=', ch, 'nx=', nx, 'ny=', ny
646                CALL message( 'init_grid', 'PA0205', 1, 2, 0, 6, 0 ) 
647             ENDIF
648          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
649             IF ( ( cys < 1 ) .OR. ( cyn > ny-1 ) .OR. ( cwy < 3 ) .OR.        &
650               ( ch < 3 ) )  THEN
651                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
652                                           '&cys=', cys, 'cyn=', cyn,          &
653                                           'cwy=', cwy,                        &
654                                           'ch=', ch, 'nx=', nx, 'ny=', ny
655                CALL message( 'init_grid', 'PA0206', 1, 2, 0, 6, 0 ) 
656             ENDIF
657          ENDIF
658          IF ( canyon_width_x /= 9999999.9_wp .AND.                            &                 
659               canyon_width_y /= 9999999.9_wp )  THEN
660             message_string = 'inconsistent canyon parameters:' //             &   
661                              '&street canyon can only be oriented' //         &
662                              '&either in x- or in y-direction'
663             CALL message( 'init_grid', 'PA0207', 1, 2, 0, 6, 0 )
664          ENDIF
665
666          nzb_local = ch
667          IF ( canyon_width_x /= 9999999.9_wp )  THEN
668             nzb_local(:,cxl+1:cxr-1) = 0
669          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
670             nzb_local(cys+1:cyn-1,:) = 0
671          ENDIF
672
673       CASE ( 'read_from_file' )
674
675          ALLOCATE ( topo_height(0:ny,0:nx) )
676
677          DO  ii = 0, io_blocks-1
678             IF ( ii == io_group )  THEN
679
680!
681!--             Arbitrary irregular topography data in PALM format (exactly
682!--             matching the grid size and total domain size)
683                OPEN( 90, FILE='TOPOGRAPHY_DATA'//TRIM( coupling_char ),       &
684                          STATUS='OLD', FORM='FORMATTED', ERR=10 )
685                DO  j = ny, 0, -1
686                   READ( 90, *, ERR=11, END=11 )  ( topo_height(j,i), i = 0,nx )
687                ENDDO
688
689                GOTO 12
690         
691 10             message_string = 'file TOPOGRAPHY'//TRIM( coupling_char )//    &
692                                 ' does not exist'
693                CALL message( 'init_grid', 'PA0208', 1, 2, 0, 6, 0 )
694
695 11             message_string = 'errors in file TOPOGRAPHY_DATA'//            &
696                                 TRIM( coupling_char )
697                CALL message( 'init_grid', 'PA0209', 1, 2, 0, 6, 0 )
698
699 12             CLOSE( 90 )
700
701             ENDIF
702#if defined( __parallel ) && ! defined ( __check )
703             CALL MPI_BARRIER( comm2d, ierr )
704#endif
705          ENDDO
706
707!
708!--       Calculate the index height of the topography
709          DO  i = 0, nx
710             DO  j = 0, ny
711                nzb_local(j,i) = MINLOC( ABS( zw - topo_height(j,i) ), 1 ) - 1
712                IF ( ABS( zw(nzb_local(j,i)  ) - topo_height(j,i) ) == &
713                     ABS( zw(nzb_local(j,i)+1) - topo_height(j,i) )    )  &
714                   nzb_local(j,i) = nzb_local(j,i) + 1
715             ENDDO
716          ENDDO
717
718          DEALLOCATE ( topo_height )
719!
720!--       Add cyclic boundaries (additional layers are for calculating
721!--       flag arrays needed for the multigrid sover)
722          nzb_local(-gls:-1,0:nx)     = nzb_local(ny-gls+1:ny,0:nx)
723          nzb_local(ny+1:ny+gls,0:nx) = nzb_local(0:gls-1,0:nx)
724          nzb_local(:,-gls:-1)        = nzb_local(:,nx-gls+1:nx)
725          nzb_local(:,nx+1:nx+gls)    = nzb_local(:,0:gls-1)
726
727       CASE DEFAULT
728!
729!--       The DEFAULT case is reached either if the parameter topography
730!--       contains a wrong character string or if the user has defined a special
731!--       case in the user interface. There, the subroutine user_init_grid
732!--       checks which of these two conditions applies.
733          CALL user_init_grid( gls, nzb_local )
734
735    END SELECT
736!
737!-- Determine the maximum level of topography. Furthermore it is used for
738!-- steering the degradation of order of the applied advection scheme.
739!-- In case of non-cyclic lateral boundaries, the order of the advection
740!-- scheme have to be reduced up to nzt (required at the lateral boundaries).
741    nzb_max = MAXVAL( nzb_local ) + 1
742    IF ( inflow_l .OR. outflow_l .OR. inflow_r .OR. outflow_r .OR.             &
743         inflow_n .OR. outflow_n .OR. inflow_s .OR. outflow_s .OR.             &
744         nest_domain )                                                         &
745    THEN
746       nzb_max = nzt
747    ENDIF
748
749!
750!-- Consistency checks and index array initialization are only required for
751!-- non-flat topography, also the initialization of topography height arrays
752!-- zu_s_inner and zw_w_inner
753    IF ( TRIM( topography ) /= 'flat' )  THEN
754
755!
756!--    Consistency checks
757       IF ( MINVAL( nzb_local ) < 0  .OR.  MAXVAL( nzb_local ) > nz + 1 )  THEN
758          WRITE( message_string, * ) 'nzb_local values are outside the',       &
759                                'model domain',                                &
760                                '&MINVAL( nzb_local ) = ', MINVAL(nzb_local),  &
761                                '&MAXVAL( nzb_local ) = ', MAXVAL(nzb_local)
762          CALL message( 'init_grid', 'PA0210', 1, 2, 0, 6, 0 )
763       ENDIF
764
765       IF ( bc_lr == 'cyclic' )  THEN
766          IF ( ANY( nzb_local(:,-1) /= nzb_local(:,nx)   )  .OR.               &
767               ANY( nzb_local(:,0)  /= nzb_local(:,nx+1) ) )  THEN
768             message_string = 'nzb_local does not fulfill cyclic' //           &
769                              ' boundary condition in x-direction'
770             CALL message( 'init_grid', 'PA0211', 1, 2, 0, 6, 0 )
771          ENDIF
772       ENDIF
773       IF ( bc_ns == 'cyclic' )  THEN
774          IF ( ANY( nzb_local(-1,:) /= nzb_local(ny,:)   )  .OR.               &
775               ANY( nzb_local(0,:)  /= nzb_local(ny+1,:) ) )  THEN
776             message_string = 'nzb_local does not fulfill cyclic' //           &
777                              ' boundary condition in y-direction'
778             CALL message( 'init_grid', 'PA0212', 1, 2, 0, 6, 0 )
779          ENDIF
780       ENDIF
781
782       IF ( topography_grid_convention == 'cell_edge' )  THEN
783!
784!--       The array nzb_local as defined using the 'cell_edge' convention
785!--       describes the actual total size of topography which is defined at the
786!--       cell edges where u=0 on the topography walls in x-direction and v=0
787!--       on the topography walls in y-direction. However, PALM uses individual
788!--       arrays nzb_u|v|w|s_inner|outer that are based on nzb_s_inner.
789!--       Therefore, the extent of topography in nzb_local is now reduced by
790!--       1dx at the E topography walls and by 1dy at the N topography walls
791!--       to form the basis for nzb_s_inner.
792          DO  j = -gls, ny + gls
793             DO  i = -gls, nx
794                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j,i+1) )
795             ENDDO
796          ENDDO
797!--       apply cyclic boundary conditions in x-direction
798!(ist das erforderlich? Ursache von Seung Bus Fehler?)
799          nzb_local(:,nx+1:nx+gls) = nzb_local(:,0:gls-1)
800          DO  i = -gls, nx + gls
801             DO  j = -gls, ny
802                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j+1,i) )
803             ENDDO
804          ENDDO
805!--       apply cyclic boundary conditions in y-direction
806!(ist das erforderlich? Ursache von Seung Bus Fehler?)
807          nzb_local(ny+1:ny+gls,:) = nzb_local(0:gls-1,:)
808       ENDIF
809
810!
811!--    Initialize index arrays nzb_s_inner and nzb_w_inner
812       nzb_s_inner = nzb_local(nysg:nyng,nxlg:nxrg)
813       nzb_w_inner = nzb_local(nysg:nyng,nxlg:nxrg)
814
815!
816!--    Initialize remaining index arrays:
817!--    first pre-initialize them with nzb_s_inner...
818       nzb_u_inner = nzb_s_inner
819       nzb_u_outer = nzb_s_inner
820       nzb_v_inner = nzb_s_inner
821       nzb_v_outer = nzb_s_inner
822       nzb_w_outer = nzb_s_inner
823       nzb_s_outer = nzb_s_inner
824
825!
826!--    ...then extend pre-initialized arrays in their according directions
827!--    based on nzb_local using nzb_tmp as a temporary global index array
828
829!
830!--    nzb_s_outer:
831!--    extend nzb_local east-/westwards first, then north-/southwards
832       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
833       DO  j = -1, ny + 1
834          DO  i = 0, nx
835             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i),             &
836                                 nzb_local(j,i+1) )
837          ENDDO
838       ENDDO
839       DO  i = nxl, nxr
840          DO  j = nys, nyn
841             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
842                                     nzb_tmp(j+1,i) )
843          ENDDO
844!
845!--       non-cyclic boundary conditions (overwritten by call of
846!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
847          IF ( nys == 0 )  THEN
848             j = -1
849             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
850          ENDIF
851          IF ( nyn == ny )  THEN
852             j = ny + 1
853             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
854          ENDIF
855       ENDDO
856!
857!--    nzb_w_outer:
858!--    identical to nzb_s_outer
859       nzb_w_outer = nzb_s_outer
860
861!
862!--    nzb_u_inner:
863!--    extend nzb_local rightwards only
864       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
865       DO  j = -1, ny + 1
866          DO  i = 0, nx + 1
867             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
868          ENDDO
869       ENDDO
870       nzb_u_inner = nzb_tmp(nysg:nyng,nxlg:nxrg)
871
872!
873!--    nzb_u_outer:
874!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
875       DO  i = nxl, nxr
876          DO  j = nys, nyn
877             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
878                                     nzb_tmp(j+1,i) )
879          ENDDO
880!
881!--       non-cyclic boundary conditions (overwritten by call of
882!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
883          IF ( nys == 0 )  THEN
884             j = -1
885             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
886          ENDIF
887          IF ( nyn == ny )  THEN
888             j = ny + 1
889             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
890          ENDIF
891       ENDDO
892
893!
894!--    nzb_v_inner:
895!--    extend nzb_local northwards only
896       nzb_tmp = nzb_local(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
897       DO  i = -1, nx + 1
898          DO  j = 0, ny + 1
899             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
900          ENDDO
901       ENDDO
902       nzb_v_inner = nzb_tmp(nys-nbgp:nyn+nbgp,nxl-nbgp:nxr+nbgp)
903
904!
905!--    nzb_v_outer:
906!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
907       DO  j = nys, nyn
908          DO  i = nxl, nxr
909             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i),             &
910                                     nzb_tmp(j,i+1) )
911          ENDDO
912!
913!--       non-cyclic boundary conditions (overwritten by call of
914!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
915          IF ( nxl == 0 )  THEN
916             i = -1
917             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
918          ENDIF
919          IF ( nxr == nx )  THEN
920             i = nx + 1
921             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
922          ENDIF
923       ENDDO
924#if ! defined ( __check )
925!
926!--    Exchange of lateral boundary values (parallel computers) and cyclic
927!--    boundary conditions, if applicable.
928!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
929!--    they do not require exchange and are not included here.
930       CALL exchange_horiz_2d_int( nzb_u_inner )
931       CALL exchange_horiz_2d_int( nzb_u_outer )
932       CALL exchange_horiz_2d_int( nzb_v_inner )
933       CALL exchange_horiz_2d_int( nzb_v_outer )
934       CALL exchange_horiz_2d_int( nzb_w_outer )
935       CALL exchange_horiz_2d_int( nzb_s_outer )
936
937!
938!--    Allocate and set the arrays containing the topography height
939       IF ( myid == 0 )  THEN
940
941          ALLOCATE( zu_s_inner(0:nx+1,0:ny+1), zw_w_inner(0:nx+1,0:ny+1) )
942
943          DO  i = 0, nx + 1
944             DO  j = 0, ny + 1
945                zu_s_inner(i,j) = zu(nzb_local(j,i))
946                zw_w_inner(i,j) = zw(nzb_local(j,i))
947             ENDDO
948          ENDDO
949         
950       ENDIF
951!
952!--    Set flag arrays to be used for masking of grid points
953       DO  i = nxlg, nxrg
954          DO  j = nysg, nyng
955             DO  k = nzb, nzt+1
956                IF ( k <= nzb_s_inner(j,i) )  rflags_s_inner(k,j,i) = 0.0_wp
957                IF ( k <= nzb_s_inner(j,i) )  rflags_invers(j,i,k)  = 0.0_wp
958             ENDDO
959          ENDDO
960       ENDDO
961#endif
962    ENDIF
963
964#if ! defined ( __check )
965!
966!-- Preliminary: to be removed after completion of the topography code!
967!-- Set the former default k index arrays nzb_2d
968    nzb_2d      = nzb
969
970!
971!-- Set the individual index arrays which define the k index from which on
972!-- the usual finite difference form (which does not use surface fluxes) is
973!-- applied
974    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
975       nzb_diff_u         = nzb_u_inner + 2
976       nzb_diff_v         = nzb_v_inner + 2
977       nzb_diff_s_inner   = nzb_s_inner + 2
978       nzb_diff_s_outer   = nzb_s_outer + 2
979    ELSE
980       nzb_diff_u         = nzb_u_inner + 1
981       nzb_diff_v         = nzb_v_inner + 1
982       nzb_diff_s_inner   = nzb_s_inner + 1
983       nzb_diff_s_outer   = nzb_s_outer + 1
984    ENDIF
985
986!
987!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
988!-- for limitation of near-wall mixing length l_wall further below
989    corner_nl = 0
990    corner_nr = 0
991    corner_sl = 0
992    corner_sr = 0
993    wall_l    = 0
994    wall_n    = 0
995    wall_r    = 0
996    wall_s    = 0
997
998    DO  i = nxl, nxr
999       DO  j = nys, nyn
1000!
1001!--       u-component
1002          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
1003             wall_u(j,i) = 1.0_wp   ! north wall (location of adjacent fluid)
1004             fym(j,i)    = 0.0_wp
1005             fyp(j,i)    = 1.0_wp
1006          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
1007             wall_u(j,i) = 1.0_wp   ! south wall (location of adjacent fluid)
1008             fym(j,i)    = 1.0_wp
1009             fyp(j,i)    = 0.0_wp
1010          ENDIF
1011!
1012!--       v-component
1013          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
1014             wall_v(j,i) = 1.0_wp   ! rigth wall (location of adjacent fluid)
1015             fxm(j,i)    = 0.0_wp
1016             fxp(j,i)    = 1.0_wp
1017          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
1018             wall_v(j,i) = 1.0_wp   ! left wall (location of adjacent fluid)
1019             fxm(j,i)    = 1.0_wp
1020             fxp(j,i)    = 0.0_wp
1021          ENDIF
1022!
1023!--       w-component, also used for scalars, separate arrays for shear
1024!--       production of tke
1025          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
1026             wall_e_y(j,i) =  1.0_wp   ! north wall (location of adjacent fluid)
1027             wall_w_y(j,i) =  1.0_wp
1028             fwym(j,i)     =  0.0_wp
1029             fwyp(j,i)     =  1.0_wp
1030          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
1031             wall_e_y(j,i) = -1.0_wp   ! south wall (location of adjacent fluid)
1032             wall_w_y(j,i) =  1.0_wp
1033             fwym(j,i)     =  1.0_wp
1034             fwyp(j,i)     =  0.0_wp
1035          ENDIF
1036          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
1037             wall_e_x(j,i) =  1.0_wp   ! right wall (location of adjacent fluid)
1038             wall_w_x(j,i) =  1.0_wp
1039             fwxm(j,i)     =  0.0_wp
1040             fwxp(j,i)     =  1.0_wp
1041          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
1042             wall_e_x(j,i) = -1.0_wp   ! left wall (location of adjacent fluid)
1043             wall_w_x(j,i) =  1.0_wp
1044             fwxm(j,i)     =  1.0_wp
1045             fwxp(j,i)     =  0.0_wp
1046          ENDIF
1047!
1048!--       Wall and corner locations inside buildings for limitation of
1049!--       near-wall mixing length l_wall
1050          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
1051
1052             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
1053
1054             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1055                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
1056                                      nzb_s_inner(j,i-1) ) + 1
1057             ENDIF
1058
1059             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1060                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
1061                                      nzb_s_inner(j,i+1) ) + 1
1062             ENDIF
1063
1064          ENDIF
1065
1066          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
1067
1068             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
1069             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1070                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
1071                                      nzb_s_inner(j,i-1) ) + 1
1072             ENDIF
1073
1074             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1075                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
1076                                      nzb_s_inner(j,i+1) ) + 1
1077             ENDIF
1078
1079          ENDIF
1080
1081          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1082             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
1083          ENDIF
1084
1085          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1086             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
1087          ENDIF
1088
1089       ENDDO
1090    ENDDO
1091
1092!
1093!-- Calculate wall flag arrays for the multigrid method
1094    IF ( psolver(1:9) == 'multigrid' )  THEN
1095!
1096!--    Gridpoint increment of the current level
1097       inc = 1
1098
1099       DO  l = maximum_grid_level, 1 , -1
1100
1101          nxl_l = nxl_mg(l)
1102          nxr_l = nxr_mg(l)
1103          nys_l = nys_mg(l)
1104          nyn_l = nyn_mg(l)
1105          nzt_l = nzt_mg(l)
1106
1107!
1108!--       Assign the flag level to be calculated
1109          SELECT CASE ( l )
1110             CASE ( 1 )
1111                flags => wall_flags_1
1112             CASE ( 2 )
1113                flags => wall_flags_2
1114             CASE ( 3 )
1115                flags => wall_flags_3
1116             CASE ( 4 )
1117                flags => wall_flags_4
1118             CASE ( 5 )
1119                flags => wall_flags_5
1120             CASE ( 6 )
1121                flags => wall_flags_6
1122             CASE ( 7 )
1123                flags => wall_flags_7
1124             CASE ( 8 )
1125                flags => wall_flags_8
1126             CASE ( 9 )
1127                flags => wall_flags_9
1128             CASE ( 10 )
1129                flags => wall_flags_10
1130          END SELECT
1131
1132!
1133!--       Depending on the grid level, set the respective bits in case of
1134!--       neighbouring walls
1135!--       Bit 0:  wall to the bottom
1136!--       Bit 1:  wall to the top (not realized in remaining PALM code so far)
1137!--       Bit 2:  wall to the south
1138!--       Bit 3:  wall to the north
1139!--       Bit 4:  wall to the left
1140!--       Bit 5:  wall to the right
1141!--       Bit 6:  inside building
1142
1143          flags = 0
1144
1145!
1146!--       In case of masking method, flags are not set and multigrid method
1147!--       works like FFT-solver
1148          IF ( .NOT. masking_method )  THEN
1149
1150             DO  i = nxl_l-1, nxr_l+1
1151                DO  j = nys_l-1, nyn_l+1
1152                   DO  k = nzb, nzt_l+1
1153                         
1154!
1155!--                   Inside/outside building (inside building does not need
1156!--                   further tests for walls)
1157                      IF ( k*inc <= nzb_local(j*inc,i*inc) )  THEN
1158
1159                         flags(k,j,i) = IBSET( flags(k,j,i), 6 )
1160
1161                      ELSE
1162!
1163!--                      Bottom wall
1164                         IF ( (k-1)*inc <= nzb_local(j*inc,i*inc) )  THEN
1165                            flags(k,j,i) = IBSET( flags(k,j,i), 0 )
1166                         ENDIF
1167!
1168!--                      South wall
1169                         IF ( k*inc <= nzb_local((j-1)*inc,i*inc) )  THEN
1170                            flags(k,j,i) = IBSET( flags(k,j,i), 2 )
1171                         ENDIF
1172!
1173!--                      North wall
1174                         IF ( k*inc <= nzb_local((j+1)*inc,i*inc) )  THEN
1175                            flags(k,j,i) = IBSET( flags(k,j,i), 3 )
1176                         ENDIF
1177!
1178!--                      Left wall
1179                         IF ( k*inc <= nzb_local(j*inc,(i-1)*inc) )  THEN
1180                            flags(k,j,i) = IBSET( flags(k,j,i), 4 )
1181                         ENDIF
1182!
1183!--                      Right wall
1184                         IF ( k*inc <= nzb_local(j*inc,(i+1)*inc) )  THEN
1185                            flags(k,j,i) = IBSET( flags(k,j,i), 5 )
1186                         ENDIF
1187
1188                      ENDIF
1189                           
1190                   ENDDO
1191                ENDDO
1192             ENDDO
1193
1194          ENDIF
1195
1196!
1197!--       Test output of flag arrays
1198!          i = nxl_l
1199!          WRITE (9,*)  ' '
1200!          WRITE (9,*)  '*** mg level ', l, ' ***', mg_switch_to_pe0_level
1201!          WRITE (9,*)  '    inc=', inc, '  i =', nxl_l
1202!          WRITE (9,*)  '    nxl_l',nxl_l,' nxr_l=',nxr_l,' nys_l=',nys_l,' nyn_l=',nyn_l
1203!          DO  k = nzt_l+1, nzb, -1
1204!             WRITE (9,'(194(1X,I2))')  ( flags(k,j,i), j = nys_l-1, nyn_l+1 )
1205!          ENDDO
1206
1207          inc = inc * 2
1208
1209       ENDDO
1210
1211    ENDIF
1212!
1213!-- Allocate flags needed for masking walls.
1214    ALLOCATE( wall_flags_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
1215              wall_flags_00(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
1216    wall_flags_0  = 0
1217    wall_flags_00 = 0
1218
1219    IF ( scalar_advec == 'ws-scheme' .OR.                                     &
1220         scalar_advec == 'ws-scheme-mono' )  THEN
1221!
1222!--    Set flags to steer the degradation of the advection scheme in advec_ws
1223!--    near topography, inflow- and outflow boundaries as well as bottom and
1224!--    top of model domain. wall_flags_0 remains zero for all non-prognostic
1225!--    grid points.
1226       DO  i = nxl, nxr
1227          DO  j = nys, nyn
1228             DO  k = nzb_s_inner(j,i)+1, nzt
1229!
1230!--             scalar - x-direction
1231!--             WS1 (0), WS3 (1), WS5 (2)
1232                IF ( k <= nzb_s_inner(j,i+1) .OR.                              &
1233                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1234                           .AND. i == nxl   )    .OR.                          &
1235                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1236                           .AND. i == nxr   ) )                                &
1237                THEN
1238                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 0 )
1239                ELSEIF ( k <= nzb_s_inner(j,i+2) .OR. k <= nzb_s_inner(j,i-1)  &
1240                                                 .OR.                          &
1241                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1242                           .AND. i == nxr-1 )    .OR.                          &
1243                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1244                           .AND. i == nxlu  ) )                                &
1245                THEN
1246                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 1 )
1247                ELSE
1248                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 2 )
1249                ENDIF
1250!
1251!--             scalar - y-direction
1252!--             WS1 (3), WS3 (4), WS5 (5)
1253                IF ( k <= nzb_s_inner(j+1,i) .OR.                              &
1254                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1255                           .AND. j == nys   )    .OR.                          &
1256                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1257                           .AND. j == nyn   ) )                                &
1258                THEN
1259                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 3 )
1260!
1261!--             WS3
1262                ELSEIF ( k <= nzb_s_inner(j+2,i) .OR. k <= nzb_s_inner(j-1,i)  &
1263                                                 .OR.                          &
1264                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1265                           .AND. j == nysv  )    .OR.                          &
1266                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1267                           .AND. j == nyn-1 ) )                                &
1268                THEN
1269                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 4 )
1270!
1271!--             WS5
1272                ELSE
1273                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 5 )
1274                ENDIF
1275!
1276!--             scalar - z-direction
1277!--             WS1 (6), WS3 (7), WS5 (8)
1278                flag_set = .FALSE.
1279                IF ( k == nzb_s_inner(j,i) + 1 .OR. k == nzt )  THEN
1280                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 6 )
1281                   flag_set = .TRUE.
1282                ELSEIF ( k == nzb_s_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1283                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 7 )
1284                   flag_set = .TRUE.
1285                ELSEIF ( k > nzb_s_inner(j,i) .AND. .NOT. flag_set )  THEN
1286                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 8 )
1287                ENDIF
1288             ENDDO
1289          ENDDO
1290       ENDDO
1291    ENDIF
1292
1293    IF ( momentum_advec == 'ws-scheme' )  THEN
1294!
1295!--    Set wall_flags_0 to steer the degradation of the advection scheme in advec_ws
1296!--    near topography, inflow- and outflow boundaries as well as bottom and
1297!--    top of model domain. wall_flags_0 remains zero for all non-prognostic
1298!--    grid points.
1299       DO  i = nxl, nxr
1300          DO  j = nys, nyn
1301             DO  k = nzb+1, nzt
1302!
1303!--             At first, set flags to WS1.
1304!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1305!--             in order to handle the left/south flux.
1306!--             near vertical walls.
1307                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 9 )
1308                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 12 )
1309!
1310!--             u component - x-direction
1311!--             WS1 (9), WS3 (10), WS5 (11)
1312                IF ( k <= nzb_u_inner(j,i+1) .OR.                              &
1313                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1314                           .AND. i <= nxlu  )    .OR.                          &
1315                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1316                           .AND. i == nxr   ) )                                &
1317                THEN
1318                    wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 9 )
1319                ELSEIF ( k <= nzb_u_inner(j,i+2) .OR. k <= nzb_u_inner(j,i-1)  &
1320                                                 .OR.                          &
1321                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1322                           .AND. i == nxr-1 )    .OR.                          &
1323                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1324                           .AND. i == nxlu+1) )                                &
1325                THEN
1326                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 10 )
1327!
1328!--                Clear flag for WS1
1329                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 9 )
1330                ELSE
1331                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 11 )
1332!
1333!--                Clear flag for WS1
1334                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 9 )
1335                ENDIF
1336!
1337!--             u component - y-direction
1338!--             WS1 (12), WS3 (13), WS5 (14)
1339                IF ( k <= nzb_u_inner(j+1,i) .OR.                              &
1340                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1341                           .AND. j == nys   )    .OR.                          &
1342                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1343                           .AND. j == nyn   ) )                                &
1344                THEN
1345                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 12 )
1346                ELSEIF ( k <= nzb_u_inner(j+2,i) .OR. k <= nzb_u_inner(j-1,i)  &
1347                                                 .OR.                          &
1348                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1349                           .AND. j == nysv  )    .OR.                          &
1350                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1351                           .AND. j == nyn-1 ) )                                &
1352                THEN
1353                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 13 )
1354!
1355!--                Clear flag for WS1
1356                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 12 )
1357                ELSE
1358                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 14 )
1359!
1360!--                Clear flag for WS1
1361                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 12 )
1362                ENDIF
1363!
1364!--             u component - z-direction
1365!--             WS1 (15), WS3 (16), WS5 (17)
1366                flag_set = .FALSE.
1367                IF ( k == nzb_u_inner(j,i) + 1 .OR. k == nzt )  THEN
1368                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 15 )
1369                   flag_set = .TRUE.
1370                ELSEIF ( k == nzb_u_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1371                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 16 )
1372                   flag_set = .TRUE.
1373                ELSEIF ( k > nzb_u_inner(j,i) .AND. .NOT. flag_set )  THEN
1374                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 17 )
1375                ENDIF
1376
1377             ENDDO
1378          ENDDO
1379       ENDDO
1380
1381       DO  i = nxl, nxr
1382          DO  j = nys, nyn
1383             DO  k = nzb+1, nzt
1384!
1385!--             At first, set flags to WS1.
1386!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1387!--             in order to handle the left/south flux.
1388                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 18 )
1389                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 21 )
1390!
1391!--             v component - x-direction
1392!--             WS1 (18), WS3 (19), WS5 (20)
1393                IF ( k <= nzb_v_inner(j,i+1) .OR.                              &
1394                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1395                           .AND. i == nxl   )    .OR.                          &
1396                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1397                           .AND. i == nxr   ) )                                &
1398                THEN
1399                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 18 )
1400!
1401!--             WS3
1402                ELSEIF ( k <= nzb_v_inner(j,i+2) .OR. k <= nzb_v_inner(j,i-1)  &
1403                                                 .OR.                          &
1404                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1405                           .AND. i == nxr-1 )    .OR.                          &
1406                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1407                           .AND. i == nxlu  ) )                                &
1408                THEN
1409                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 19 )
1410!
1411!--                Clear flag for WS1
1412                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 18 )
1413                ELSE
1414                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 20 )
1415!
1416!--                Clear flag for WS1
1417                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 18 )
1418                ENDIF
1419!
1420!--             v component - y-direction
1421!--             WS1 (21), WS3 (22), WS5 (23)
1422                IF ( k <= nzb_v_inner(j+1,i) .OR.                              &
1423                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1424                           .AND. j <= nysv  )    .OR.                          &
1425                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1426                           .AND. j == nyn   ) )                                &
1427                THEN
1428                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 21 )
1429                ELSEIF ( k <= nzb_v_inner(j+2,i) .OR. k <= nzb_v_inner(j-1,i)  &
1430                                                 .OR.                          &
1431                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1432                           .AND. j == nysv+1)    .OR.                          &
1433                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1434                           .AND. j == nyn-1 ) )                                &
1435                THEN
1436                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 22 )
1437!
1438!--                Clear flag for WS1
1439                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 21 )
1440                ELSE
1441                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 23 )
1442!
1443!--                Clear flag for WS1
1444                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 21 )
1445                ENDIF
1446!
1447!--             v component - z-direction
1448!--             WS1 (24), WS3 (25), WS5 (26)
1449                flag_set = .FALSE.
1450                IF ( k == nzb_v_inner(j,i) + 1 .OR. k == nzt )  THEN
1451                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 24 )
1452                   flag_set = .TRUE.
1453                ELSEIF ( k == nzb_v_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1454                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 25 )
1455                   flag_set = .TRUE.
1456                ELSEIF ( k > nzb_v_inner(j,i) .AND. .NOT. flag_set )  THEN
1457                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 26 )
1458                ENDIF
1459
1460             ENDDO
1461          ENDDO
1462       ENDDO
1463       DO  i = nxl, nxr
1464          DO  j = nys, nyn
1465             DO  k = nzb+1, nzt
1466!
1467!--             At first, set flags to WS1.
1468!--             Since fluxes are swapped in advec_ws.f90, this is necessary to
1469!--             in order to handle the left/south flux.
1470                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 27 )
1471                wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 30 )
1472!
1473!--             w component - x-direction
1474!--             WS1 (27), WS3 (28), WS5 (29)
1475                IF ( k <= nzb_w_inner(j,i+1) .OR.                              &
1476                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1477                           .AND. i == nxl   )    .OR.                          &
1478                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1479                           .AND. i == nxr   ) )                                &
1480                THEN
1481                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 27 )
1482                ELSEIF ( k <= nzb_w_inner(j,i+2) .OR. k <= nzb_w_inner(j,i-1)  &
1483                                                 .OR.                          &
1484                         ( ( inflow_r .OR. outflow_r .OR. nest_bound_r )       &
1485                           .AND. i == nxr-1 )    .OR.                          &
1486                         ( ( inflow_l .OR. outflow_l .OR. nest_bound_l )       &
1487                           .AND. i == nxlu  ) )                                &
1488                THEN
1489                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 28 )
1490!
1491!--                Clear flag for WS1
1492                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 27 )
1493                ELSE
1494                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i),29 )
1495!
1496!--                Clear flag for WS1
1497                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 27 )
1498                ENDIF
1499!
1500!--             w component - y-direction
1501!--             WS1 (30), WS3 (31), WS5 (32)
1502                IF ( k <= nzb_w_inner(j+1,i) .OR.                              &
1503                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1504                           .AND. j == nys   )    .OR.                          &
1505                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1506                           .AND. j == nyn   ) )                                &
1507                THEN
1508                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 30 )
1509                ELSEIF ( k <= nzb_w_inner(j+2,i) .OR. k <= nzb_w_inner(j-1,i)  &
1510                                                 .OR.                          &
1511                         ( ( inflow_s .OR. outflow_s .OR. nest_bound_s )       &
1512                           .AND. j == nysv  )    .OR.                          &
1513                         ( ( inflow_n .OR. outflow_n .OR. nest_bound_n )       &
1514                           .AND. j == nyn-1 ) )                                &
1515                THEN
1516                   wall_flags_0(k,j,i) = IBSET( wall_flags_0(k,j,i), 31 )
1517!
1518!--                Clear flag for WS1
1519                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 30 )
1520                ELSE
1521                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 0 )
1522!
1523!--                Clear flag for WS1
1524                   wall_flags_0(k,j,i) = IBCLR( wall_flags_0(k,j,i), 30 )
1525                ENDIF
1526!
1527!--             w component - z-direction
1528!--             WS1 (33), WS3 (34), WS5 (35)
1529                flag_set = .FALSE.
1530                IF ( k == nzb_w_inner(j,i) .OR. k == nzb_w_inner(j,i) + 1      &
1531                                           .OR. k == nzt )  THEN
1532!
1533!--                Please note, at k == nzb_w_inner(j,i) a flag is explictely
1534!--                set, although this is not a prognostic level. However,
1535!--                contrary to the advection of u,v and s this is necessary
1536!--                because flux_t(nzb_w_inner(j,i)) is used for the tendency
1537!--                at k == nzb_w_inner(j,i)+1.
1538                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 1 )
1539                   flag_set = .TRUE.
1540                ELSEIF ( k == nzb_w_inner(j,i) + 2 .OR. k == nzt - 1 )  THEN
1541                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 2 )
1542                   flag_set = .TRUE.
1543                ELSEIF ( k > nzb_w_inner(j,i) .AND. .NOT. flag_set )  THEN
1544                   wall_flags_00(k,j,i) = IBSET( wall_flags_00(k,j,i), 3 )
1545                ENDIF
1546
1547             ENDDO
1548          ENDDO
1549       ENDDO
1550
1551    ENDIF
1552
1553!
1554!-- Exchange 3D integer wall_flags.
1555    IF ( momentum_advec == 'ws-scheme' .OR. scalar_advec == 'ws-scheme'     &
1556    .OR. scalar_advec == 'ws-scheme-mono' )  THEN 
1557!
1558!--    Exchange ghost points for advection flags
1559       CALL exchange_horiz_int( wall_flags_0,  nbgp )
1560       CALL exchange_horiz_int( wall_flags_00, nbgp )
1561!
1562!--    Set boundary flags at inflow and outflow boundary in case of
1563!--    non-cyclic boundary conditions.
1564       IF ( inflow_l .OR. outflow_l .OR. nest_bound_l )  THEN
1565          wall_flags_0(:,:,nxl-1)  = wall_flags_0(:,:,nxl)
1566          wall_flags_00(:,:,nxl-1) = wall_flags_00(:,:,nxl)
1567       ENDIF
1568
1569       IF ( inflow_r .OR. outflow_r .OR. nest_bound_r )  THEN
1570          wall_flags_0(:,:,nxr+1)  = wall_flags_0(:,:,nxr)
1571          wall_flags_00(:,:,nxr+1) = wall_flags_00(:,:,nxr)
1572       ENDIF
1573
1574       IF ( inflow_n .OR. outflow_n .OR. nest_bound_n )  THEN
1575          wall_flags_0(:,nyn+1,:)  = wall_flags_0(:,nyn,:)
1576          wall_flags_00(:,nyn+1,:) = wall_flags_00(:,nyn,:)
1577       ENDIF
1578
1579       IF ( inflow_s .OR. outflow_s  .OR. nest_bound_s )  THEN
1580          wall_flags_0(:,nys-1,:)  = wall_flags_0(:,nys,:)
1581          wall_flags_00(:,nys-1,:) = wall_flags_00(:,nys,:)
1582       ENDIF
1583
1584    ENDIF
1585
1586!
1587!-- In case of topography: limit near-wall mixing length l_wall further:
1588!-- Go through all points of the subdomain one by one and look for the closest
1589!-- surface
1590    IF ( TRIM(topography) /= 'flat' )  THEN
1591       DO  i = nxl, nxr
1592          DO  j = nys, nyn
1593
1594             nzb_si = nzb_s_inner(j,i)
1595             vi     = vertical_influence(nzb_si)
1596
1597             IF ( wall_n(j,i) > 0 )  THEN
1598!
1599!--             North wall (y distance)
1600                DO  k = wall_n(j,i), nzb_si
1601                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5_wp * dy )
1602                ENDDO
1603!
1604!--             Above North wall (yz distance)
1605                DO  k = nzb_si + 1, nzb_si + vi
1606                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),                     &
1607                                          SQRT( 0.25_wp * dy**2 +              &
1608                                          ( zu(k) - zw(nzb_si) )**2 ) )
1609                ENDDO
1610!
1611!--             Northleft corner (xy distance)
1612                IF ( corner_nl(j,i) > 0 )  THEN
1613                   DO  k = corner_nl(j,i), nzb_si
1614                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
1615                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1616                   ENDDO
1617!
1618!--                Above Northleft corner (xyz distance)
1619                   DO  k = nzb_si + 1, nzb_si + vi
1620                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),              &
1621                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1622                                            ( zu(k) - zw(nzb_si) )**2 ) )
1623                   ENDDO
1624                ENDIF
1625!
1626!--             Northright corner (xy distance)
1627                IF ( corner_nr(j,i) > 0 )  THEN
1628                   DO  k = corner_nr(j,i), nzb_si
1629                       l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),             &
1630                                                0.5_wp * SQRT( dx**2 + dy**2 ) )
1631                   ENDDO
1632!
1633!--                Above northright corner (xyz distance)
1634                   DO  k = nzb_si + 1, nzb_si + vi
1635                      l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),              &
1636                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1637                                            ( zu(k) - zw(nzb_si) )**2 ) )
1638                   ENDDO
1639                ENDIF
1640             ENDIF
1641
1642             IF ( wall_s(j,i) > 0 )  THEN
1643!
1644!--             South wall (y distance)
1645                DO  k = wall_s(j,i), nzb_si
1646                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5_wp * dy )
1647                ENDDO
1648!
1649!--             Above south wall (yz distance)
1650                DO  k = nzb_si + 1, nzb_si + vi
1651                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),                     &
1652                                          SQRT( 0.25_wp * dy**2 +              &
1653                                          ( zu(k) - zw(nzb_si) )**2 ) )
1654                ENDDO
1655!
1656!--             Southleft corner (xy distance)
1657                IF ( corner_sl(j,i) > 0 )  THEN
1658                   DO  k = corner_sl(j,i), nzb_si
1659                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1660                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1661                   ENDDO
1662!
1663!--                Above southleft corner (xyz distance)
1664                   DO  k = nzb_si + 1, nzb_si + vi
1665                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1666                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1667                                            ( zu(k) - zw(nzb_si) )**2 ) )
1668                   ENDDO
1669                ENDIF
1670!
1671!--             Southright corner (xy distance)
1672                IF ( corner_sr(j,i) > 0 )  THEN
1673                   DO  k = corner_sr(j,i), nzb_si
1674                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1675                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1676                   ENDDO
1677!
1678!--                Above southright corner (xyz distance)
1679                   DO  k = nzb_si + 1, nzb_si + vi
1680                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1681                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1682                                            ( zu(k) - zw(nzb_si) )**2 ) )
1683                   ENDDO
1684                ENDIF
1685
1686             ENDIF
1687
1688             IF ( wall_l(j,i) > 0 )  THEN
1689!
1690!--             Left wall (x distance)
1691                DO  k = wall_l(j,i), nzb_si
1692                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5_wp * dx )
1693                ENDDO
1694!
1695!--             Above left wall (xz distance)
1696                DO  k = nzb_si + 1, nzb_si + vi
1697                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),                     &
1698                                       SQRT( 0.25_wp * dx**2 +                 &
1699                                       ( zu(k) - zw(nzb_si) )**2 ) )
1700                ENDDO
1701             ENDIF
1702
1703             IF ( wall_r(j,i) > 0 )  THEN
1704!
1705!--             Right wall (x distance)
1706                DO  k = wall_r(j,i), nzb_si
1707                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5_wp * dx )
1708                ENDDO
1709!
1710!--             Above right wall (xz distance)
1711                DO  k = nzb_si + 1, nzb_si + vi
1712                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),                     &
1713                                          SQRT( 0.25_wp * dx**2 +              &
1714                                          ( zu(k) - zw(nzb_si) )**2 ) )
1715                ENDDO
1716
1717             ENDIF
1718
1719          ENDDO
1720       ENDDO
1721
1722    ENDIF
1723
1724!
1725!-- Multiplication with wall_adjustment_factor
1726    l_wall = wall_adjustment_factor * l_wall
1727
1728!
1729!-- Set lateral boundary conditions for l_wall
1730    CALL exchange_horiz( l_wall, nbgp )
1731
1732    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
1733                nzb_tmp, vertical_influence, wall_l, wall_n, wall_r, wall_s )
1734
1735#endif
1736
1737 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.