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

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

last commit documented

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