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

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

PE-wise reading of topography file

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