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

Last change on this file since 4701 was 4691, checked in by suehring, 4 years ago

Reference model topography to the lowest grid point also in ASCII input case

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