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

Last change on this file since 4407 was 4386, checked in by Giersch, 5 years ago

Allocation statements, comments, naming of variables revised and _wp added to real type values

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