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

Last change on this file since 1978 was 1969, checked in by suehring, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 61.7 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 terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2016 Leibniz Universitaet Hannover
17!-------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: init_grid.f90 1969 2016-07-18 12:03:12Z maronga $
26!
27! 1968 2016-07-18 12:01:49Z suehring
28! Changed: PE-wise reading of topography file in order to avoid global definition
29! of arrays nzb_local and nzb_tmp. Thereby, topography definition for single
30! buildings and street canyons has changed, as well as flag setting for
31! multigrid scheme.
32!
33! Bugfix in checking l_grid anisotropy.
34! Simplify initial computation of lwall and vertical_influence, i.e. remove
35! nzb_s_inner as it is still zero at this point.
36!
37! 1942 2016-06-14 12:18:18Z suehring
38! Topography filter implemented to fill holes resolved by only one grid point.
39! Initialization of flags for ws-scheme moved to advec_ws. 
40!
41! 1931 2016-06-10 12:06:59Z suehring
42! Rename multigrid into multigrid_noopt and multigrid_fast into multigrid
43!
44! 1910 2016-05-26 06:49:46Z raasch
45! Bugfix: if topography is read from file, Neumann conditions are used for the
46! nzb_local array (instead of cyclic conditions) in case that non-cyclic
47! boundary conditions are switched on for the run
48!
49! 1902 2016-05-09 11:18:56Z suehring
50! Set topography flags for multigrid solver only (not for multigrid_fast)
51!
52! 1886 2016-04-21 11:20:47Z suehring
53! Bugfix: setting advection flags near walls
54! reformulated index values for nzb_v_inner
55! variable discriptions added in declaration block
56!
57! 1845 2016-04-08 08:29:13Z raasch
58! nzb_2d removed
59!
60! 1804 2016-04-05 16:30:18Z maronga
61! Removed code for parameter file check (__check)
62!
63! 1779 2016-03-03 08:01:28Z raasch
64! coupling_char is trimmed at every place it occurs, because it can have
65! different length now
66!
67! 1762 2016-02-25 12:31:13Z hellstea
68! Introduction of nested domain feature
69!
70! 1743 2016-01-13 10:23:51Z raasch
71! Bugfix for calculation of nzb_s_outer and nzb_u_outer at north boundary of
72! total domain
73!
74! 1691 2015-10-26 16:17:44Z maronga
75! Renamed prandtl_layer to constant_flux_layer.
76!
77! 1682 2015-10-07 23:56:08Z knoop
78! Code annotations made doxygen readable
79!
80! 1677 2015-10-02 13:25:23Z boeske
81! Bugfix: Ghost points are included in wall_flags_0 and wall_flags_00
82!
83! 1675 2015-10-02 08:28:59Z gronemeier
84! Bugfix: Definition of topography grid levels
85!
86! 1660 2015-09-21 08:15:16Z gronemeier
87! Bugfix: Definition of topography grid levels if vertical grid stretching
88!         starts below the maximum topography height.
89!
90! 1580 2015-04-10 13:43:49Z suehring
91! Bugfix: setting flags for 5th order scheme near buildings
92!
93! 1575 2015-03-27 09:56:27Z raasch
94! adjustments for psolver-queries
95!
96! 1557 2015-03-05 16:43:04Z suehring
97! Adjustment for monotoinic limiter
98!
99! 1418 2014-06-06 13:05:08Z fricke
100! Bugfix: Change if-condition for stretched grid in the ocean, with the old
101!          condition and a negative value for dz_stretch_level the condition
102!          was always true for the whole model domain
103!
104! 1409 2014-05-23 12:11:32Z suehring
105! Bugfix: set wall_flags_0 at inflow and outflow boundary also for i <= nxlu
106! j <= nysv
107!
108! 1353 2014-04-08 15:21:23Z heinze
109! REAL constants provided with KIND-attribute
110!
111! 1322 2014-03-20 16:38:49Z raasch
112! REAL constants defined as wp-kind
113!
114! 1320 2014-03-20 08:40:49Z raasch
115! ONLY-attribute added to USE-statements,
116! kind-parameters added to all INTEGER and REAL declaration statements,
117! kinds are defined in new module kinds,
118! revision history before 2012 removed,
119! comment fields (!:) to be used for variable explanations added to
120! all variable declaration statements
121!
122! 1221 2013-09-10 08:59:13Z raasch
123! wall_flags_00 introduced to hold bits 32-63,
124! additional 3D-flag arrays for replacing the 2D-index array nzb_s_inner in
125! loops optimized for openACC (pres + flow_statistics)
126!
127! 1092 2013-02-02 11:24:22Z raasch
128! unused variables removed
129!
130! 1069 2012-11-28 16:18:43Z maronga
131! bugfix: added coupling_char to TOPOGRAPHY_DATA to allow topography in the
132!         ocean model in case of coupled runs
133!
134! 1036 2012-10-22 13:43:42Z raasch
135! code put under GPL (PALM 3.9)
136!
137! 1015 2012-09-27 09:23:24Z raasch
138! lower index for calculating wall_flags_0 set to nzb_w_inner instead of
139! nzb_w_inner+1
140!
141! 996 2012-09-07 10:41:47Z raasch
142! little reformatting
143!
144! 978 2012-08-09 08:28:32Z fricke
145! Bugfix: nzb_max is set to nzt at non-cyclic lateral boundaries
146! Bugfix: Set wall_flags_0 for inflow boundary
147!
148! 927 2012-06-06 19:15:04Z raasch
149! Wall flags are not set for multigrid method in case of masking method
150!
151! 864 2012-03-27 15:10:33Z gryschka
152! In case of ocean and Dirichlet bottom bc for u and v dzu_mg and ddzu_pres
153! were not correctly defined for k=1.
154!
155! 861 2012-03-26 14:18:34Z suehring
156! Set wall_flags_0. The array is needed for degradation in ws-scheme near walls,
157! inflow and outflow boundaries as well as near the bottom and the top of the
158! model domain.!
159! Initialization of nzb_s_inner and nzb_w_inner.
160! gls has to be at least nbgp to do not exceed the array bounds of nzb_local
161! while setting wall_flags_0
162!
163! 843 2012-02-29 15:16:21Z gryschka
164! In case of ocean and dirichlet bc for u and v at the bottom
165! the first u-level ist defined at same height as the first w-level
166!
167! 818 2012-02-08 16:11:23Z maronga
168! Bugfix: topo_height is only required if topography is used. It is thus now
169! allocated in the topography branch
170!
171! 809 2012-01-30 13:32:58Z maronga
172! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
173!
174! 807 2012-01-25 11:53:51Z maronga
175! New cpp directive "__check" implemented which is used by check_namelist_files
176!
177! Revision 1.1  1997/08/11 06:17:45  raasch
178! Initial revision (Testversion)
179!
180!
181! Description:
182! ------------
183!> Creating grid depending constants
184!------------------------------------------------------------------------------!
185 SUBROUTINE init_grid
186 
187    USE advec_ws,                                                              &
188        ONLY:  ws_init_flags
189
190    USE arrays_3d,                                                             &
191        ONLY:  dd2zu, ddzu, ddzu_pres, ddzw, dzu, dzu_mg, dzw, dzw_mg, f1_mg,  &
192               f2_mg, f3_mg, l_grid, l_wall, zu, zw
193       
194    USE control_parameters,                                                    &
195        ONLY:  bc_lr_cyc, bc_ns_cyc, building_height, building_length_x,       &
196               building_length_y, building_wall_left, building_wall_south,     &
197               canyon_height, canyon_wall_left, canyon_wall_south,             &
198               canyon_width_x, canyon_width_y, constant_flux_layer,            &
199               coupling_char, dp_level_ind_b, dz, dz_max, dz_stretch_factor,   &
200               dz_stretch_level, dz_stretch_level_index, grid_level, ibc_uv_b, &
201               io_blocks, io_group, inflow_l, inflow_n, inflow_r, inflow_s,    &
202               masking_method, maximum_grid_level, message_string,             &
203               momentum_advec, nest_domain, ocean, outflow_l, outflow_n,       &
204               outflow_r, outflow_s, psolver, scalar_advec, topography,        &
205               topography_grid_convention, use_surface_fluxes, use_top_fluxes, &
206               wall_adjustment_factor 
207       
208    USE grid_variables,                                                        &
209        ONLY:  ddx, ddx2, ddx2_mg, ddy, ddy2, ddy2_mg, dx, dx2, dy, dy2, fwxm, &
210               fwxp, fwym, fwyp, fxm, fxp, fym, fyp, wall_e_x, wall_e_y,       &
211               wall_u, wall_v, wall_w_x, wall_w_y, zu_s_inner, zw_w_inner
212       
213    USE indices,                                                               &
214        ONLY:  flags, nbgp, nx, nxl, nxlg, nxl_mg, nxr, nxrg, nxr_mg,          &
215               ny, nyn, nyng, nyn_mg, nys, nys_mg, nysg, nz, nzb,              &
216               nzb_diff, nzb_diff_s_inner, nzb_diff_s_outer, nzb_diff_u,       &
217               nzb_diff_v, nzb_max, nzb_s_inner, nzb_s_outer, nzb_u_inner,     &
218               nzb_u_outer, nzb_v_inner, nzb_v_outer, nzb_w_inner,             &
219               nzb_w_outer, nzt, nzt_diff, nzt_mg, rflags_invers,              &
220               rflags_s_inner, wall_flags_0, wall_flags_00, wall_flags_1,      &
221               wall_flags_10, wall_flags_2, wall_flags_3,  wall_flags_4,       &
222               wall_flags_5, wall_flags_6, wall_flags_7, wall_flags_8,         &
223               wall_flags_9
224   
225    USE kinds
226   
227    USE pegrid
228
229    IMPLICIT NONE
230
231    INTEGER(iwp) ::  bh            !< temporary vertical index of building height
232    INTEGER(iwp) ::  blx           !< grid point number of building size along x
233    INTEGER(iwp) ::  bly           !< grid point number of building size along y
234    INTEGER(iwp) ::  bxl           !< index for left building wall
235    INTEGER(iwp) ::  bxr           !< index for right building wall
236    INTEGER(iwp) ::  byn           !< index for north building wall
237    INTEGER(iwp) ::  bys           !< index for south building wall
238    INTEGER(iwp) ::  ch            !< temporary vertical index for canyon height
239    INTEGER(iwp) ::  cwx           !< grid point number of canyon size along x
240    INTEGER(iwp) ::  cwy           !< grid point number of canyon size along y
241    INTEGER(iwp) ::  cxl           !< index for left canyon wall
242    INTEGER(iwp) ::  cxr           !< index for right canyon wall
243    INTEGER(iwp) ::  cyn           !< index for north canyon wall
244    INTEGER(iwp) ::  cys           !< index for south canyon wall
245    INTEGER(iwp) ::  i             !< index variable along x
246    INTEGER(iwp) ::  ii            !< loop variable for reading topography file
247    INTEGER(iwp) ::  inc           !< incremental parameter for coarsening grid level
248    INTEGER(iwp) ::  j             !< index variable along y
249    INTEGER(iwp) ::  k             !< index variable along z
250    INTEGER(iwp) ::  l             !< loop variable
251    INTEGER(iwp) ::  nxl_l         !< index of left PE boundary for multigrid level
252    INTEGER(iwp) ::  nxr_l         !< index of right PE boundary for multigrid level
253    INTEGER(iwp) ::  nyn_l         !< index of north PE boundary for multigrid level
254    INTEGER(iwp) ::  nys_l         !< index of south PE boundary for multigrid level
255    INTEGER(iwp) ::  nzb_local_max !< vertical grid index of maximum topography height
256    INTEGER(iwp) ::  nzb_local_min !< vertical grid index of minimum topography height
257    INTEGER(iwp) ::  nzb_si        !< dummy index for local nzb_s_inner
258    INTEGER(iwp) ::  nzt_l         !< index of top PE boundary for multigrid level
259    INTEGER(iwp) ::  num_hole      !< number of holes (in topography) resolved by only one grid point
260    INTEGER(iwp) ::  num_hole_l    !< number of holes (in topography) resolved by only one grid point on local PE     
261    INTEGER(iwp) ::  num_wall      !< number of surrounding vertical walls for a single grid point
262    INTEGER(iwp) ::  skip_n_rows   !< counting variable to skip rows while reading topography file   
263    INTEGER(iwp) ::  vi            !< dummy for vertical influence
264
265    INTEGER(iwp), DIMENSION(:), ALLOCATABLE   ::                               &
266                     vertical_influence  !< number of vertical grid points above obstacle where adjustment of near-wall mixing length is required
267                                         
268    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nl      !< index of north-left corner location to limit near-wall mixing length
269    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_nr      !< north-right
270    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sl      !< south-left
271    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  corner_sr      !< south-right
272    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_local      !< index for topography top at cell-center
273    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_tmp        !< dummy to calculate topography indices on u- and v-grid
274    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_l         !< distance to adjacent left-facing wall
275    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_n         !< north-facing
276    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_r         !< right-facing
277    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  wall_s         !< right-facing
278
279    REAL(wp) ::  dum           !< dummy variable to skip columns while reading topography file   
280    REAL(wp) ::  dx_l          !< grid spacing along x on different multigrid level
281    REAL(wp) ::  dy_l          !< grid spacing along y on different multigrid level
282    REAL(wp) ::  dz_stretched  !< stretched vertical grid spacing
283
284    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  topo_height   !< input variable for topography height
285    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  zu_s_inner_l  !< dummy array on global scale to write topography output array
286    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  zw_w_inner_l  !< dummy array on global scale to write topography output array
287
288   
289!
290!-- Calculation of horizontal array bounds including ghost layers
291    nxlg = nxl - nbgp
292    nxrg = nxr + nbgp
293    nysg = nys - nbgp
294    nyng = nyn + nbgp
295
296!
297!-- Allocate grid arrays
298    ALLOCATE( ddzu(1:nzt+1), ddzw(1:nzt+1), dd2zu(1:nzt), dzu(1:nzt+1),        &
299              dzw(1:nzt+1), l_grid(1:nzt), zu(nzb:nzt+1), zw(nzb:nzt+1) )
300
301!
302!-- Compute height of u-levels from constant grid length and dz stretch factors
303    IF ( dz == -1.0_wp )  THEN
304       message_string = 'missing dz'
305       CALL message( 'init_grid', 'PA0200', 1, 2, 0, 6, 0 ) 
306    ELSEIF ( dz <= 0.0_wp )  THEN
307       WRITE( message_string, * ) 'dz=',dz,' <= 0.0'
308       CALL message( 'init_grid', 'PA0201', 1, 2, 0, 6, 0 )
309    ENDIF
310
311!
312!-- Define the vertical grid levels
313    IF ( .NOT. ocean )  THEN
314!
315!--    Grid for atmosphere with surface at z=0 (k=0, w-grid).
316!--    The second u-level (k=1) corresponds to the top of the
317!--    Prandtl-layer.
318
319       IF ( ibc_uv_b == 0 .OR. ibc_uv_b == 2 ) THEN
320          zu(0) = 0.0_wp
321      !    zu(0) = - dz * 0.5_wp
322       ELSE
323          zu(0) = - dz * 0.5_wp
324       ENDIF
325       zu(1) =   dz * 0.5_wp
326
327       dz_stretch_level_index = nzt+1
328       dz_stretched = dz
329       DO  k = 2, nzt+1
330          IF ( dz_stretch_level <= zu(k-1)  .AND.  dz_stretched < dz_max )  THEN
331             dz_stretched = dz_stretched * dz_stretch_factor
332             dz_stretched = MIN( dz_stretched, dz_max )
333             IF ( dz_stretch_level_index == nzt+1 ) dz_stretch_level_index = k-1
334          ENDIF
335          zu(k) = zu(k-1) + dz_stretched
336       ENDDO
337
338!
339!--    Compute the w-levels. They are always staggered half-way between the
340!--    corresponding u-levels. In case of dirichlet bc for u and v at the
341!--    ground the first u- and w-level (k=0) are defined at same height (z=0).
342!--    The top w-level is extrapolated linearly.
343       zw(0) = 0.0_wp
344       DO  k = 1, nzt
345          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
346       ENDDO
347       zw(nzt+1) = zw(nzt) + 2.0_wp * ( zu(nzt+1) - zw(nzt) )
348
349    ELSE
350!
351!--    Grid for ocean with free water surface is at k=nzt (w-grid).
352!--    In case of neumann bc at the ground the first first u-level (k=0) lies
353!--    below the first w-level (k=0). In case of dirichlet bc the first u- and
354!--    w-level are defined at same height, but staggered from the second level.
355!--    The second u-level (k=1) corresponds to the top of the Prandtl-layer.
356       zu(nzt+1) =   dz * 0.5_wp
357       zu(nzt)   = - dz * 0.5_wp
358
359       dz_stretch_level_index = 0
360       dz_stretched = dz
361       DO  k = nzt-1, 0, -1
362!
363!--       The default value of dz_stretch_level is positive, thus the first
364!--       condition is always true. Hence, the second condition is necessary.
365          IF ( dz_stretch_level >= zu(k+1)  .AND.  dz_stretch_level <= 0.0  &
366               .AND.  dz_stretched < dz_max )  THEN
367             dz_stretched = dz_stretched * dz_stretch_factor
368             dz_stretched = MIN( dz_stretched, dz_max )
369             IF ( dz_stretch_level_index == 0 ) dz_stretch_level_index = k+1
370          ENDIF
371          zu(k) = zu(k+1) - dz_stretched
372       ENDDO
373
374!
375!--    Compute the w-levels. They are always staggered half-way between the
376!--    corresponding u-levels, except in case of dirichlet bc for u and v
377!--    at the ground. In this case the first u- and w-level are defined at
378!--    same height. The top w-level (nzt+1) is not used but set for
379!--    consistency, since w and all scalar variables are defined up tp nzt+1.
380       zw(nzt+1) = dz
381       zw(nzt)   = 0.0_wp
382       DO  k = 0, nzt
383          zw(k) = ( zu(k) + zu(k+1) ) * 0.5_wp
384       ENDDO
385
386!
387!--    In case of dirichlet bc for u and v the first u- and w-level are defined
388!--    at same height.
389       IF ( ibc_uv_b == 0 ) THEN
390          zu(0) = zw(0)
391       ENDIF
392
393    ENDIF
394
395!
396!-- Compute grid lengths.
397    DO  k = 1, nzt+1
398       dzu(k)  = zu(k) - zu(k-1)
399       ddzu(k) = 1.0_wp / dzu(k)
400       dzw(k)  = zw(k) - zw(k-1)
401       ddzw(k) = 1.0_wp / dzw(k)
402    ENDDO
403
404    DO  k = 1, nzt
405       dd2zu(k) = 1.0_wp / ( dzu(k) + dzu(k+1) )
406    ENDDO
407   
408!   
409!-- The FFT- SOR-pressure solvers assume grid spacings of a staggered grid
410!-- everywhere. For the actual grid, the grid spacing at the lowest level
411!-- is only dz/2, but should be dz. Therefore, an additional array
412!-- containing with appropriate grid information is created for these
413!-- solvers.
414    IF ( psolver(1:9) /= 'multigrid' )  THEN
415       ALLOCATE( ddzu_pres(1:nzt+1) )
416       ddzu_pres = ddzu
417       ddzu_pres(1) = ddzu_pres(2)  ! change for lowest level
418    ENDIF   
419
420!
421!-- In case of multigrid method, compute grid lengths and grid factors for the
422!-- grid levels
423    IF ( psolver(1:9) == 'multigrid' )  THEN
424
425       ALLOCATE( ddx2_mg(maximum_grid_level), ddy2_mg(maximum_grid_level), &
426                 dzu_mg(nzb+1:nzt+1,maximum_grid_level),                   &
427                 dzw_mg(nzb+1:nzt+1,maximum_grid_level),                   &
428                 f1_mg(nzb+1:nzt,maximum_grid_level),                      &
429                 f2_mg(nzb+1:nzt,maximum_grid_level),                      &
430                 f3_mg(nzb+1:nzt,maximum_grid_level) )
431
432       dzu_mg(:,maximum_grid_level) = dzu
433!       
434!--    Next line to ensure an equally spaced grid.
435       dzu_mg(1,maximum_grid_level) = dzu(2)
436
437       dzw_mg(:,maximum_grid_level) = dzw
438       nzt_l = nzt
439       DO  l = maximum_grid_level-1, 1, -1
440           dzu_mg(nzb+1,l) = 2.0_wp * dzu_mg(nzb+1,l+1)
441           dzw_mg(nzb+1,l) = 2.0_wp * dzw_mg(nzb+1,l+1)
442           nzt_l = nzt_l / 2
443           DO  k = 2, nzt_l+1
444              dzu_mg(k,l) = dzu_mg(2*k-2,l+1) + dzu_mg(2*k-1,l+1)
445              dzw_mg(k,l) = dzw_mg(2*k-2,l+1) + dzw_mg(2*k-1,l+1)
446           ENDDO
447       ENDDO
448
449       nzt_l = nzt
450       dx_l  = dx
451       dy_l  = dy
452       DO  l = maximum_grid_level, 1, -1
453          ddx2_mg(l) = 1.0_wp / dx_l**2
454          ddy2_mg(l) = 1.0_wp / dy_l**2
455          DO  k = nzb+1, nzt_l
456             f2_mg(k,l) = 1.0_wp / ( dzu_mg(k+1,l) * dzw_mg(k,l) )
457             f3_mg(k,l) = 1.0_wp / ( dzu_mg(k,l)   * dzw_mg(k,l) )
458             f1_mg(k,l) = 2.0_wp * ( ddx2_mg(l) + ddy2_mg(l) ) + &
459                          f2_mg(k,l) + f3_mg(k,l)
460          ENDDO
461          nzt_l = nzt_l / 2
462          dx_l  = dx_l * 2.0_wp
463          dy_l  = dy_l * 2.0_wp
464       ENDDO
465
466    ENDIF
467
468!
469!-- Compute the reciprocal values of the horizontal grid lengths.
470    ddx = 1.0_wp / dx
471    ddy = 1.0_wp / dy
472    dx2 = dx * dx
473    dy2 = dy * dy
474    ddx2 = 1.0_wp / dx2
475    ddy2 = 1.0_wp / dy2
476
477!
478!-- Compute the grid-dependent mixing length.
479    DO  k = 1, nzt
480       l_grid(k)  = ( dx * dy * dzw(k) )**0.33333333333333_wp
481    ENDDO
482
483!
484!-- Allocate outer and inner index arrays for topography and set
485!-- defaults.
486
487    ALLOCATE( corner_nl(nys:nyn,nxl:nxr), corner_nr(nys:nyn,nxl:nxr),       &
488              corner_sl(nys:nyn,nxl:nxr), corner_sr(nys:nyn,nxl:nxr),       &
489              wall_l(nys:nyn,nxl:nxr), wall_n(nys:nyn,nxl:nxr),             &
490              wall_r(nys:nyn,nxl:nxr), wall_s(nys:nyn,nxl:nxr) )                     
491     
492    ALLOCATE( fwxm(nysg:nyng,nxlg:nxrg), fwxp(nysg:nyng,nxlg:nxrg),         &
493              fwym(nysg:nyng,nxlg:nxrg), fwyp(nysg:nyng,nxlg:nxrg),         &
494              fxm(nysg:nyng,nxlg:nxrg), fxp(nysg:nyng,nxlg:nxrg),           &
495              fym(nysg:nyng,nxlg:nxrg), fyp(nysg:nyng,nxlg:nxrg),           &
496              nzb_s_inner(nysg:nyng,nxlg:nxrg),                             &
497              nzb_s_outer(nysg:nyng,nxlg:nxrg),                             &
498              nzb_u_inner(nysg:nyng,nxlg:nxrg),                             &
499              nzb_u_outer(nysg:nyng,nxlg:nxrg),                             &
500              nzb_v_inner(nysg:nyng,nxlg:nxrg),                             &
501              nzb_v_outer(nysg:nyng,nxlg:nxrg),                             &
502              nzb_w_inner(nysg:nyng,nxlg:nxrg),                             &
503              nzb_w_outer(nysg:nyng,nxlg:nxrg),                             &
504              nzb_diff_s_inner(nysg:nyng,nxlg:nxrg),                        &
505              nzb_diff_s_outer(nysg:nyng,nxlg:nxrg),                        &
506              nzb_diff_u(nysg:nyng,nxlg:nxrg),                              &
507              nzb_diff_v(nysg:nyng,nxlg:nxrg),                              &
508              nzb_local(nysg:nyng,nxlg:nxrg),                               &
509              nzb_tmp(nysg:nyng,nxlg:nxrg),                                 &
510              rflags_s_inner(nzb:nzt+2,nysg:nyng,nxlg:nxrg),                &
511              rflags_invers(nysg:nyng,nxlg:nxrg,nzb:nzt+2),                 &
512              wall_e_x(nysg:nyng,nxlg:nxrg),                                &
513              wall_e_y(nysg:nyng,nxlg:nxrg),                                &
514              wall_u(nysg:nyng,nxlg:nxrg),                                  &
515              wall_v(nysg:nyng,nxlg:nxrg),                                  &
516              wall_w_x(nysg:nyng,nxlg:nxrg),                                &
517              wall_w_y(nysg:nyng,nxlg:nxrg) )
518
519
520
521    ALLOCATE( l_wall(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
522
523
524    nzb_s_inner = nzb;  nzb_s_outer = nzb
525    nzb_u_inner = nzb;  nzb_u_outer = nzb
526    nzb_v_inner = nzb;  nzb_v_outer = nzb
527    nzb_w_inner = nzb;  nzb_w_outer = nzb
528
529    rflags_s_inner = 1.0_wp
530    rflags_invers  = 1.0_wp
531
532!
533!-- Define vertical gridpoint from (or to) which on the usual finite difference
534!-- form (which does not use surface fluxes) is applied
535    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
536       nzb_diff = nzb + 2
537    ELSE
538       nzb_diff = nzb + 1
539    ENDIF
540    IF ( use_top_fluxes )  THEN
541       nzt_diff = nzt - 1
542    ELSE
543       nzt_diff = nzt
544    ENDIF
545
546    nzb_diff_s_inner = nzb_diff;  nzb_diff_s_outer = nzb_diff
547    nzb_diff_u = nzb_diff;  nzb_diff_v = nzb_diff
548
549    wall_e_x = 0.0_wp;  wall_e_y = 0.0_wp;  wall_u = 0.0_wp;  wall_v = 0.0_wp
550    wall_w_x = 0.0_wp;  wall_w_y = 0.0_wp
551    fwxp = 1.0_wp;  fwxm = 1.0_wp;  fwyp = 1.0_wp;  fwym = 1.0_wp
552    fxp  = 1.0_wp;  fxm  = 1.0_wp;  fyp  = 1.0_wp;  fym  = 1.0_wp
553
554!
555!-- Initialize near-wall mixing length l_wall only in the vertical direction
556!-- for the moment,
557!-- multiplication with wall_adjustment_factor near the end of this routine
558    l_wall(nzb,:,:)   = l_grid(1)
559    DO  k = nzb+1, nzt
560       l_wall(k,:,:)  = l_grid(k)
561    ENDDO
562    l_wall(nzt+1,:,:) = l_grid(nzt)
563
564    ALLOCATE ( vertical_influence(nzb:nzt) )
565    DO  k = 1, nzt
566       vertical_influence(k) = MIN ( INT( l_grid(k) / &
567                     ( wall_adjustment_factor * dzw(k) ) + 0.5_wp ), nzt - k )
568    ENDDO
569
570    DO  k = 1, nzt
571       IF ( l_grid(k) > 1.5_wp * dx * wall_adjustment_factor .OR.  &
572            l_grid(k) > 1.5_wp * dy * wall_adjustment_factor )  THEN
573          WRITE( message_string, * ) 'grid anisotropy exceeds ', &
574                                     'threshold given by only local', &
575                                     ' &horizontal reduction of near_wall ', &
576                                     'mixing length l_wall', &
577                                     ' &starting from height level k = ', k, '.'
578          CALL message( 'init_grid', 'PA0202', 0, 1, 0, 6, 0 )
579          EXIT
580       ENDIF
581    ENDDO
582    vertical_influence(0) = vertical_influence(1)
583
584    DO  k = nzb + 1, nzb + vertical_influence(nzb)
585       l_wall(k,:,:) = zu(k) - zw(nzb)
586    ENDDO
587
588!
589!-- Set outer and inner index arrays for non-flat topography.
590!-- Here consistency checks concerning domain size and periodicity are
591!-- necessary.
592!-- Within this SELECT CASE structure only nzb_local is initialized
593!-- individually depending on the chosen topography type, all other index
594!-- arrays are initialized further below.
595    SELECT CASE ( TRIM( topography ) )
596
597       CASE ( 'flat' )
598!
599!--       nzb_local is required for the multigrid solver
600          nzb_local = 0
601
602       CASE ( 'single_building' )
603!
604!--       Single rectangular building, by default centered in the middle of the
605!--       total domain
606          blx = NINT( building_length_x / dx )
607          bly = NINT( building_length_y / dy )
608          bh  = MINLOC( ABS( zw - building_height ), 1 ) - 1
609          IF ( ABS( zw(bh  ) - building_height ) == &
610               ABS( zw(bh+1) - building_height )    )  bh = bh + 1
611
612          IF ( building_wall_left == 9999999.9_wp )  THEN
613             building_wall_left = ( nx + 1 - blx ) / 2 * dx
614          ENDIF
615          bxl = NINT( building_wall_left / dx )
616          bxr = bxl + blx
617
618          IF ( building_wall_south == 9999999.9_wp )  THEN
619             building_wall_south = ( ny + 1 - bly ) / 2 * dy
620          ENDIF
621          bys = NINT( building_wall_south / dy )
622          byn = bys + bly
623
624!
625!--       Building size has to meet some requirements
626          IF ( ( bxl < 1 ) .OR. ( bxr > nx-1 ) .OR. ( bxr < bxl+3 ) .OR.  &
627               ( bys < 1 ) .OR. ( byn > ny-1 ) .OR. ( byn < bys+3 ) )  THEN
628             WRITE( message_string, * ) 'inconsistent building parameters:',   &
629                                      '& bxl=', bxl, 'bxr=', bxr, 'bys=', bys, &
630                                      'byn=', byn, 'nx=', nx, 'ny=', ny
631             CALL message( 'init_grid', 'PA0203', 1, 2, 0, 6, 0 )
632          ENDIF
633
634!
635!--       Define the building.
636          nzb_local = 0
637          IF ( bxl <= nxr  .AND.  bxr >= nxl  .AND.                            &
638               bys <= nyn  .AND.  byn >= nys )                                 &       
639             nzb_local(MAX(nys,bys):MIN(nyn,byn),MAX(nxl,bxl):MIN(nxr,bxr)) = bh
640
641       CASE ( 'single_street_canyon' )
642!
643!--       Single quasi-2D street canyon of infinite length in x or y direction.
644!--       The canyon is centered in the other direction by default.
645          IF ( canyon_width_x /= 9999999.9_wp )  THEN
646!
647!--          Street canyon in y direction
648             cwx = NINT( canyon_width_x / dx )
649             IF ( canyon_wall_left == 9999999.9_wp )  THEN
650                canyon_wall_left = ( nx + 1 - cwx ) / 2 * dx
651             ENDIF
652             cxl = NINT( canyon_wall_left / dx )
653             cxr = cxl + cwx
654
655          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
656!
657!--          Street canyon in x direction
658             cwy = NINT( canyon_width_y / dy )
659             IF ( canyon_wall_south == 9999999.9_wp )  THEN
660                canyon_wall_south = ( ny + 1 - cwy ) / 2 * dy
661             ENDIF
662             cys = NINT( canyon_wall_south / dy )
663             cyn = cys + cwy
664
665          ELSE
666             
667             message_string = 'no street canyon width given'
668             CALL message( 'init_grid', 'PA0204', 1, 2, 0, 6, 0 )
669 
670          ENDIF
671
672          ch  = MINLOC( ABS( zw - canyon_height ), 1 ) - 1
673          IF ( ABS( zw(ch  ) - canyon_height ) == &
674               ABS( zw(ch+1) - canyon_height )    )  ch = ch + 1
675
676          dp_level_ind_b = ch
677!
678!--       Street canyon size has to meet some requirements
679          IF ( canyon_width_x /= 9999999.9_wp )  THEN
680             IF ( ( cxl < 1 ) .OR. ( cxr > nx-1 ) .OR. ( cwx < 3 ) .OR.        &
681               ( ch < 3 ) )  THEN
682                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
683                                           '&cxl=', cxl, 'cxr=', cxr,          &
684                                           'cwx=', cwx,                        &
685                                           'ch=', ch, 'nx=', nx, 'ny=', ny
686                CALL message( 'init_grid', 'PA0205', 1, 2, 0, 6, 0 ) 
687             ENDIF
688          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
689             IF ( ( cys < 1 ) .OR. ( cyn > ny-1 ) .OR. ( cwy < 3 ) .OR.        &
690               ( ch < 3 ) )  THEN
691                WRITE( message_string, * ) 'inconsistent canyon parameters:',  &
692                                           '&cys=', cys, 'cyn=', cyn,          &
693                                           'cwy=', cwy,                        &
694                                           'ch=', ch, 'nx=', nx, 'ny=', ny
695                CALL message( 'init_grid', 'PA0206', 1, 2, 0, 6, 0 ) 
696             ENDIF
697          ENDIF
698          IF ( canyon_width_x /= 9999999.9_wp .AND.                            &                 
699               canyon_width_y /= 9999999.9_wp )  THEN
700             message_string = 'inconsistent canyon parameters:' //             &   
701                              '&street canyon can only be oriented' //         &
702                              '&either in x- or in y-direction'
703             CALL message( 'init_grid', 'PA0207', 1, 2, 0, 6, 0 )
704          ENDIF
705
706          nzb_local = ch
707          IF ( canyon_width_x /= 9999999.9_wp )  THEN
708             IF ( cxl <= nxr  .AND.  cxr >= nxl )                              &
709                nzb_local(:,MAX(nxl,cxl+1):MIN(nxr,cxr-1)) = 0
710          ELSEIF ( canyon_width_y /= 9999999.9_wp )  THEN
711             IF ( cys <= nyn  .AND.  cyn >= nys )                              &         
712                nzb_local(MAX(nys,cys+1):MIN(nyn,cyn-1),:) = 0
713          ENDIF
714
715       CASE ( 'read_from_file' )
716
717          ALLOCATE ( topo_height(nys:nyn,nxl:nxr) )
718
719          DO  ii = 0, io_blocks-1
720             IF ( ii == io_group )  THEN
721
722!
723!--             Arbitrary irregular topography data in PALM format (exactly
724!--             matching the grid size and total domain size)
725                OPEN( 90, FILE='TOPOGRAPHY_DATA'//TRIM( coupling_char ),       &
726                          STATUS='OLD', FORM='FORMATTED', ERR=10 )
727!
728!--             Read topography PE-wise. Rows are read from nyn to nys, columns
729!--             are read from nxl to nxr. At first, ny-nyn rows need to be skipped.
730                skip_n_rows = 0
731                DO WHILE ( skip_n_rows < ny - nyn )
732                   READ( 90, * ) 
733                   skip_n_rows = skip_n_rows + 1
734                ENDDO
735!
736!--             Read data from nyn to nys and nxl to nxr. Therefore, skip
737!--             column until nxl-1 is reached
738                DO  j = nyn, nys, -1
739                   READ( 90, *, ERR=11, END=11 )                               &
740                                              ( dum, i = 0, nxl-1 ),           &
741                                              ( topo_height(j,i), i = nxl, nxr )
742                ENDDO
743
744                GOTO 12
745         
746 10             message_string = 'file TOPOGRAPHY'//TRIM( coupling_char )//    &
747                                 ' does not exist'
748                CALL message( 'init_grid', 'PA0208', 1, 2, 0, 6, 0 )
749
750 11             message_string = 'errors in file TOPOGRAPHY_DATA'//            &
751                                 TRIM( coupling_char )
752                CALL message( 'init_grid', 'PA0209', 1, 2, 0, 6, 0 )
753
754 12             CLOSE( 90 )
755
756             ENDIF
757#if defined( __parallel )
758             CALL MPI_BARRIER( comm2d, ierr )
759#endif
760          ENDDO
761
762!
763!--       Calculate the index height of the topography
764          nzb_local = 0
765          DO  i = nxl, nxr
766             DO  j = nys, nyn
767                nzb_local(j,i) = MINLOC( ABS( zw - topo_height(j,i) ), 1 ) - 1
768                IF ( ABS( zw(nzb_local(j,i)  ) - topo_height(j,i) ) == &
769                     ABS( zw(nzb_local(j,i)+1) - topo_height(j,i) )    )  &
770                   nzb_local(j,i) = nzb_local(j,i) + 1
771             ENDDO
772          ENDDO
773
774          DEALLOCATE ( topo_height )
775!
776!--       Filter topography, i.e. fill holes resolved by only one grid point. 
777!--       Such holes are suspected to lead to velocity blow-ups as continuity
778!--       equation on discrete grid cannot be fulfilled in such case.
779!--       For now, check only for holes and fill them to the lowest height level
780!--       of the directly adjoining grid points along x- and y- direction.
781!--       Before checking for holes, set lateral boundary conditions for
782!--       topography. After hole-filling, boundary conditions must be set again!
783          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
784         
785          IF ( .NOT. bc_ns_cyc )  THEN
786             IF ( nys == 0  )  nzb_local(-1,:)   = nzb_local(0,:)
787             IF ( nyn == ny )  nzb_local(ny+1,:) = nzb_local(ny,:)
788          ENDIF
789
790          IF ( .NOT. bc_lr_cyc )  THEN
791             IF ( nxl == 0  )  nzb_local(:,-1)   = nzb_local(:,0)
792             IF ( nxr == nx )  nzb_local(:,nx+1) = nzb_local(:,nx)         
793          ENDIF
794
795          num_hole_l = 0
796          DO i = nxl, nxr
797             DO j = nys, nyn
798
799                num_wall = 0
800
801                IF ( nzb_local(j-1,i) > nzb_local(j,i) )                       &
802                   num_wall = num_wall + 1
803                IF ( nzb_local(j+1,i) > nzb_local(j,i) )                       &
804                   num_wall = num_wall + 1
805                IF ( nzb_local(j,i-1) > nzb_local(j,i) )                       &
806                   num_wall = num_wall + 1
807                IF ( nzb_local(j,i+1) > nzb_local(j,i) )                       &
808                   num_wall = num_wall + 1
809
810                IF ( num_wall == 4 )  THEN
811                   nzb_local(j,i) = MIN( nzb_local(j-1,i), nzb_local(j+1,i),   &
812                                         nzb_local(j,i-1), nzb_local(j,i+1) )
813                   num_hole_l     = num_hole_l + 1
814                ENDIF
815             ENDDO
816          ENDDO
817!
818!--       Count the total number of holes, required for informative message.
819#if defined( __parallel )
820          CALL MPI_ALLREDUCE( num_hole_l, num_hole, 1, MPI_INTEGER, MPI_SUM,   &
821                              comm2d, ierr )
822#else
823          num_hole = num_hole_l
824#endif   
825!
826!--       Create an informative message if any hole was removed.
827          IF ( num_hole > 0 )  THEN
828             WRITE( message_string, * ) num_hole, 'hole(s) resolved by only '//&
829                                                  'one grid point were filled'
830             CALL message( 'init_grid', 'PA0430', 0, 0, 0, 6, 0 )
831          ENDIF
832!
833!--       Exchange ghost-points, as well as add cyclic or Neumann boundary
834!--       conditions.
835          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
836         
837          IF ( .NOT. bc_ns_cyc )  THEN
838             IF ( nys == 0  )  nzb_local(-1,:)   = nzb_local(0,:)
839             IF ( nyn == ny )  nzb_local(ny+1,:) = nzb_local(ny,:)
840          ENDIF
841
842          IF ( .NOT. bc_lr_cyc )  THEN
843             IF ( nxl == 0  )  nzb_local(:,-1)   = nzb_local(:,0)
844             IF ( nxr == nx )  nzb_local(:,nx+1) = nzb_local(:,nx)         
845          ENDIF
846
847       CASE DEFAULT
848!
849!--       The DEFAULT case is reached either if the parameter topography
850!--       contains a wrong character string or if the user has defined a special
851!--       case in the user interface. There, the subroutine user_init_grid
852!--       checks which of these two conditions applies.
853          CALL user_init_grid( nzb_local )
854
855    END SELECT
856!
857!-- Determine the maximum level of topography. Furthermore it is used for
858!-- steering the degradation of order of the applied advection scheme.
859!-- In case of non-cyclic lateral boundaries, the order of the advection
860!-- scheme has to be reduced up to nzt (required at the lateral boundaries).
861#if defined( __parallel )
862    CALL MPI_ALLREDUCE( MAXVAL( nzb_local ) + 1, nzb_max, 1, MPI_INTEGER,      &
863                        MPI_MAX, comm2d, ierr )
864#else
865    nzb_max = MAXVAL( nzb_local ) + 1
866#endif
867    IF ( inflow_l .OR. outflow_l .OR. inflow_r .OR. outflow_r .OR.             &
868         inflow_n .OR. outflow_n .OR. inflow_s .OR. outflow_s .OR.             &
869         nest_domain )                                                         &
870    THEN
871       nzb_max = nzt
872    ENDIF
873
874!
875!-- Consistency checks and index array initialization are only required for
876!-- non-flat topography, also the initialization of topography height arrays
877!-- zu_s_inner and zw_w_inner
878    IF ( TRIM( topography ) /= 'flat' )  THEN
879#if defined( __parallel )
880       CALL MPI_ALLREDUCE( MAXVAL( nzb_local ), nzb_local_max, 1, MPI_INTEGER, &
881                           MPI_MAX, comm2d, ierr )
882       CALL MPI_ALLREDUCE( MAXVAL( nzb_local ), nzb_local_min, 1, MPI_INTEGER, &
883                           MPI_MIN, comm2d, ierr )                           
884#else
885       nzb_local_max = MAXVAL( nzb_local )
886       nzb_local_min = MINVAL( nzb_local )
887#endif
888!
889!--    Consistency checks
890       IF ( nzb_local_min < 0  .OR.  nzb_local_max  > nz + 1 )  THEN
891          WRITE( message_string, * ) 'nzb_local values are outside the',       &
892                                'model domain',                                &
893                                '&MINVAL( nzb_local ) = ', nzb_local_min,      &
894                                '&MAXVAL( nzb_local ) = ', nzb_local_max
895          CALL message( 'init_grid', 'PA0210', 1, 2, 0, 6, 0 )
896       ENDIF
897
898       IF ( topography_grid_convention == 'cell_edge' )  THEN
899!
900!--       The array nzb_local as defined using the 'cell_edge' convention
901!--       describes the actual total size of topography which is defined at the
902!--       cell edges where u=0 on the topography walls in x-direction and v=0
903!--       on the topography walls in y-direction. However, PALM uses individual
904!--       arrays nzb_u|v|w|s_inner|outer that are based on nzb_s_inner.
905!--       Therefore, the extent of topography in nzb_local is now reduced by
906!--       1dx at the E topography walls and by 1dy at the N topography walls
907!--       to form the basis for nzb_s_inner.
908!--       Note, the reverse memory access (i-j instead of j-i) is absolutely
909!--       required at this point.
910          DO  j = nys+1, nyn+1
911             DO  i = nxl-1, nxr
912                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j,i+1) )
913             ENDDO
914          ENDDO
915!
916!--       Exchange ghost points
917          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
918
919          DO  i = nxl, nxr+1
920             DO  j = nys-1, nyn
921                nzb_local(j,i) = MIN( nzb_local(j,i), nzb_local(j+1,i) )
922             ENDDO
923          ENDDO
924!
925!--       Exchange ghost points         
926          CALL exchange_horiz_2d_int( nzb_local, nys, nyn, nxl, nxr, nbgp )
927       ENDIF
928!
929!--    Initialize index arrays nzb_s_inner and nzb_w_inner
930       nzb_s_inner = nzb_local
931       nzb_w_inner = nzb_local
932
933!
934!--    Initialize remaining index arrays:
935!--    first pre-initialize them with nzb_s_inner...
936       nzb_u_inner = nzb_s_inner
937       nzb_u_outer = nzb_s_inner
938       nzb_v_inner = nzb_s_inner
939       nzb_v_outer = nzb_s_inner
940       nzb_w_outer = nzb_s_inner
941       nzb_s_outer = nzb_s_inner
942
943!
944!--    ...then extend pre-initialized arrays in their according directions
945!--    based on nzb_local using nzb_tmp as a temporary global index array
946
947!
948!--    nzb_s_outer:
949!--    extend nzb_local east-/westwards first, then north-/southwards
950       nzb_tmp = nzb_local
951       DO  j = nys, nyn
952          DO  i = nxl, nxr
953             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i),             &
954                                 nzb_local(j,i+1) )
955          ENDDO
956       ENDDO
957       
958       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )
959       
960       DO  i = nxl, nxr
961          DO  j = nys, nyn
962             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
963                                     nzb_tmp(j+1,i) )
964          ENDDO
965!
966!--       non-cyclic boundary conditions (overwritten by call of
967!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
968          IF ( nys == 0 )  THEN
969             j = -1
970             nzb_s_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
971          ENDIF
972          IF ( nyn == ny )  THEN
973             j = ny + 1
974             nzb_s_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
975          ENDIF
976       ENDDO
977!
978!--    nzb_w_outer:
979!--    identical to nzb_s_outer
980       nzb_w_outer = nzb_s_outer
981
982!
983!--    nzb_u_inner:
984!--    extend nzb_local rightwards only
985       nzb_tmp = nzb_local
986       DO  j = nys, nyn
987          DO  i = nxl, nxr
988             nzb_tmp(j,i) = MAX( nzb_local(j,i-1), nzb_local(j,i) )
989          ENDDO
990       ENDDO
991       
992       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )
993       
994       nzb_u_inner = nzb_tmp
995!
996!--    nzb_u_outer:
997!--    extend current nzb_tmp (nzb_u_inner) north-/southwards
998       DO  i = nxl, nxr
999          DO  j = nys, nyn
1000             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i),             &
1001                                     nzb_tmp(j+1,i) )
1002          ENDDO
1003!
1004!--       non-cyclic boundary conditions (overwritten by call of
1005!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
1006          IF ( nys == 0 )  THEN
1007             j = -1
1008             nzb_u_outer(j,i) = MAX( nzb_tmp(j+1,i), nzb_tmp(j,i) )
1009          ENDIF
1010          IF ( nyn == ny )  THEN
1011             j = ny + 1
1012             nzb_u_outer(j,i) = MAX( nzb_tmp(j-1,i), nzb_tmp(j,i) )
1013          ENDIF
1014       ENDDO
1015
1016!
1017!--    nzb_v_inner:
1018!--    extend nzb_local northwards only
1019       nzb_tmp = nzb_local
1020       DO  i = nxl, nxr
1021          DO  j = nys, nyn
1022             nzb_tmp(j,i) = MAX( nzb_local(j-1,i), nzb_local(j,i) )
1023          ENDDO
1024       ENDDO
1025       
1026       CALL exchange_horiz_2d_int( nzb_tmp, nys, nyn, nxl, nxr, nbgp )     
1027       nzb_v_inner = nzb_tmp
1028
1029!
1030!--    nzb_v_outer:
1031!--    extend current nzb_tmp (nzb_v_inner) right-/leftwards
1032       DO  j = nys, nyn
1033          DO  i = nxl, nxr
1034             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i),             &
1035                                     nzb_tmp(j,i+1) )
1036          ENDDO
1037!
1038!--       non-cyclic boundary conditions (overwritten by call of
1039!--       exchange_horiz_2d_int below in case of cyclic boundary conditions)
1040          IF ( nxl == 0 )  THEN
1041             i = -1
1042             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i+1), nzb_tmp(j,i) )
1043          ENDIF
1044          IF ( nxr == nx )  THEN
1045             i = nx + 1
1046             nzb_v_outer(j,i) = MAX( nzb_tmp(j,i-1), nzb_tmp(j,i) )
1047          ENDIF
1048       ENDDO
1049
1050!
1051!--    Exchange of lateral boundary values (parallel computers) and cyclic
1052!--    boundary conditions, if applicable.
1053!--    Since nzb_s_inner and nzb_w_inner are derived directly from nzb_local
1054!--    they do not require exchange and are not included here.
1055       CALL exchange_horiz_2d_int( nzb_u_inner, nys, nyn, nxl, nxr, nbgp )
1056       CALL exchange_horiz_2d_int( nzb_u_outer, nys, nyn, nxl, nxr, nbgp )
1057       CALL exchange_horiz_2d_int( nzb_v_inner, nys, nyn, nxl, nxr, nbgp )
1058       CALL exchange_horiz_2d_int( nzb_v_outer, nys, nyn, nxl, nxr, nbgp )
1059       CALL exchange_horiz_2d_int( nzb_w_outer, nys, nyn, nxl, nxr, nbgp )
1060       CALL exchange_horiz_2d_int( nzb_s_outer, nys, nyn, nxl, nxr, nbgp )
1061
1062!
1063!--    Allocate and set the arrays containing the topography height
1064       ALLOCATE( zu_s_inner(0:nx+1,0:ny+1), zw_w_inner(0:nx+1,0:ny+1),         &
1065                 zu_s_inner_l(0:nx+1,0:ny+1), zw_w_inner_l(0:nx+1,0:ny+1) )
1066                 
1067       zu_s_inner   = 0.0_wp
1068       zw_w_inner   = 0.0_wp
1069       zu_s_inner_l = 0.0_wp
1070       zw_w_inner_l = 0.0_wp
1071       
1072       DO  i = nxl, nxr
1073          DO  j = nys, nyn
1074             zu_s_inner_l(i,j) = zu(nzb_local(j,i))
1075             zw_w_inner_l(i,j) = zw(nzb_local(j,i))
1076          ENDDO
1077       ENDDO
1078       
1079#if defined( __parallel )
1080       CALL MPI_REDUCE( zu_s_inner_l, zu_s_inner, (nx+2)*(ny+2),         &
1081                           MPI_REAL, MPI_SUM, 0, comm2d, ierr )       
1082       CALL MPI_REDUCE( zw_w_inner_l, zw_w_inner, (nx+2)*(ny+2),         &
1083                           MPI_REAL, MPI_SUM, 0, comm2d, ierr ) 
1084#else
1085       zu_s_inner = zu_s_inner_l
1086       zw_w_inner = zw_w_inner_l
1087#endif
1088
1089      DEALLOCATE( zu_s_inner_l, zw_w_inner_l )
1090      IF ( myid /= 0 )  DEALLOCATE( zu_s_inner, zw_w_inner )
1091!
1092!--   Set south and left ghost points, required for netcdf output
1093      IF ( myid == 0 )  THEN
1094         IF( bc_lr_cyc )  THEN
1095            zu_s_inner(nx+1,:) = zu_s_inner(0,:)
1096            zw_w_inner(nx+1,:) = zw_w_inner(0,:)
1097         ELSE
1098            zu_s_inner(nx+1,:) = zu_s_inner(nx,:)
1099            zw_w_inner(nx+1,:) = zw_w_inner(nx,:)
1100         ENDIF
1101         IF( bc_ns_cyc )  THEN
1102            zu_s_inner(:,ny+1) = zu_s_inner(:,0)
1103            zw_w_inner(:,ny+1) = zw_w_inner(:,0)
1104         ELSE
1105            zu_s_inner(:,ny+1) = zu_s_inner(:,ny)
1106            zw_w_inner(:,ny+1) = zw_w_inner(:,ny)
1107         ENDIF
1108      ENDIF
1109!
1110!--    Set flag arrays to be used for masking of grid points
1111       DO  i = nxlg, nxrg
1112          DO  j = nysg, nyng
1113             DO  k = nzb, nzt+1
1114                IF ( k <= nzb_s_inner(j,i) )  rflags_s_inner(k,j,i) = 0.0_wp
1115                IF ( k <= nzb_s_inner(j,i) )  rflags_invers(j,i,k)  = 0.0_wp
1116             ENDDO
1117          ENDDO
1118       ENDDO
1119
1120    ENDIF
1121!
1122!-- Deallocate temporary array, as it might be reused for different
1123!-- grid-levels further below.
1124    DEALLOCATE( nzb_tmp )
1125
1126!
1127!-- Set the individual index arrays which define the k index from which on
1128!-- the usual finite difference form (which does not use surface fluxes) is
1129!-- applied
1130    IF ( constant_flux_layer  .OR.  use_surface_fluxes )  THEN
1131       nzb_diff_u         = nzb_u_inner + 2
1132       nzb_diff_v         = nzb_v_inner + 2
1133       nzb_diff_s_inner   = nzb_s_inner + 2
1134       nzb_diff_s_outer   = nzb_s_outer + 2
1135    ELSE
1136       nzb_diff_u         = nzb_u_inner + 1
1137       nzb_diff_v         = nzb_v_inner + 1
1138       nzb_diff_s_inner   = nzb_s_inner + 1
1139       nzb_diff_s_outer   = nzb_s_outer + 1
1140    ENDIF
1141
1142!
1143!-- Calculation of wall switches and factors required by diffusion_u/v.f90 and
1144!-- for limitation of near-wall mixing length l_wall further below
1145    corner_nl = 0
1146    corner_nr = 0
1147    corner_sl = 0
1148    corner_sr = 0
1149    wall_l    = 0
1150    wall_n    = 0
1151    wall_r    = 0
1152    wall_s    = 0
1153
1154    DO  i = nxl, nxr
1155       DO  j = nys, nyn
1156!
1157!--       u-component
1158          IF ( nzb_u_outer(j,i) > nzb_u_outer(j+1,i) )  THEN
1159             wall_u(j,i) = 1.0_wp   ! north wall (location of adjacent fluid)
1160             fym(j,i)    = 0.0_wp
1161             fyp(j,i)    = 1.0_wp
1162          ELSEIF ( nzb_u_outer(j,i) > nzb_u_outer(j-1,i) )  THEN
1163             wall_u(j,i) = 1.0_wp   ! south wall (location of adjacent fluid)
1164             fym(j,i)    = 1.0_wp
1165             fyp(j,i)    = 0.0_wp
1166          ENDIF
1167!
1168!--       v-component
1169          IF ( nzb_v_outer(j,i) > nzb_v_outer(j,i+1) )  THEN
1170             wall_v(j,i) = 1.0_wp   ! rigth wall (location of adjacent fluid)
1171             fxm(j,i)    = 0.0_wp
1172             fxp(j,i)    = 1.0_wp
1173          ELSEIF ( nzb_v_outer(j,i) > nzb_v_outer(j,i-1) )  THEN
1174             wall_v(j,i) = 1.0_wp   ! left wall (location of adjacent fluid)
1175             fxm(j,i)    = 1.0_wp
1176             fxp(j,i)    = 0.0_wp
1177          ENDIF
1178!
1179!--       w-component, also used for scalars, separate arrays for shear
1180!--       production of tke
1181          IF ( nzb_w_outer(j,i) > nzb_w_outer(j+1,i) )  THEN
1182             wall_e_y(j,i) =  1.0_wp   ! north wall (location of adjacent fluid)
1183             wall_w_y(j,i) =  1.0_wp
1184             fwym(j,i)     =  0.0_wp
1185             fwyp(j,i)     =  1.0_wp
1186          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j-1,i) )  THEN
1187             wall_e_y(j,i) = -1.0_wp   ! south wall (location of adjacent fluid)
1188             wall_w_y(j,i) =  1.0_wp
1189             fwym(j,i)     =  1.0_wp
1190             fwyp(j,i)     =  0.0_wp
1191          ENDIF
1192          IF ( nzb_w_outer(j,i) > nzb_w_outer(j,i+1) )  THEN
1193             wall_e_x(j,i) =  1.0_wp   ! right wall (location of adjacent fluid)
1194             wall_w_x(j,i) =  1.0_wp
1195             fwxm(j,i)     =  0.0_wp
1196             fwxp(j,i)     =  1.0_wp
1197          ELSEIF ( nzb_w_outer(j,i) > nzb_w_outer(j,i-1) )  THEN
1198             wall_e_x(j,i) = -1.0_wp   ! left wall (location of adjacent fluid)
1199             wall_w_x(j,i) =  1.0_wp
1200             fwxm(j,i)     =  1.0_wp
1201             fwxp(j,i)     =  0.0_wp
1202          ENDIF
1203!
1204!--       Wall and corner locations inside buildings for limitation of
1205!--       near-wall mixing length l_wall
1206          IF ( nzb_s_inner(j,i) > nzb_s_inner(j+1,i) )  THEN
1207
1208             wall_n(j,i) = nzb_s_inner(j+1,i) + 1            ! North wall
1209
1210             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1211                corner_nl(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northleft corner
1212                                      nzb_s_inner(j,i-1) ) + 1
1213             ENDIF
1214
1215             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1216                corner_nr(j,i) = MAX( nzb_s_inner(j+1,i),  & ! Northright corner
1217                                      nzb_s_inner(j,i+1) ) + 1
1218             ENDIF
1219
1220          ENDIF
1221
1222          IF ( nzb_s_inner(j,i) > nzb_s_inner(j-1,i) )  THEN
1223
1224             wall_s(j,i) = nzb_s_inner(j-1,i) + 1            ! South wall
1225             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1226                corner_sl(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southleft corner
1227                                      nzb_s_inner(j,i-1) ) + 1
1228             ENDIF
1229
1230             IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1231                corner_sr(j,i) = MAX( nzb_s_inner(j-1,i),  & ! Southright corner
1232                                      nzb_s_inner(j,i+1) ) + 1
1233             ENDIF
1234
1235          ENDIF
1236
1237          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i-1) )  THEN
1238             wall_l(j,i) = nzb_s_inner(j,i-1) + 1            ! Left wall
1239          ENDIF
1240
1241          IF ( nzb_s_inner(j,i) > nzb_s_inner(j,i+1) )  THEN
1242             wall_r(j,i) = nzb_s_inner(j,i+1) + 1            ! Right wall
1243          ENDIF
1244
1245       ENDDO
1246    ENDDO
1247!
1248!-- Calculate wall flag arrays for the multigrid method.
1249!-- Please note, wall flags are only applied in the not cache-optimized
1250!-- version.
1251    IF ( psolver == 'multigrid_noopt' )  THEN
1252
1253!
1254!--    Gridpoint increment of the current level.
1255       inc = 1
1256       DO  l = maximum_grid_level, 1 , -1
1257!
1258!--       Set grid_level as it is required for exchange_horiz_2d_int
1259          grid_level = l
1260
1261          nxl_l = nxl_mg(l)
1262          nxr_l = nxr_mg(l)
1263          nys_l = nys_mg(l)
1264          nyn_l = nyn_mg(l)
1265          nzt_l = nzt_mg(l)
1266!
1267!--       Assign the flag level to be calculated
1268          SELECT CASE ( l )
1269             CASE ( 1 )
1270                flags => wall_flags_1
1271             CASE ( 2 )
1272                flags => wall_flags_2
1273             CASE ( 3 )
1274                flags => wall_flags_3
1275             CASE ( 4 )
1276                flags => wall_flags_4
1277             CASE ( 5 )
1278                flags => wall_flags_5
1279             CASE ( 6 )
1280                flags => wall_flags_6
1281             CASE ( 7 )
1282                flags => wall_flags_7
1283             CASE ( 8 )
1284                flags => wall_flags_8
1285             CASE ( 9 )
1286                flags => wall_flags_9
1287             CASE ( 10 )
1288                flags => wall_flags_10
1289          END SELECT
1290
1291!
1292!--       Depending on the grid level, set the respective bits in case of
1293!--       neighbouring walls
1294!--       Bit 0:  wall to the bottom
1295!--       Bit 1:  wall to the top (not realized in remaining PALM code so far)
1296!--       Bit 2:  wall to the south
1297!--       Bit 3:  wall to the north
1298!--       Bit 4:  wall to the left
1299!--       Bit 5:  wall to the right
1300!--       Bit 6:  inside building
1301
1302          flags = 0
1303
1304!
1305!--       In case of masking method, flags are not set and multigrid method
1306!--       works like FFT-solver
1307          IF ( .NOT. masking_method )  THEN
1308
1309!
1310!--          Allocate temporary array for topography heights on coarser grid
1311!--          level. Please note, 2 ghoist points are required, in order to
1312!--          calculate flags() on the interior ghost point.
1313             ALLOCATE( nzb_tmp(nys_l-2:nyn_l+2,nxl_l-2:nxr_l+2) )
1314             nzb_tmp = 0
1315             
1316             DO  i = nxl_l, nxr_l
1317                DO  j = nys_l, nyn_l
1318                   nzb_tmp(j,i) = nzb_local(j*inc,i*inc)
1319                ENDDO
1320             ENDDO
1321!
1322!--          Exchange ghost points on respective multigrid level. 2 ghost points
1323!--          are required, in order to calculate flags on
1324!--          nys_l-1 / nyn_l+1 / nxl_l-1 / nxr_l+1. The alternative would be to
1325!--          exchange 3D-INTEGER array flags on the respective multigrid level.
1326             CALL exchange_horiz_2d_int( nzb_tmp, nys_l, nyn_l, nxl_l, nxr_l, 2 )
1327!
1328!--          Set non-cyclic boundary conditions on respective multigrid level
1329             IF ( .NOT. bc_ns_cyc )  THEN
1330                IF ( nys == 0  )  THEN
1331                   nzb_tmp(-2,:) = nzb_tmp(0,:)
1332                   nzb_tmp(-1,:) = nzb_tmp(0,:)
1333                ENDIF
1334                IF ( nyn == ny )  THEN
1335                   nzb_tmp(ny+2,:) = nzb_tmp(ny,:)
1336                   nzb_tmp(ny+1,:) = nzb_tmp(ny,:)
1337                ENDIF
1338             ENDIF
1339             IF ( .NOT. bc_lr_cyc )  THEN
1340                IF ( nxl == 0  )  THEN
1341                   nzb_tmp(:,-2) = nzb_tmp(:,0)
1342                   nzb_tmp(:,-1) = nzb_tmp(:,0)
1343                ENDIF
1344                IF ( nxr == nx )  THEN
1345                   nzb_tmp(:,nx+1) = nzb_tmp(:,nx)   
1346                   nzb_tmp(:,nx+2) = nzb_tmp(:,nx)     
1347                ENDIF       
1348             ENDIF
1349                       
1350             DO  i = nxl_l-1, nxr_l+1
1351                DO  j = nys_l-1, nyn_l+1
1352                   DO  k = nzb, nzt_l+1     
1353!
1354!--                   Inside/outside building (inside building does not need
1355!--                   further tests for walls)
1356                      IF ( k*inc <= nzb_tmp(j,i) )  THEN
1357
1358                         flags(k,j,i) = IBSET( flags(k,j,i), 6 )
1359
1360                      ELSE
1361!
1362!--                      Bottom wall
1363                         IF ( (k-1)*inc <= nzb_tmp(j,i) )  THEN
1364                            flags(k,j,i) = IBSET( flags(k,j,i), 0 )
1365                         ENDIF
1366!
1367!--                      South wall
1368                         IF ( k*inc <= nzb_tmp(j-1,i) )  THEN
1369                            flags(k,j,i) = IBSET( flags(k,j,i), 2 )
1370                         ENDIF
1371!
1372!--                      North wall
1373                         IF ( k*inc <= nzb_tmp(j+1,i) )  THEN
1374                            flags(k,j,i) = IBSET( flags(k,j,i), 3 )
1375                         ENDIF
1376!
1377!--                      Left wall
1378                         IF ( k*inc <= nzb_tmp(j,i-1) )  THEN
1379                            flags(k,j,i) = IBSET( flags(k,j,i), 4 )
1380                         ENDIF
1381!
1382!--                      Right wall
1383                         IF ( k*inc <= nzb_tmp(j,i+1) )  THEN
1384                            flags(k,j,i) = IBSET( flags(k,j,i), 5 )
1385                         ENDIF
1386
1387                      ENDIF
1388                           
1389                   ENDDO
1390                ENDDO
1391             ENDDO
1392
1393             DEALLOCATE( nzb_tmp )
1394
1395          ENDIF
1396
1397          inc = inc * 2
1398
1399       ENDDO
1400!
1401!--    Reset grid_level to "normal" grid
1402       grid_level = 0
1403       
1404    ENDIF
1405!
1406!-- Allocate flags needed for masking walls. Even though these flags are only
1407!-- required in the ws-scheme, the arrays need to be allocated as they are
1408!-- used in OpenACC directives.
1409    ALLOCATE( wall_flags_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
1410              wall_flags_00(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
1411    wall_flags_0  = 0
1412    wall_flags_00 = 0
1413!
1414!-- Init flags for ws-scheme to degrade order near walls
1415    IF ( momentum_advec == 'ws-scheme'  .OR.  scalar_advec == 'ws-scheme'  .OR.&
1416         scalar_advec   == 'ws-scheme-mono' )  THEN
1417       CALL ws_init_flags
1418    ENDIF
1419
1420!
1421!-- In case of topography: limit near-wall mixing length l_wall further:
1422!-- Go through all points of the subdomain one by one and look for the closest
1423!-- surface
1424    IF ( TRIM(topography) /= 'flat' )  THEN
1425       DO  i = nxl, nxr
1426          DO  j = nys, nyn
1427
1428             nzb_si = nzb_s_inner(j,i)
1429             vi     = vertical_influence(nzb_si)
1430
1431             IF ( wall_n(j,i) > 0 )  THEN
1432!
1433!--             North wall (y distance)
1434                DO  k = wall_n(j,i), nzb_si
1435                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i), 0.5_wp * dy )
1436                ENDDO
1437!
1438!--             Above North wall (yz distance)
1439                DO  k = nzb_si + 1, nzb_si + vi
1440                   l_wall(k,j+1,i) = MIN( l_wall(k,j+1,i),                     &
1441                                          SQRT( 0.25_wp * dy**2 +              &
1442                                          ( zu(k) - zw(nzb_si) )**2 ) )
1443                ENDDO
1444!
1445!--             Northleft corner (xy distance)
1446                IF ( corner_nl(j,i) > 0 )  THEN
1447                   DO  k = corner_nl(j,i), nzb_si
1448                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1), &
1449                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1450                   ENDDO
1451!
1452!--                Above Northleft corner (xyz distance)
1453                   DO  k = nzb_si + 1, nzb_si + vi
1454                      l_wall(k,j+1,i-1) = MIN( l_wall(k,j+1,i-1),              &
1455                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1456                                            ( zu(k) - zw(nzb_si) )**2 ) )
1457                   ENDDO
1458                ENDIF
1459!
1460!--             Northright corner (xy distance)
1461                IF ( corner_nr(j,i) > 0 )  THEN
1462                   DO  k = corner_nr(j,i), nzb_si
1463                       l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),             &
1464                                                0.5_wp * SQRT( dx**2 + dy**2 ) )
1465                   ENDDO
1466!
1467!--                Above northright corner (xyz distance)
1468                   DO  k = nzb_si + 1, nzb_si + vi
1469                      l_wall(k,j+1,i+1) = MIN( l_wall(k,j+1,i+1),              &
1470                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1471                                            ( zu(k) - zw(nzb_si) )**2 ) )
1472                   ENDDO
1473                ENDIF
1474             ENDIF
1475
1476             IF ( wall_s(j,i) > 0 )  THEN
1477!
1478!--             South wall (y distance)
1479                DO  k = wall_s(j,i), nzb_si
1480                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i), 0.5_wp * dy )
1481                ENDDO
1482!
1483!--             Above south wall (yz distance)
1484                DO  k = nzb_si + 1, nzb_si + vi
1485                   l_wall(k,j-1,i) = MIN( l_wall(k,j-1,i),                     &
1486                                          SQRT( 0.25_wp * dy**2 +              &
1487                                          ( zu(k) - zw(nzb_si) )**2 ) )
1488                ENDDO
1489!
1490!--             Southleft corner (xy distance)
1491                IF ( corner_sl(j,i) > 0 )  THEN
1492                   DO  k = corner_sl(j,i), nzb_si
1493                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1494                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1495                   ENDDO
1496!
1497!--                Above southleft corner (xyz distance)
1498                   DO  k = nzb_si + 1, nzb_si + vi
1499                      l_wall(k,j-1,i-1) = MIN( l_wall(k,j-1,i-1),              &
1500                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1501                                            ( zu(k) - zw(nzb_si) )**2 ) )
1502                   ENDDO
1503                ENDIF
1504!
1505!--             Southright corner (xy distance)
1506                IF ( corner_sr(j,i) > 0 )  THEN
1507                   DO  k = corner_sr(j,i), nzb_si
1508                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1509                                               0.5_wp * SQRT( dx**2 + dy**2 ) )
1510                   ENDDO
1511!
1512!--                Above southright corner (xyz distance)
1513                   DO  k = nzb_si + 1, nzb_si + vi
1514                      l_wall(k,j-1,i+1) = MIN( l_wall(k,j-1,i+1),              &
1515                                            SQRT( 0.25_wp * (dx**2 + dy**2) +  &
1516                                            ( zu(k) - zw(nzb_si) )**2 ) )
1517                   ENDDO
1518                ENDIF
1519
1520             ENDIF
1521
1522             IF ( wall_l(j,i) > 0 )  THEN
1523!
1524!--             Left wall (x distance)
1525                DO  k = wall_l(j,i), nzb_si
1526                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1), 0.5_wp * dx )
1527                ENDDO
1528!
1529!--             Above left wall (xz distance)
1530                DO  k = nzb_si + 1, nzb_si + vi
1531                   l_wall(k,j,i-1) = MIN( l_wall(k,j,i-1),                     &
1532                                       SQRT( 0.25_wp * dx**2 +                 &
1533                                       ( zu(k) - zw(nzb_si) )**2 ) )
1534                ENDDO
1535             ENDIF
1536
1537             IF ( wall_r(j,i) > 0 )  THEN
1538!
1539!--             Right wall (x distance)
1540                DO  k = wall_r(j,i), nzb_si
1541                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1), 0.5_wp * dx )
1542                ENDDO
1543!
1544!--             Above right wall (xz distance)
1545                DO  k = nzb_si + 1, nzb_si + vi
1546                   l_wall(k,j,i+1) = MIN( l_wall(k,j,i+1),                     &
1547                                          SQRT( 0.25_wp * dx**2 +              &
1548                                          ( zu(k) - zw(nzb_si) )**2 ) )
1549                ENDDO
1550
1551             ENDIF
1552
1553          ENDDO
1554       ENDDO
1555
1556    ENDIF
1557
1558!
1559!-- Multiplication with wall_adjustment_factor
1560    l_wall = wall_adjustment_factor * l_wall
1561
1562!
1563!-- Set lateral boundary conditions for l_wall
1564    CALL exchange_horiz( l_wall, nbgp )
1565
1566    DEALLOCATE( corner_nl, corner_nr, corner_sl, corner_sr, nzb_local, &
1567                vertical_influence, wall_l, wall_n, wall_r, wall_s )
1568
1569
1570 END SUBROUTINE init_grid
Note: See TracBrowser for help on using the repository browser.