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

Last change on this file since 4327 was 4314, checked in by suehring, 5 years ago

plant-canopy: Bugfix, plant canopy was still considered at building edges on for the u- and v-component; Relax restriction of LAD on building tops. LAD is only omitted at locations where building grid points emerged artificially by the topography filtering

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