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

Last change on this file since 3026 was 2968, checked in by suehring, 6 years ago

Bugfixes in initalization of land-surface model as well as timeseries data output in case of elevated model surfaces.

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