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

Last change on this file since 2774 was 2747, checked in by suehring, 7 years ago

Bugfix, topography is rounded to the nearest discrete grid level

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