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

Last change on this file since 2604 was 2550, checked in by boeske, 7 years ago

enable simulations with complex terrain

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