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

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

Relax restrictions for topography input via static input file, terrain and building heights, as well as building IDs can be input separately and are not mandatory any more.

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