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

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

last commit documented

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