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

Last change on this file since 2723 was 2718, checked in by maronga, 7 years ago

deleting of deprecated files; headers updated where needed

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