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

Last change on this file since 2696 was 2696, checked in by kanani, 6 years ago

Merge of branch palm4u into trunk

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