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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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