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

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

Bugfix in 3D building initialization

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