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

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

last commit documented

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