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

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

last commit documented

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