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

Last change on this file since 2343 was 2319, checked in by suehring, 7 years ago

remove print statements

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