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

Last change on this file since 1320 was 1320, checked in by raasch, 10 years ago

ONLY-attribute added to USE-statements,
kind-parameters added to all INTEGER and REAL declaration statements,
kinds are defined in new module kinds,
old module precision_kind is removed,
revision history before 2012 removed,
comment fields (!:) to be used for variable explanations added to all variable declaration statements

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