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

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

Bugfixes concerning top fluxes and TKE production

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