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

Last change on this file since 1910 was 1910, checked in by raasch, 8 years ago

Bugfix: if topography is read from file, Neumann conditions are used for the nzb_local array (instead of cyclic conditions) in case that non-cyclic boundary conditions are used for the run

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