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

Last change on this file since 4598 was 4564, checked in by raasch, 4 years ago

Vertical nesting method of Huq et al. (2019) removed

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