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

Last change on this file since 3083 was 3068, checked in by Giersch, 6 years ago

Warning message related to vertical grid stretching has been added, old example_cbl_rc file restored but with the new header, example_cbl_p3d changed to get the old run control output

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