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

Last change on this file since 4230 was 4189, checked in by suehring, 5 years ago

init_grid: Add check for proper setting of namelist parameter topography; Set bit flag to indicate land surfaces in case no topography is provided

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