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

Last change on this file since 3051 was 3051, checked in by suehring, 6 years ago

Speed-up NetCDF input; Revise NetCDF-input routines and remove input via io-blocks; Temporarily revoke renaming of input variables in dynamic driver; More detailed error messages created; Bugfix in mapping 3D buildings; Bugfix in land-surface model at pavement surfaces; Bugfix in initialization with inifor

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