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

Last change on this file since 1354 was 1354, checked in by heinze, 10 years ago

last commit documented

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