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

Last change on this file since 2232 was 2232, checked in by suehring, 7 years ago

Adjustments according new topography and surface-modelling concept implemented

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