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

Last change on this file since 4543 was 4543, checked in by gronemeier, 4 years ago

Remove non-required check for canyon height (init_grid)

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