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

Last change on this file since 2290 was 2274, checked in by Giersch, 7 years ago

Changed error messages

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