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

Last change on this file since 1801 was 1780, checked in by raasch, 9 years ago

last commit documented

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