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

Last change on this file since 4873 was 4868, checked in by raasch, 4 years ago

Height of level k=0 for the u,v-grid is always set 0.0

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