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

Last change on this file since 2953 was 2927, checked in by suehring, 7 years ago

Bugfix in setting boundary condition for topography index array in case of non-cyclic boundary conditions and Piascek-Williams advection scheme is used.

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