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

Last change on this file since 3206 was 3200, checked in by suehring, 6 years ago

Missing pre-processor directive

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