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

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