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

Last change on this file since 2412 was 2365, checked in by kanani, 7 years ago

Vertical nesting implemented (SadiqHuq?)

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