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

Last change on this file since 2918 was 2918, checked in by gronemeier, 6 years ago

moved initialization of mixing length to turbulence_closure_mod; consider walls in calculation of ml in rans-mode

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