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

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

Code adjusted according to coding standards, renamed namelists, error messages revised until PA0347, output CASE 108 disabled

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