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

Last change on this file since 4835 was 4828, checked in by Giersch, 4 years ago

Copyright updated to year 2021, interface pmc_sort removed to accelarate the nesting code

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