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

Last change on this file since 2241 was 2233, checked in by suehring, 8 years ago

last commit documented

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